Description of the functions
Submodules
pybott.bott module
Code to compute the Bott index following the definition given by T. A. Loring and M. B. Hastings in https://iopscience.iop.org/article/10.1209/0295-5075/92/67004/meta
The Bott index measures the commutativity of projected position operators, providing a topological invariant that helps distinguish topological insulators from trivial insulators.
- pybott.bott.all_bott(lattice, ham, orb=1, dagger=False, energy_max=0)
Compute the Bott index for a given Hamiltonian and lattice for all energy levels or up to a specified limit.
This function calculates the Bott index for each energy state in the system, from the lowest to the highest energy state, unless a stopping point is specified via the energy_max parameter.
The Bott index is computed for each eigenstate, and its evolution can be tracked across the energy spectrum of the system.
- Parameters:
lattice (ndarray) – Array of shape
(N_sites, 2)
containing the coordinates of the lattice sites.ham (ndarray) – Hamiltonian matrix of shape
(orb * N_sites, orb * N_sites)
. The Hamiltonian must be Hermitian.orb (int, optional) – Number of orbitals considered per lattice site. Default is
1
.dagger (bool, optional) – If
True
, computes the Bott index using the Hermitian conjugate (dagger) of the projected position operators. IfFalse
, computes using the inverse of the position operators. Default isFalse
.energy_max (float, optional) – The maximum energy to consider. If not
0
, the function will compute the Bott index only for eigenstates with energy less thanenergy_max
. Default is0
, meaning the function computes the Bott index for all energy levels.
- Returns:
A dictionary where the keys are the energy values and the values are the corresponding Bott index calculated for each energy level.
- Return type:
dict
- Raises:
ValueError – If the Hamiltonian is not Hermitian.
Note
The function iterates over all the eigenstates (or up to the specified limit) and computes the Bott index for each state. This allows one to track the evolution of the topological properties of the system across its entire energy spectrum. This is particularly useful in systems with energy-dependent topological transitions.
- pybott.bott.all_bott_vect(lattice, evects, energies, orb=1, dagger=False, energy_max=inf)
Compute the Bott index for all energy levels or up to a specified limit.
This function calculates the Bott index for each energy state in the system, sequentially from the lowest to the highest energy state, unless a stopping point is specified via the energy_max parameter.
- Parameters:
lattice (ndarray) – Array of shape
(N_sites, 2)
containing the coordinates of the lattice sites.evects (ndarray) – Array of shape
(orb * N_sites, orb * N_sites)
containing the eigenvectors of the system.energies (ndarray) – Array of shape
(orb * N_sites,)
containing the energy values corresponding to the eigenstates. These energies may differ from the eigenvalues of the Hamiltonian for more complex systems.orb (int, optional) – Number of orbitals considered per lattice site. Default is
1
.dagger (bool, optional) – If
True
, computes the Bott index using the Hermitian conjugate (dagger) of the projected position operators. IfFalse
, computes using the inverse of the position operators. Default isFalse
.energy_max (float, optional) – The maximum energy to consider. If not
np.inf
, the calculation will only be performed for eigenstates with energy less thanenergy_max
. Default isnp.inf
, meaning the function computes the Bott index for all energy levels.
- Returns:
A dictionary where the keys are the energy values and the values are the corresponding Bott index calculated for each energy level.
- Return type:
dict
- Raises:
ValueError – If the Hamiltonian is not Hermitian.
Note
The function iterates over all the eigenstates (or up to the specified limit) and computes the Bott index for each state. This allows one to track the evolution of the topological properties of the system across its entire energy spectrum. This is particularly useful in systems with energy-dependent topological transitions.
- pybott.bott.bott(lattice, ham, fermi_energy=0, gap=None, orb=1, dagger=False)
Calculate the Bott index of a system described by a given Hamiltonian and lattice.
This function calculates the Bott index, which is a topological invariant used to distinguish topological phases in a system described by the Hamiltonian. If the Hamiltonian is not Hermitian, you must compute the eigenvalues and eigenvectors independently and use bott_vect instead. If the theoretical width of the gap is provided and the Hamiltonian is large, eigenvalues and eigenvectors will be computed in a restricted Hilbert space to save computation time.
- Parameters:
lattice (ndarray) – Array of shape
(N_sites, 2)
containing the coordinates of the lattice sites.ham (ndarray) – Hamiltonian matrix of shape
(orb * N_sites, orb * N_sites)
. Must be Hermitian.fermi_energy (float, optional) – Value of energy for which the Bott index is computed, must be in the bulk gap to match the Chern number. Not defined outside of the bulk gap but usually gives 0. Optional.
gap (tuple of float, optional) – Energy gap used for filtering eigenvalues when calculating the Bott index. Must be a tuple of two ordered real numbers. If
None
, the entire spectrum is computed. Optional.orb (int, optional) – Number of orbitals considered per lattice site. Default is
1
.dagger (bool) – Specifies the method to compute the Bott index. If
True
, uses the dagger of the projected position operator; otherwise, it computes the inverse of the operator.
- Returns:
The computed Bott index.
- Return type:
float
- Raises:
ValueError – If the Hamiltonian is not Hermitian, or if gap is not a valid tuple of floats.
- pybott.bott.bott_matrix(u_mat, v_mat, dagger=False)
This function computes the Bott index for two invertible matrices, U and V. The Bott index is a topological invariant used to distinguish different topological phases. The function either computes the standard Bott index or uses the dagger of the projected position operator depending on the value of the dagger parameter.
The Bott index is mathematically defined as:
\[\text{Bott}(U, V) = \frac{1}{2 \pi i} \text{Tr} \log (UVU^{- 1} V^{- 1})\]- Parameters:
u_mat (ndarray) – The matrix
U
.v_mat (ndarray) – The matrix
V
.dagger (bool, optional) – If
True
, the method uses the conjugate transpose (dagger) of the matrices in the Bott index computation. Default isFalse
, in which case the inverse matrices are used.
- Returns:
The computed Bott index.
- Return type:
float
- Raises:
np.linalg.LinAlgError – If either of the matrices
U
orV
is not invertible.
- pybott.bott.bott_vect(lattice, evects, energies, fermi_energy=0, orb=1, dagger=False)
Compute the Bott index for a given set of eigenvectors and energies.
This function computes the Bott index, which is useful when the Hamiltonian is not Hermitian, and there is a need for additional preparation of the eigenvalues and eigenvectors before sending them to the Bott routine. This can be particularly useful in systems beyond tight-binding models. See the example on topological photonic systems provided in the documentation.
- Parameters:
lattice (ndarray) – Array of shape
(N_sites, 2)
containing the coordinates of the lattice sites.evects (ndarray) – Array of shape
(orb * N_sites, orb * N_sites)
containing the eigenvectors.energies (ndarray) – Array of shape
(orb * N_sites,)
containing the energies. These energies may differ from the eigenvalues of the Hamiltonian in more complex systems beyond tight-binding models.fermi_energy (float, optional) – Energy value at which the Bott index is computed. It must be within the bulk gap to match the Chern number. Outside of the bulk gap, it usually returns 0. Default is
0
.orb (int, optional) – Number of orbitals considered per lattice site. Default is
1
.dagger (bool, optional) – Specifies which method to use for computing the Bott index. If
True
, the method uses the dagger of the projected position operator; otherwise, it computes the inverse of the operator. Default isFalse
.
- Returns:
The Bott index value.
- Return type:
float
- pybott.bott.compute_uv(lattice, eigenvectors, pos_omega, orb, lx=None, ly=None)
Compute U and V matrices, which are the projected position operators.
This function computes the U and V matrices (projected position operators) based on the given lattice coordinates and eigenvectors.
- Parameters:
lattice (ndarray) – Array of shape
(N_sites, 2)
containing the coordinates of the lattice sites.eigenvectors (ndarray) – Array of shape
(orb * N_sites, orb * N_sites)
containing the eigenvectors.pos_omega (int) – Position of the frequency in the ordered list of frequencies.
orb (int) – Number of orbitals.
lx (float, optional) – Size of the sample along the x-axis. If
None
, the function will determine it automatically.ly (float, optional) – Size of the sample along the y-axis. If
None
, the function will determine it automatically.
- Returns:
u_proj
andv_proj
, which are the projected position operators on the x and y coordinates, respectively.- Return type:
tuple of (ndarray, ndarray)
- pybott.bott.phase_diagram(lattice, ham_function, p1, p2, fermi_energy=0, name_of_file='phase_diagram.csv')
Generate a phase diagram by calculating the Bott index for each pair of parameters in p1 and p2.
This function calculates the Bott index for each pair of parameter values from p1 and p2 and generates a phase diagram. The results are saved in a CSV file, with columns for p1, p2, and the corresponding Bott index.
- Parameters:
lattice (ndarray) – Array of shape
(N_sites, 2)
containing the coordinates of the lattice sites.ham_function (callable) – Callable function that generates the Hamiltonian matrix given the parameters. It should have the signature
ham_function(param1, param2)
and return the Hamiltonian matrix.p1 (list) – List of values for the first parameter to vary in the phase diagram.
p2 (list) – List of values for the second parameter to vary in the phase diagram.
fermi_energy (float, optional) – The Fermi energy at which to calculate the Bott index. Default is
0
.name_of_file (str, optional) – Name of the output CSV file where the phase diagram will be saved. Default is
"phase_diagram.csv"
.
- Returns:
None
- Return type:
None
- Output:
A CSV file containing the phase diagram with columns for p1, p2, and the corresponding Bott index.
Note
This function iterates over each combination of values from p1 and p2, calculates the Bott index for each pair, and writes the results to the specified CSV file. The resulting file can be used to visualize the topological phases as a function of the two parameters.
- pybott.bott.phase_diagram_disorder(ham_lattice_function, disorder, energies, name_of_file='phase_diagram_disorder.csv', n_realisations=1)
Generate a phase diagram by calculating the averaged Bott index over multiple disorder realizations for a range of energy levels.
This function computes the Bott index for a series of energy levels and disorder strengths, averaging the results over multiple disorder realizations. The phase diagram is saved to a CSV file.
- Parameters:
ham_lattice_function (callable) – Callable function that generates the lattice and Hamiltonian matrix given a disorder parameter. It should have the signature
ham_lattice_function(disorder_value)
and return a tuple(lattice, hamiltonian)
.disorder (list) – A list of disorder strength values to use in generating the Hamiltonian.
energies (list) – A list of energy levels at which the Bott index will be calculated.
name_of_file (str, optional) – The name of the output CSV file where the phase diagram will be saved. Default is
"phase_diagram_disorder.csv"
.n_realisations (int, optional) – Number of disorder realizations to compute for each pair of (disorder, energy) values. The average Bott index over all realizations will be saved in the CSV file. Default is 1.
- Returns:
None
- Return type:
None
- Output:
Writes a CSV file with columns “energy”, “disorder”, and the averaged Bott index over all realizations for each combination of disorder and energy.
Note
The function outputs a progress bar showing the progress of the calculations across disorder values.
- pybott.bott.plot_phase_diagram(filename='phase_diagram.csv', title_fig='Phase Diagram', save_fig='phase_diagram.pdf', xkey='p2', ykey='p1', xlabel='p2', ylabel='p1', colorbar_label='Bott Index', fontsize=20, cmap='coolwarm')
Plot a phase diagram from a CSV file generated by the phase_diagram function.
This function reads a CSV file containing the phase diagram data (with columns ‘p1’, ‘p2’, and ‘Bott Index’) and generates a heatmap plot. The phase diagram can be saved as a figure in PDF format.
- Parameters:
filename (str, optional) – The name of the CSV file to read, which contains columns ‘p1’, ‘p2’, and ‘Bott Index’. Default is
"phase_diagram.csv"
.title_fig (str, optional) – The title of the plot. Default is
"Phase Diagram"
.save_fig (str, optional) – The name of the file to save the plot. Default is
"phase_diagram.pdf"
.xkey (str, optional) – The key for the x-axis data in the CSV file. Default is
"p2"
.ykey (str, optional) – The key for the y-axis data in the CSV file. Default is
"p1"
.xlabel (str, optional) – Label for the x-axis. Default is
"p2"
.ylabel (str, optional) – Label for the y-axis. Default is
"p1"
.colorbar_label (str, optional) – Label for the colorbar. Default is
"Bott Index"
.fontsize (int, optional) – Size of all the fonts in the plot. Default is 20.
cmap (str, optional) – The colormap to use for the heatmap. Default is
"coolwarm"
.
- Returns:
None
- Return type:
None
- Output:
Displays a heatmap plot of the phase diagram, with an option to save the plot as a PDF file.
Note
The function assumes that the input CSV file contains the following columns: ‘p1’, ‘p2’, and ‘Bott Index’.
- pybott.bott.sort_evals_evects(evals, evects, rev=False)
This function sorts the eigenvalues and eigenvectors in ascending or descending order depending on the rev flag. Useful when the energy are not exactly the eigenvalues, typically in photonics systems or non Hermitian systems in general.
- Parameters:
evals (ndarray) – Array containing the eigenvalues to be sorted.
evects (ndarray) – Array containing the corresponding eigenvectors.
rev (bool, optional) – If
True
, the eigenvalues and eigenvectors are sorted in descending order; otherwise, they are sorted in ascending order. Default isFalse
.
- Returns:
A tuple containing the sorted eigenvalues and the corresponding sorted eigenvectors.
- Return type:
tuple of (ndarray, ndarray)
pybott.spin_bott module
Code to compute the spin Bott index following the definition given by Huaqing Huang and Feng Liu in https://journals.aps.org/prb/abstract/10.1103/PhysRevB.98.125130
The Spin Bott index is used to characterize quantum spin Hall (QSH) states in both periodic and non-periodic systems.
The Spin Bott index is an extension of the Bott index that incorporates spin, allowing it to identify quantum spin Hall states. It involves projecting the spin operator onto the occupied states and calculating the Bott index separately for the spin-up and spin-down sectors. The Spin Bott index is defined as half the difference between the Bott indices of these two spin sectors.
Key steps:
Construct the projector onto the occupied states.
Calculate the projected position operators.
Compute the Bott index from the commutativity of the position operators.
For the Spin Bott index, introduce the projected spin operator and compute Bott indices for spin-up and spin-down components.
Functions included:
make_projector: Creates the projector matrix.
get_p_sigma_p_bott: Computes the PσP operator for the spin Bott index.
plot_psp_spectrum: Visualizes the spectrum of the PσP operator.
spin_bott: Calculates the spin Bott index for a given system configuration.
all_spin_bott: Computes the spin Bott index for all eigenvalues of PσP.
This method is applicable to both periodic and non-periodic systems, including disordered systems.
- pybott.spin_bott.all_spin_bott(grid, evals, evects, sigma, threshold_psp=0)
Calculate the spin Bott index for all eigenvalues of the PσP operator.
- Parameters:
grid (array) – The grid defining the system.
evals (array) – Array of eigenvalues of the Hamiltonian.
evects (array) – Array of eigenvectors of the Hamiltonian.
sigma (numpy.ndarray) – The spin operator matrix (σ).
threshold_psp (float, optional) – Threshold for the PσP projection. Defaults to 0.
- Returns:
A dictionary where each eigenvalue of PσP maps to its corresponding spin Bott index.
- Return type:
dict
- pybott.spin_bott.get_p_sigma_p_bott(w, vl, vr, sigma, omega)
Calculate the projected spin operator PσP using the eigenvectors and eigenvalues.
- Parameters:
w (array) – Array of eigenvalues.
vl (array) – Array of left eigenvectors.
vr (array) – Array of right eigenvectors.
sigma (numpy.ndarray) – The spin operator matrix (σ).
omega (float) – Threshold energy for the projection.
- Returns:
The PσP projected spin operator.
- Return type:
numpy.ndarray
- pybott.spin_bott.make_projector(vl, vr)
Create a projector from the left and right eigenvectors.
- Parameters:
vl (array) – Array of left eigenvectors.
vr (array) – Array of right eigenvectors.
- Returns:
The projector matrix constructed from the eigenvectors.
- Return type:
numpy.ndarray
- pybott.spin_bott.plot_psp_spectrum(threshold_energy, w_psp, name_psp)
Plot the spectrum of the PσP operator.
- Parameters:
threshold_energy (float) – Threshold value for energy levels to plot.
w_psp (array) – Array of eigenvalues of the PσP operator.
name_psp (str) – Filename to save the plot as a PDF.
- Returns:
None
- pybott.spin_bott.spin_bott(grid, evals, evects, sigma, fermi_energy=0, threshold_bott=-0.1, plot_psp=False, name_psp='spectrum_psp')
Calculate the spin Bott index for a given set of eigenvalues and eigenvectors.
- Parameters:
grid (array) – The grid defining the system.
evals (array) – Array of eigenvalues of the Hamiltonian.
evects (array) – Array of eigenvectors of the Hamiltonian.
sigma (numpy.ndarray) – The spin operator matrix (σ).
threshold_psp (float, optional) – Threshold for the PσP projection. Defaults to 0.
threshold_energy (float, optional) – Threshold for energy levels. Defaults to 0.
plot_psp (bool, optional) – Whether to plot the PσP spectrum. Defaults to False.
name_psp (str, optional) – Filename to save the PσP plot as a PDF.
- Returns:
The spin Bott index of the system.
- Return type:
float