API documentation of kitaev_clusters

symmetry_functions module

kitaev_clusters.symmetry_functions.get_neighbors(Lx, Ly, bond)

For a given bond-direction in the honeycomb lattice (x, y, or z), this finds all the pairs of sites along this type of bond.

Parameters:
  • Lx (int) – Number of units cells in the x-direction.

  • Ly (int) – Number of units cells in the y-direction.

  • bond (str) – A string which is either “x”, “y”, or “z” corresponding to the desired honeycomb bond direction.

Returns:

neighbors – A list of tuples of site pairs which each constitute a bond of the chosen type.

Return type:

list

kitaev_clusters.symmetry_functions.get_representative_states(L, translation_order_list, ints_filename='kept_ints', states_filename='state_map', nunique_filename='n_unique_list')

For a lattice with L sites and therefore 3^L possible configurations, this function finds all the representative states we need to keep, thus reducing the Hilbert space dimension from 3^L to a much smaller number ndim.

Parameters:
  • L (int) – The total number of lattice sites.

  • translation_order_list (list) – A list of translationally equivalent site orderings. This is returned by the function lattice_translations.

  • ints_filename (str) – Filename for the .npy file created which stores the ndim representative states (in base-10).

  • states_filename (str) – Filename for the .npy file created which stores the representative states (in base-10) for each of the 3^L possible lattice configurations, i.e. a mapping from a state to its representative state.

  • nunique_filename (str) – Filename for the .npy file created which stores the number of unique mirror states for each representative state.

Returns:

  • kept_ints (ndarray) – An array of the base-10 integers corresponding to the representative states.

  • state_map (ndarray) – An array containing the representative state (in base-10) for each of the 3^L possible states.

  • n_unique_list (ndarray) – An array containing the number of unique mirror states for each representative state.

  • ndim (int) – The total number of representative states we keep, i.e. the Hilbert space dimension is reduced from 3^L to ndim.

kitaev_clusters.symmetry_functions.lattice_translations(Lx, Ly)

For a specified honeycomb lattice (with Lx unit cells in the x-direction, and Ly units cells in y-direction), returns a list of site orderings which are equivalent via translational symmetry (with periodic boundary conditions). For example, with Lx=3 and Ly=2 there are 6 unit cells and L=12 total sites, thus for any ordered labelling of sites, there are 5 equivalent orders which are generated by the translations (1,0), (2,0), (0,1), (1,1), and (2,1).

Parameters:
  • Lx (int) – Number of units cells in the x-direction.

  • Ly (int) – Number of unit cells in the y-direction.

Returns:

translation_order_list – A list of all equivalent site orderings generated via non-zero translations (with periodic boundary conditions).

Return type:

list

kitaev_clusters.symmetry_functions.mirror_states(state, translation_order_list)

For a single state represented by a ternary string, this functions finds all equivalent ‘mirror’ states related by symmetry, and their base-10 integer representation. The state with the smallest base-10 integer representation is chosen to be the ‘representative state’. The number of unique mirror states for a given state, n_unique, is also found.

Parameters:
  • state (str) – A state represented by a ternary string.

  • translation_order_list (list) – A list of translationally equivalent site orderings. This is returned by the function lattice_translations.

Returns:

  • sorted_ints (list) – A list of the base-10 integers (in ascending order) representing the equivalent mirror states.

  • n_unique (int) – The number of unique mirror states.

  • rep_state (str) – A ternary string corresponding to the state’s representative state.

  • rep_int (int) – A base-10 integer corresponding to the state’s representative state.

kitaev_clusters.symmetry_functions.tern_to_base10(state)

Converts a ternary string to a base-10 integer.

Parameters:

state (str) – A string representing a number in ternary.

Returns:

integer_rep – The base-10 integer corresponding to the ternary string.

Return type:

int

kitaev_clusters.symmetry_functions.ternary(n)

Converts a base-10 integer to a ternary string.

Parameters:

n (int) – A base-10 integer.

Returns:

A string representing an integer in ternary.

Return type:

str

kitaev_clusters.symmetry_functions.ternary_pad(n, L)

Converts a base-10 integer to a ternary string of length L. Will left-pad ternary string with zeros such that the string has length L.

Parameters:
  • n (int) – A base-10 integer.

  • L (int) – Desired total length of ternary string. The string will be left-padded with zeros to have length L.

Returns:

A string of length L representing an integer in ternary.

Return type:

str

hamiltonian_functions module

kitaev_clusters.hamiltonian_functions.SmSm(state, j, k)

Applies spin lowering operator S- to sites j and k, i.e. acts S-(j) S-(k) on a state.

Parameters:
  • state (str) – A string representing a state in ternary.

  • j (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site j to lower.

  • k (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site k to lower.

Returns:

new_state – The ternary representation of the new state after the spin lowering operator S- has acted on sites j and k.

Return type:

str

kitaev_clusters.hamiltonian_functions.SmSp(state, j, k)

Applies spin lowering operator S- to site j, and spin raising operator S+ to site k, i.e. acts S-(j) S+(k) on a state.

Parameters:
  • state (str) – A string representing a state in ternary.

  • j (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site j to lower.

  • k (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site k to raise.

Returns:

new_state – The ternary representation of the new state after the spin lowering operator S- has acted on sites j, and the spin raising operator S+ has acted on site k.

Return type:

str

kitaev_clusters.hamiltonian_functions.SpSm(state, j, k)

Applies spin raising operator S+ to site j, and spin lowering operator S- to site k, i.e. acts S+(j) S-(k) on a state.

Parameters:
  • state (str) – A string representing a state in ternary.

  • j (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site j to raise.

  • k (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site k to lower.

Returns:

new_state – The ternary representation of the new state after the spin raising operator S+ has acted on site j, and the spin lowering operator S- has acted on site k.

Return type:

str

kitaev_clusters.hamiltonian_functions.SpSp(state, j, k)

Applies spin raising operator S+ to sites j and k, i.e. acts S+(j) S+(k) on a state.

Parameters:
  • state (str) – A string representing a state in ternary.

  • j (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site j to raise.

  • k (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site k to raise.

Returns:

new_state – The ternary representation of the new state after the spin raising operator S+ has acted on sites j and k.

Return type:

str

kitaev_clusters.hamiltonian_functions.Ss(state, i)

Applies squared spin operator (Sx + Sy + Sz)^2 to a site j, and returns the resultant states and their coefficients.

Parameters:
  • state (str) – A string representing a state in ternary.

  • i (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site i where we apply the operator (Sx + Sy + Sz)^2.

Returns:

A list of tuples of the form (state, coefficient), which result from applying the operator (Sx + Sy + Sz)^2 to our initial state at site i.

Return type:

list

kitaev_clusters.hamiltonian_functions.SzSz(state, j, k)

Applies Sz operator to sites j and k, and returns the coefficient multiplying the resultant (identical) state.

Parameters:
  • state (str) – A string representing a state in ternary.

  • j (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site j where we apply the Sz operator.

  • k (int) – An integer in the range [0, L-1] (where L is the total number of lattice sites) denoting a site k where we apply the Sz operator.

Returns:

new_state – The ternary representation of the new state after the Sz operator has acted on sites j and k.

Return type:

str

kitaev_clusters.hamiltonian_functions.construct_HD(L, kept_ints, state_map, n_unique_list, filename, as_csr=True)

Constructs HD, i.e. the component of the Hamiltonian which describes a single-ion anisotropy in the [1, 1, 1] direction. The HD operator involves a sum over all sites of D * (Sx + Sy + Sz)^2, where D is the magnitude of the single-ion anisotropy.

Parameters:
  • L (int) – The total number of lattice sites.

  • kept_ints (ndarray) – An array of the base-10 integers corresponding to the representative states. This is the first value returned by the function get_representative_states.

  • state_map (ndarray) – An array containing the representative state (in base-10) for each of the 3^L possible states. This is the second value returned by the function get_representative_states.

  • n_unique_list (ndarray) – An array containing the number of unique mirror states for each representative state. This is the third value returned by the function get_representative_states.

  • filename (str) – Filename for the HD matrix which will be saved as a .npz file.

  • as_csr (bool, optional) – If True (default), the HD sparse matrix is converted from List of Lists (LIL) format to Compressed Sparse Row (CSR) format before saving.

Returns:

HD – The matrix HD in Compressed Sparse Row (CSR) format, i.e. the component of the Hamiltonian which describes a single-ion anisotropy in the [1, 1, 1] direction.

Return type:

csr_matrix

kitaev_clusters.hamiltonian_functions.construct_Hxx(L, x_neighbors, kept_ints, state_map, n_unique_list, filename, as_csr=True)

Constructs Hxx, i.e. the component of the Kitaev Hamiltonian which sums over all x-direction bonds of the honeycomb lattice.

Parameters:
  • L (int) – The total number of lattice sites.

  • x_neighbors (list) – A list of tuples of site pairs which each constitute an x-bond. This is returned by the function get_neighbors with the parameter bond=”x”.

  • kept_ints (ndarray) – An array of the base-10 integers corresponding to the representative states. This is the first value returned by the function get_representative_states.

  • state_map (ndarray) – An array containing the representative state (in base-10) for each of the 3^L possible states. This is the second value returned by the function get_representative_states.

  • n_unique_list (ndarray) – An array containing the number of unique mirror states for each representative state. This is the third value returned by the function get_representative_states.

  • filename (str) – Filename for the Hxx matrix which will be saved as a .npz file.

  • as_csr (bool, optional) – If True (default), the Hxx sparse matrix is converted from List of Lists (LIL) format to Compressed Sparse Row (CSR) format before saving.

Returns:

Hxx – The matrix Hxx in Compressed Sparse Row (CSR) format, i.e. the component of the Kitaev Hamiltonian which sums over all x-direction bonds of the honeycomb lattice. It will have dimensions ndim x ndim, where ndim is the number of representative states.

Return type:

csr_matrix

kitaev_clusters.hamiltonian_functions.construct_Hyy(L, y_neighbors, kept_ints, state_map, n_unique_list, filename, as_csr=True)

Constructs Hyy, i.e. the component of the Kitaev Hamiltonian which sums over all y-direction bonds of the honeycomb lattice.

Parameters:
  • L (int) – The total number of lattice sites.

  • y_neighbors (list) – A list of tuples of site pairs which each constitute a y-bond. This is returned by the function get_neighbors with the parameter bond=”y”.

  • kept_ints (ndarray) – An array of the base-10 integers corresponding to the representative states. This is the first value returned by the function get_representative_states.

  • state_map (ndarray) – An array containing the representative state (in base-10) for each of the 3^L possible states. This is the second value returned by the function get_representative_states.

  • n_unique_list (ndarray) – An array containing the number of unique mirror states for each representative state. This is the third value returned by the function get_representative_states.

  • filename (str) – Filename for the Hyy matrix which will be saved as a .npz file.

  • as_csr (bool, optional) – If True (default), the Hxx sparse matrix is converted from List of Lists (LIL) format to Compressed Sparse Row (CSR) format before saving.

Returns:

Hyy – The matrix Hxx in Compressed Sparse Row (CSR) format, i.e. the component of the Kitaev Hamiltonian which sums over all y-direction bonds of the honeycomb lattice. It will have dimensions ndim x ndim, where ndim is the number of representative states.

Return type:

csr_matrix

kitaev_clusters.hamiltonian_functions.construct_Hzz(L, z_neighbors, kept_ints, filename, as_csr=True)

Constructs Hzz, i.e. the component of the Kitaev Hamiltonian which sums over all z-direction bonds of the honeycomb lattice.

Parameters:
  • L (int) – The total number of lattice sites.

  • z_neighbors (list) – A list of tuples of site pairs which each constitute a z-bond. This is returned by the function get_neighbors with the parameter bond=”z”.

  • kept_ints (ndarray) – An array of the base-10 integers corresponding to the representative states. This is the first value returned by the function get_representative_states.

  • filename (str) – Filename for the Hzz matrix which will be saved as a .npz file.

  • as_csr (bool, optional) – If True (default), the Hzz sparse matrix is converted from List of Lists (LIL) format to Compressed Sparse Row (CSR) format before saving.

Returns:

Hzz – The matrix Hxx in Compressed Sparse Row (CSR) format, i.e. the component of the Kitaev Hamiltonian which sums over all z-direction bonds of the honeycomb lattice. It will have dimensions ndim x ndim, where ndim is the number of representative states.

Return type:

csr_matrix

kitaev_clusters.hamiltonian_functions.construct_hamiltonian(Kx, Ky, Kz, D, Hxx_file, Hyy_file, Hzz_file, HD_file)

Loads from file separate components of the Hamiltonian (Hxx, Hyy, Hzz, HD) and constructs the full Hamiltonian matrix.

Parameters:
  • Kx (float or int) – The Kitaev coupling along x-bonds.

  • Ky (float or int) – The Kitaev coupling along y-bonds.

  • Kz (float or int) – The Kitaev coupling along z-bonds.

  • D (float or int) – The magnitude of the single-ion anisotropy term, i.e. the coefficient which multiplies HD.

  • Hxx_file (.npz file) – Saved file (in .npz format) for the Hxx matrix.

  • Hyy_file (.npz file) – Saved file (in .npz format) for the Hyy matrix.

  • Hzz_file (.npz file) – Saved file (in .npz format) for the Hzz matrix.

  • HD_file (.npz file) – Saved file (in .npz format) for the HD matrix.

Returns:

H – The full sparse matrix Kitaev Hamiltonian corresponding to H = (Kx * Hxx) + (Ky * Hyy) + (Kz * Hzz) + (D * HD).

Return type:

csr_matrix

ground_state_functions module

kitaev_clusters.ground_state_functions.entanglement_entropy(L, psi, state_map, n_unique_list)

Finds the bipartite entanglement entropy when the system is divided into two equal halves. Performs a singular value decomposition (SVD) to obtain the singular values.

Parameters:
  • L (int) – The total number of lattice sites. Used to calculate <(Sx + Sy + Sz)^2> per site.

  • psi (ndarray) – The state of the system, i.e. a 1d array with ndim complex coefficients. Typically the ground state psi_0.

  • state_map (ndarray) – An array containing the representative state (in base-10) for each of the 3^L possible states. This is the second value returned by the function get_representative_states.

  • n_unique_list (ndarray) – An array containing the number of unique mirror states for each representative state. This is the third value returned by the function get_representative_states.

Returns:

entropy – The bipartite entanglement entropy when the system is divided into two equal halves.

Return type:

float

kitaev_clusters.ground_state_functions.ground_state(H, L, k=1, ncv=20)

Finds the ground state psi_0 and the ground state energy E_0, by applying Lanczos algorithm to matrix Hamiltonian H.

Parameters:
  • H (csr_matrix or lil_matrix) – The full sparse matrix Hamiltonian.

  • L (int) – The total number of lattice sites. Used to calculate the energy per site.

  • k (int, optional) – Number of eigenvalues and eigenvectors desired. By default k=1 to find state only.

  • ncv (int, optional) – Number of Lanczos vectors generated.

Returns:

  • E_0 (float) – The ground state energy per site.

  • psi_0 (ndarray) – The ground state, i.e. a 1d array with ndim complex coefficients.

kitaev_clusters.ground_state_functions.s_squared(HD, L, psi_0)

Finds the expectation value of (Sx + Sy + Sz)^2 in the ground state (per site).

Parameters:
  • HD (csr_matrix or lil_matrix) – The sparse matrix HD corresponding to the single-ion anisotropy term in the Hamiltonian.

  • L (int) – The total number of lattice sites. Used to calculate <(Sx + Sy + Sz)^2> per site.

  • psi_0 (ndarray) – The ground state, i.e. a 1d array with ndim complex coefficients.

Returns:

Ss – The expectation value of (Sx + Sy + Sz)^2 in the ground state (per site).

Return type:

float