Electron
Tight-binding models from Wannier90.
- class elphmod.el.Model(seedname=None, N=None, a=None, r=None, divide_ndegen=True, read_xsf=False, normalize_wf=False, buffer_wf=False, check_ortho=False, rydberg=False, shared_memory=False)[source]
Tight-binding model for the electrons.
- Parameters:
- seednamestr
Common prefix of Wannier90 output files: seedname_hr.dat with the Hamiltonian in the Wannier basis and, optionally, seedname_wsvec.dat with superlattice vectors used to symmetrize the long-range hopping. Alternatively, dat.h_mat_r from RESPACK can be used.
- Ntuple of int, optional
Numbers of unit cells per direction on which RESPACK data is defined. This can be omitted if all numbers are even.
- andarray, optional
Bravais lattice vectors used to map RESPACK data to Wigner-Seitz cell. By default, a cubic cell is assumed.
- rndarray, optional
Positions of orbital centers used to map RESPACK data to Wigner-Seitz cell. By default, all orbitals are assumed to be located at the origin of the unit cell.
- divide_ndegenbool, default True
Divide hopping by degeneracy of Wigner-Seitz point and apply the abovementioned correction? Only
True
yields correct bands.False
is only used in combination withdecayH()
.- read_xsfbool, default False
Read Wannier functions in position representation in the XCrySDen structure format (XSF)?
- normalize_wfbool, default False
Normalize Wannier functions? This is only recommended if a single and complete image of each Wannier function is contained in the supercell written to the XSF file. This may not be the case for a small number of k points along one direction (e.g., a single k point in the out-of-plane direction in the case of 2D materials), especially in combination with a large supercell.
- buffer_wfbool, default False
After reading XSF files, store Wannier functions and corresponding real-space mesh in binary files, which will be read next time? Since reading the XSF files is slow, this can save a lot of time for large systems. Make sure to delete the binary files seedname_wf.npy and seedname_xyz.npy whensoever the XSF files change.
- check_orthobool, default False
Check if Wannier functions are orthogonal?
- rydbergbool, default False
Convert energies from eV to Ry?
- shared_memorybool, default False
Store Wannier functions in shared memory?
- Attributes:
- Rndarray
Lattice vectors \(\vec R\) of Wigner-Seitz supercell.
- datandarray
Corresponding on-site energies and hoppings in eV.
\[H_{\vec R \alpha \beta} = \bra{0 \alpha} H \ket{\vec R \beta}\]If
rydberg
isTrue
, the units are Ry instead.- sizeint
Number of Wannier functions/bands.
- nktuple of int
Guessed shape of original k-point mesh.
- cellslist of tuple of int, optional
Lattice vectors of unit cells if the model describes a supercell.
- Nlist of tuple of int, optional
Primitive vectors of supercell if the model describes a supercell.
- Wndarray, optional
Wannier functions in position representation if read_xsf.
- rndarray, optional
Cartesian positions belonging to Wannier functions if read_xsf.
- taundarray, optional
Positions of basis atoms if read_xsf.
- atom_orderlist of str, optional
Ordered list of atoms if read_xsf.
- dVfloat, optional
Volume element/voxel volume belonging to r if read_xsf.
- divide_ndegenbool
Have hoppings been divided by degeneracy of Wigner-Seitz point?
- rydbergbool
Have energies been converted from eV to Ry?
- H(k1=0, k2=0, k3=0)[source]
Set up Hamilton operator for arbitrary k point.
- Parameters:
- k1, k2, k3float, default 0.0
k point in crystal coordinates with period \(2 \pi\).
- Returns:
- ndarray
Fourier transform of
data
.
- order_orbitals(*order)[source]
Reorder Wannier functions.
Warning: Wannier functions in position representation and related attributes remain unchanged!
Together with
shift_orbitals()
, this function helps reconcile inconsistent definitions of the basis/motif of the Bravais lattice.- Parameters:
- *orderint
New order of Wannier functions.
- shift_orbitals(s, S)[source]
Move selected Wannier functions across unit-cell boundary.
Warning: Wannier functions in position representation and related attributes remain unchanged!
Together with
order_orbitals()
, this function helps reconcile inconsistent definitions of the basis/motif of the Bravais lattice.- Parameters:
- sslice
Slice of orbital indices corresponding to selected basis atom(s).
- Stuple of int
Shift of as multiple of primitive lattice vectors.
- standardize(eps=0.0, symmetrize=False)[source]
Standardize tight-binding data.
Keep only nonzero hopping matrices.
Sum over repeated lattice vectors.
Sort lattice vectors.
Optionally symmetrize hopping:
\[H_{\vec k} = H_{\vec k}^\dagger, H_{\vec R} = H_{-\vec R}^\dagger\]- Parameters:
- epsfloat
Threshold for “nonzero” matrix elements in units of the maximum matrix element.
- symmetrizebool
Symmetrize hopping?
- supercell(N1=1, N2=1, N3=1, sparse=False)[source]
Map tight-binding model onto supercell.
- Parameters:
- N1, N2, N3tuple of int or int, default 1
Supercell lattice vectors in units of primitive lattice vectors.
- sparsebool, default False
Only calculate k = 0 Hamiltonian as a sparse matrix to save memory? The result, which is assumed to be real, is stored in the attribute
Hs
. Consider usingstandardize()
with nonzero eps and symmetrize before.
- Returns:
- object
Tight-binding model for supercell.
See also
- t(R1=0, R2=0, R3=0)[source]
Get on-site or hopping energy for arbitrary lattice vector.
- Parameters:
- R1, R2, R3int, default 0
Lattice vector in units of primitive vectors.
- Returns:
- ndarray
Element of
data
or zero.
- to_hrdat(seedname)[source]
Save tight-binding model to _hr.dat file.
- Parameters:
- seednamestr
Common prefix of Wannier90 input and output files.
- v(k1=0, k2=0, k3=0)[source]
Set up band-velocity operator for arbitrary k point.
- Parameters:
- k1, k2, k3float, default 0.0
k point in crystal coordinates with period \(2 \pi\).
- Returns:
- ndarray
k derivative of Fourier transform of
data
in the basis of primitive lattice vectors and orbitals. It can be transformed into the Cartesian and band basis afterward (Hellmann-Feynman theorem).
- elphmod.el.decayH(seedname, **kwargs)[source]
Calculate the decay of the Hamiltonian.
This function should yield the same data as
read_decayH()
.- Parameters:
- seednamestr
Prefix of Wannier90 output file.
- **kwargs
Arguments for
elphmod.bravais.primitives()
: Choose the right Bravais lattice (ibrav
) and lattice constants (a, b, c, ...
).For a simple cubic lattice:
decayH(seedname, ibrav=1, a=a)
.
- Returns:
- Rndarray
The distance of every Wigner-Seitz grid point measured from the center in angstrom.
- Hndarray
The maximum absolute value of the elements of the Hamiltonian matrix in rydberg.
- elphmod.el.demet_from_qe_pwo(pw_scf_out, subset=None)[source]
Calculate the entropy contribution \(-TS\) to the total free energy.
In Quantum ESPRESSO
demet
is calculated ingweights.f90
.- Parameters:
- pw_scf_outstr
The name of the output file (typically
'pw.out'
).- subsetlist or array
List of indices to pick only a subset of the bands for the integration.
- Returns:
- demetfloat
The \(-TS\) contribution to the total free energy.
- elphmod.el.eband(pw_scf_out, subset=None)[source]
Calculate
eband
part of one-electron energy.The ‘one-electron contribution’ energy in the Quantum ESPRESSO PWscf output is a sum
eband + deband
. Here, we can calculate theeband
part.To compare it with the Quantum ESPRESSO result, you need to modify the
SUBROUTINE print_energies ( printout )
from electrons.f90.Change:
WRITE( stdout, 9062 ) (eband + deband), ehart, ( etxc - etxcc ), ewld
to:
WRITE( stdout, 9062 ) eband, & (eband + deband), ehart, ( etxc - etxcc ), ewld
and:
9062 FORMAT( ' one-electron contribution =',F17.8,' Ry' &
to:
9062 FORMAT( ' sum bands =',F17.8,' Ry' & /' one-electron contribution =',F17.8,' Ry' &
- Parameters:
- pw_scf_outstr
Name of the output file (typically
'pw.out'
or'scf.out'
).
- Returns:
- ebandfloat
The band energy.
- elphmod.el.k2r(el, H, a, r, fft=True, rydberg=False)[source]
Interpolate Hamilontian matrices on uniform k-point mesh.
- Parameters:
- el
Model
Tight-binding model.
- Hndarray
Hamiltonian matrices on complete uniform k-point mesh.
- andarray
Bravais lattice vectors.
- rndarray
Positions of orbital centers.
- fftbool
Perform Fourier transform? If
False
, only the mapping to the Wigner-Seitz cell is performed.- rydbergbool, default False
Is input Hamiltonian given in Ry rather than eV units? This is independent of
el.rydberg
, which is always respected.
- el
- elphmod.el.proj_sum(proj, orbitals, *groups, **kwargs)[source]
Sum over selected atomic projections.
Example:
proj = read_atomic_projections('atomic_proj.xml') orbitals = read_projwf_out('projwfc.out') proj = proj_sum(proj, orbitals, 'S-p', 'Ta-d{z2, x2-y2, xy}')
- Parameters:
- projndarray
Atomic projections from
read_atomic_projections()
.- orbitalslist of str
Names of all orbitals from
read_projwfc_out()
.- *groups
Comma-separated lists of names of selected orbitals. Omitted name components are summed over. Curly braces are expanded.
- otherbool, default False
Return remaining orbital weight too? If you just want to select the “other” orbitals from
read_atomic_projections()
, add the namex
to orbitals instead.
- Returns:
- ndarray
Summed-over atomic projections.
- elphmod.el.read_Fermi_level(pw_scf_out)[source]
Read Fermi level from output of self-consistent PW run.
- Parameters:
- pw_scf_outstr
PWscf output file.
- Returns:
- float
Fermi level if found,
None
otherwise.
- elphmod.el.read_atomic_projections(atomic_proj_xml, order=False, from_fermi=True, squared=True, other=False, **order_kwargs)[source]
Read projected bands from outdir/prefix.save/atomic_proj.xml.
- Parameters:
- atomic_proj_xmlstr
XML file with atomic projections generated by
projwfc.x
.- orderbool
Order/disentangle bands via their k-local character?
- from_fermibool
Subtract Fermi level from electronic energies?
- squaredbool
Return squared complex modulus of projection?
- otherbool
Estimate projection onto “other” orbitals not defined in pseudopotential files as difference of band weights to one? This requires squared.
- **order_kwargs
Keyword arguments passed to
band_order()
.
- Returns:
- ndarray
k points.
- ndarray
Cumulative path distance.
- ndarray
Electronic energies.
- ndarray
Projections onto (pseudo) atomic orbitals.
- elphmod.el.read_atomic_projections_old(atomic_proj_xml, order=False, from_fermi=True, **order_kwargs)[source]
Read projected bands from outdir/prefix.save/atomic_proj.xml.
- elphmod.el.read_bands(filband)[source]
Read bands from filband just like Quantum ESRESSO’s
plotband.x
.- Parameters:
- filbandstr
Filename.
- Returns:
- kndarray
k points in Cartesian coordiantes.
- bandsndarray
Band energies.
- elphmod.el.read_bands_plot(filbandgnu, bands)[source]
Read bands from filband.gnu produced by Quantum ESPRESSO’s
bands.x
.- Parameters:
- filbandgnustr
Name of file with plotted bands.
- bandsint
Number of bands.
- Returns:
- ndarray
Cumulative reciprocal distance.
- ndarray
Band energies.
- elphmod.el.read_decayH(file)[source]
Read decay.H output from EPW.
- Parameters:
- filestr
The name of the decay.H output from EPW.
- Returns:
- Rndarray
The distance of every Wigner-Seitz grid point measured from the center in angstrom.
- Hndarray
The maximum absolute value of the elements of the Hamiltonian matrix in rydberg.
- elphmod.el.read_energy_contributions_scf_out(filename)[source]
Read energy contributions to the total energy
scf
output file.- Parameters:
- filenamestr
Name of Quantum ESPRESSO’s
scf
output file.
- Returns:
- dict
Energy contributions.
- elphmod.el.read_eps_nk_from_qe_pwo(pw_scf_out)[source]
Read electronic eigenenergies and k points and calculate the occupations.
- Parameters:
- pw_scf_outstr
Name of the output file (typically
'pw.out'
or'scf.out'
).
- Returns:
- energiesndarray
Electronic eigenenergies (
eps_nk
) fromscf
output, shape(nk, nbnd)
.- kpointsndarray
k points and weights from
scf
output, shape(nk, 4)
, where the 4th column contains the weights.- f_occndarray
Occupations of
eps_nk
(same shape as energies).- smearing_typestr
Type of smearing (for example
'Fermi-Dirac'
).- mufloat
Chemical potential in eV.
- kTfloat
Value of smearing (
degauss
) in eV.
- elphmod.el.read_hrdat(seedname, divide_ndegen=True)[source]
Read _hr.dat (or _tb.dat) file from Wannier90.
- Parameters:
- seednamestr
Common prefix of Wannier90 input and output files.
- divide_ndegenbool
Divide hopping by degeneracy of Wigner-Seitz point?
- Returns:
- ndarray
Lattice vectors.
- ndarray
On-site and hopping parameters.
- elphmod.el.read_pp_density(filename)[source]
Read file filout with charge density generated by
pp.x
.Calculate total charge
\[Q = \int \rho(\vec r) dx dy dz = = \frac \Omega {n_1 n_2 n_3} \sum_{i j k} \rho_{i j k},\]where \(\Omega\) is the unit-cell volume, \(n_r\) are the FFT grid dimensions, and \(i, j, k\) run over FFT real-space grid points.
- Parameters:
- filenamestr
Filename.
- Returns:
- rhondarray
Electronic charge density \(\rho_{i j k}\) on FFT real-space grid points.
- tot_charge: float
Total charge calculated from charge density \(\rho_{i j k}\).
- elphmod.el.read_projwfc_out(projwfc_out, other=True)[source]
Identify orbitals in atomic_proj.xml via output of
projwfc.x
.- Parameters:
- projwfc_outstr
Output file of
projwfc.x
.
- Returns:
- list of str
Common names of (pseudo) atomic orbitals listed in projwfc_out (in that order). If spin-orbit coupling is considered, symmetry labels related to the magnetic quantum number are omitted.
- otherbool, default True
Add name
X-x
to list of orbitals, which corresponds to “other” orbitals fromread_atomic_projections()
?
- elphmod.el.read_rhoG_density(filename, ibrav, a=1.0, b=1.0, c=1.0)[source]
Read the charge density output from Quantum ESPRESSO.
The purpose of \(\rho(\vec G)\) is to calculate the charge density in real space \(\rho(\vec r)\) or the Hartree energy
ehart
.- Parameters:
- filenamestr
Filename.
- ibravinteger
Bravais lattice index (see
pw.x
input description).- a, b, cfloat
Bravais lattice parameters.
- Returns:
- rho_gndarray
Electronic charge density \(\rho(\vec G)\) on reciprocal-space grid points.
- g_vectndarray
Reciprocal lattice vectors \(\vec G\).
- ngm_ginteger
Number of reciprocal lattice vectors.
- uc_volumefloat
Unit-cell volume.
- mill_gndarray
Miller indices.
- elphmod.el.read_symmetry_points(bandsout)[source]
Read positions of symmetry points along path.
- Parameters:
- bandsoutstr
File with standard output from Quantum ESPRESSO’s
bands.x
.
- Returns:
- list
Positions of symmetry points.
- elphmod.el.read_wannier90_eig_file(seedname, num_bands, nkpts)[source]
Read Kohn-Sham energies (eV) from the Wannier90 output seedname.eig file.
- Parameters:
- seednamestr
For example ‘tas2’, if the file is named ‘tas2.eig’.
- num_bandsint
Number of bands in your pseudopotential.
- nkptsint
Number of k points in your Wannier90 calculations. For example 1296 for 36x36x1.
- Returns:
- ndarray
Kohn-Sham energies:
eig[num_bands, nkpts]
.
- elphmod.el.read_wfc(filename, ibrav, a=1.0, b=1.0, c=1.0)[source]
Read the wave function output from Quantum ESPRESSO.
\[\psi_{n \vec k}(\vec r) = \frac 1 {\sqrt V} \sum_{\vec G} c_{n, \vec k + \vec G} \E^{\I (\vec k + \vec G) \vec r}\]- Parameters:
- filenamestr
Filename.
- ibravinteger
Bravais lattice index (see
pw.x
input description).- a, b, cfloat
Bravais lattice parameters.
- Returns:
- evcndarray
Electronic wave function \(c_{n, \vec k + \vec G}\).
- igwxinteger
Number of reciprocal lattice vectors.
- xkndarray
k point.
- k_plus_Gndarray
k point plus reciprocal lattice vectors \(\vec G\).
- g_vectndarray
Reciprocal lattice vectors \(\vec G\).
- mill_gndarray
Miller indices.
- elphmod.el.write_bands(filband, k, bands, fmt='%8.3f', cols=10)[source]
Write bands to filband just like Quantum ESPRESSO’s
plotband.x
.- Parameters:
- filbandstr
Filename.
- kndarray
k points in Cartesian coordiantes.
- bandsndarray
Band energies.
- fmtstr
Format string for band energies.
- colsint
Number of band energies per line.
- elphmod.el.write_hrdat(seedname, R, H, ndegen=None)[source]
Write _hr.dat file as generated by Wannier90.
- Parameters:
- seednamestr
Common prefix of Wannier90 input and output files.
- Rndarray
Lattice vectors of Wigner-Seitz supercell.
- Hndarray
Corresponding on-site energies and hoppings.
- ndegenlist of int
Degeneracies of Wigner-Seitz lattice vectors. This is just what is written to the file header, and it is not used to further modify H.