kpoint
- class irrep.kpoint.Kpoint(ik=None, num_bands=None, RecLattice=None, symmetries_SG=None, calculate_traces=False, spinor=None, kpt=None, WF=None, Energy=None, ig=None, upper=None, degen_thresh=1e-08, symmetries=None, refUC=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), shiftUC=array([0., 0., 0.]), symmetries_tables=None, save_wf=True, v=0)[source]
Parses files and organizes info about the states and energy-levels of a particular k-point in attributes. Contains methods to calculate and write traces (and irreps), for the separation of the band structure in terms of a symmetry operation and for the calculation of the Zak phase.
symmetries=None, symmetries_tables=None # calculate_traces needs it
- Parameters:
ik (int) – Index of kpoint, starting count from 0.
NBin (int) – Number of bands considered at the k-point in the DFT calculation.
Ecut (float) – Plane-wave cutoff (in eV) to consider in the expansion of wave-functions.
Ecut0 (float) – Plane-wave cutoff (in eV) used in the DFT calulation. Always read from DFT files. Insignificant if `code`=’wannier90’.
RecLattice (array, shape=(3,3)) – Each row contains the cartesian coordinates of a basis vector forming the unit-cell in reciprocal space.
symmetries_SG (list, default=None) – Each element is an instance of class SymmetryOperation corresponding to a symmetry in the point group of the space-group. The value passed is the attribute symmetries of class SpaceGroup.
spinor (bool, default=None) – True if wave functions are spinors, False if they are scalars.
kpt (list or array, default=None) – Direct coordinates of the k-point.
Energy (array) – Energy levels of the states with band indices between IBstart and IBend.
ig (array) – Array returned by
sortIG()
. It contains data about the plane waves in the expansion of wave functions.upper (float) – Energy of the state `IBend`+1. Used to calculate the gap with upper bands.
degen_thresh (float, default=1e-8) – Threshold to identify degenerate energy levels.
refUC (array, default=None) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=None) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
symmetries_tables (list) – Attribute symmetries of class IrrepTable. Each component is an instance of class SymopTable corresponding to a symmetry operation in the “point-group” of the space-group.
calculate_traces (bool) – If True, traces of symmetries will be calculated. Useful to create instances faster.
save_wf (bool) – Whether wave functions should be kept as attribute after calculating traces.
v (int) – Verbosity level. Default set to minimalistic printing
- spinor
True if wave-functions are spinors, False if they are scalars.
- Type:
bool
- ik0
Index of the k-point, starting the count from 1.
- Type:
int
- num_bands
Number of bands whose traces should be calculated.
- Type:
int
- RecLattice
Each row contains the cartesian coordinates of a basis vector forming the unit-cell in reciprocal space.
- Type:
array, shape=(3,3)
- WF
Coefficients of wave-functions in the plane-wave expansion. A row for each wave-function, a column for each plane-wave.
- Type:
array
- ig
Returned by sortIG. Every column corresponds to a plane-wave of energy smaller than Ecut. The number of rows is 6: the first 3 contain direct coordinates of the plane-wave, the fourth row stores indices needed to short plane-waves based on energy (ascending order). Fitfth (sixth) row contains the index of the first (last) plane-wave with the same energy as the plane-wave of the current column.
- Type:
array
- k
Direct coordinates of the k point in the DFT cell setting.
- Type:
array, shape=(3,)
- K
Property getter for self.k. NEEDED TO PRESERVE COMPATIBILITY WITH BANDUPPY<=0.3.3. DO NOT CHANGE UNLESS NECESSARY. NOTIFY THE DEVELOPERS IF ANY CHANGES ARE MADE.
- Type:
array, shape=(3,)
- k_refUC
Direct coordinates of the k point in the reference cell setting.
- Type:
array, shape=(3,)
- Energy_raw
energy levels of each state (1 energy for every row in WF)
- Type:
array
- Energy_mean
Energy-levels of degenerate froups of bands whose traces should be calculated.
- Type:
array
- upper
Energy of the first band above the set of bands whose traces should be calculated. It will be set to numpy.NaN if the last band matches the last band in the DFT calculation (default IBend).
- Type:
float
- symmetries_SG
Each element is an instance of class SymmetryOperation corresponding to a symmetry in the point group of the space-group. The value passed is the attribute symmetries of class SpaceGroup.
- Type:
list
- symmetries
Each key is an instance of class SymmetryOperation corresponding to an operation in the little-(co)group and the attached value is an array with the traces of the operation.
- Type:
dict
- char
Each row corresponds to a set of degenerate states. Each column is the trace of a symmetry in the little cogroup in the DFT cell setting.
- Type:
array
- char_refUC
The same as char, but in the reference cell setting.
- Type:
array
- degeneracies
Degeneracies of energy levels between IBstart and IBend.
- Type:
array
- borders
Integers representing the band index of the first state in each set of degenerate states. The bounds can be obtained as for ibot, itop in zip(borders[:-1], borders[1:]).
- Type:
array
- Energy_mean
Average of energy levels within each set of degenerate states
- Type:
array
- num_bandinvs
Property getter for the number of inversion-odd states. If the k point is not inversion symmetric, None is returned.
- Type:
int
- NG
Property getter for the number of plane waves in the expansion of wave functions.
- Type:
int
- onlytraces
False if irreps have been identified and have to be written.
- Type:
bool
- property K
Getter for the redfuced coordinates of the k-point needed to keep compatibility with banduppy
ACCESSED BY BANDUPPY.NEEDED TO PRESERVE COMPATIBILITY WITH BANDUPPY<=0.3.3. AVOID CHANGING UNLESS NECESSARY.
- property NG
Getter for the number of plane-waves in current k-point
- Separate(symop, degen_thresh, groupKramers=True, v=0)[source]
Separate the band structure in a particular k-point according to the eigenvalues of a symmetry operation.
- Parameters:
isymop (int) – Index of symmetry used for the separation.
degen_thresh (float, default=1e-5) – Energy threshold used to determine degeneracy of energy-levels.
groupKramers (bool, default=True) – If True, states will be coupled by pairs of Kramers.
v (int, default=0) – Verbosity level. Default set to minimalistic printing
- Returns:
subspaces – Each key is an eigenvalue of the symmetry operation and the corresponding value is an instance of class Kpoint for the states with that eigenvalue.
- Return type:
dict
- calculate_traces(refUC, shiftUC, symmetries_tables, v=0)[source]
Calculate traces of symmetry operations
- Parameters:
refUC (array, default=None) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=None) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
symmetries (list) – Indices of symmetries whose traces will be calculated.
symmetries_tables (list) – Attribute symmetries of class IrrepTable. Each component is an instance of class SymopTable corresponding to a symmetry operation in the “point-group” of the space-group.
degen_thresh (float, default=1e-8) – Threshold to identify degenerate energy levels.
- Returns:
char (array) – Each row corresponds to a set of degenerate states. Each column is the trace of a symmetry in the little cogroup in the DFT cell setting.
char_refUC (array) – The same as char, but in the reference cell setting.
Energy_mean (array) – Average of energy levels within each set of degenerate states
- copy_sub(E, WF, inds)[source]
Create an instance of class Kpoint for a restricted set of states.
- Parameters:
E (array) – Energy-levels of states.
WF (array) – Coefficients of the plane-wave expansion of wave-functions. Each row corresponds to a wave-function, each row to a plane-wave.
- Returns:
other – Instance of Kpoints corresponding to the group of states passed. They are shorted by energy-levels.
- Return type:
class
- get_rho_spin(degen_thresh=0.0001)[source]
Evaluates the matrix <i|M|j> in every group of degenerate bands labeled by i and j, where M is \(\sigma_0\), \(\sigma_x\), \(\sigma_y\) or \(\sigma_z\) for the spinor case, \(\sigma_0\) for the spinor case.
- Parameters:
degen_thresh (float) – Bands with energy difference smaller that the threshold will be considered as one group.
- Returns:
Each tuple contains b1 ,b2, E, (M , Sx , Sy , Sz), where b1 - index of first band in the group b2 - index of the last band inthe group +1 E - average energy of the group M - <i|j> Sx, Sy, Sz - <i|sigma|j>
- Return type:
list of tuples
- identify_irreps(irreptable=None)[source]
Identify irreps based on traces. Sets attributes onlytraces and irreps.
- Parameters:
irreptable (dict) – Each key is the label of an irrep, each value another dict. Keys of every secondary dict are indices of symmetries (starting from 1 and following order of operations in tables of BCS) and values are traces of symmetries.
- json()[source]
Prepare the data to save it in JSON format.
- Returns:
json_data – Data that will be saved in JSON format.
- Return type:
dict
- k_close_mod1(kpt, prec=1e-06)[source]
Check if the k-point is close to another k-point modulo 1. (in reduced coordinates)
ACCESSED BY BANDUPPY. AVOID CHANGING UNLESS NECESSARY. NOTIFY DEVELOPERS IF ANY CHANGE IS MADE.
- Parameters:
kpt (array) – Coordinates of the k-point to compare.
prec (float, default=1e-6) – Threshold to consider the k-points as equal.
- overlap(other)[source]
Calculates the overlap matrix of elements < u_m(k) | u_n(k+g) >.
- Parameters:
other (class) – Instance of Kpoints corresponding to k+g (next k-point in path).
- Returns:
res – Matrix of complex elements < u_m(k) | u_n(k+g) >.
- Return type:
array
- unfold(supercell, kptPBZ, degen_thresh=0.0001)[source]
Unfolds a kpoint of a supercell onto the point of the primitive cell kptPBZ.
- Parameters:
supercell (array, shape=(3,3)) – Describes how the lattice vectors of the (super)cell used in the calculation are expressed in the basis vectors of the primitive cell.
kptPBZ (array, shape=(3,)) – Coordinates of the k-point in the primitive Brillouin zone (PBZ), on which the present kpoint should be unfolded.
degen_thresh (float) – Bands with energy difference smaller that the threshold will be considered as one band, and only one total weight will be given for them.
- Returns:
Array containing 2 columns (5 for spinor case): E, W, Sx, Sy, Sz. E - energy of the band or average energy of the group of bands. W - weight of the band(s) projected onto the PBZ kpoint. Sx, Sy, Sz - Spin components projected onto the PBZ kpoint.
- Return type:
array
- write_characters()[source]
Write the block of data of the k point, including energy levels, degeneracies, traces and irreps.
- write_irrepsfile(file)[source]
Write the irreps of this k point into irreps.dat file.
- Parameters:
file (File object) – File object for the irreps.dat file.
- write_trace()[source]
Write in trace.txt the block corresponding to a single k-point.
- Parameters:
degen_thresh (float, default=1e-8) – Threshold energy used to decide whether wave-functions are degenerate in energy.
symmetries (list, default=None) – Index of symmetry operations whose traces will be printed.
efermi (float, default=0.0) – Fermi-energy. Used as origin for energy-levels.
- Returns:
res – Block to write in trace.txt with description of traces in a single k-point.
- Return type:
str