Input/Output
The reading and writing functions are implemented in the WannierIO.jl package. However, here are also some convenience functions which wrap the corresponding functions in WannierIO.jl, to utilize the structs defined in Wannier90.jl, e.g. BVectors, RGrid, etc.
In most cases, the units of the function arguments and returns are in angstrom unit for lattice, and fractional w.r.t lattice for atomic positions, etc.
The following abbreviations are used throughout the code and documentation:
n_bandsfor number of bandsn_wannfor number of WFsn_kptsfor number of kpointsn_bvecsfor number of b-vectorsn_atomsfor number of atomsUforamnor the gauge matricesMformmnmatricesEforeigmatricesSforspnmatrices
In most cases, for arrays we adopt the convention that n_bands is the first index, n_wann is the second index, and n_kpts is the third index. For example, U for the gauge matrices is a 3D array of size (n_bands, n_wann, n_kpts).
Contents
Index
Brillouin.KPaths.KPathInterpolantWannier.ModelWannier._read_w90_tbWannier.compute_mudWannier.get_symm_idx_labelWannier.read_bxsfWannier.read_cubeWannier.read_nnkpWannier.read_orthonorm_amnWannier.read_w90Wannier.read_w90_bandWannier.read_w90_interpWannier.read_w90_tbWannier.read_xsfWannier.truncateWannier.truncate_mmn_eigWannier.truncate_unkWannier.truncate_w90Wannier.write_bxsfWannier.write_chkWannier.write_chkWannier.write_cubeWannier.write_nnkpWannier.write_w90Wannier.write_w90_bandWannier.write_w90_kpt_labelWannier.write_xsf
Wannier90 files
Wannier.read_orthonorm_amn — Methodread_orthonorm_amn(filename::AbstractString)Read and orthonormalize the amn file.
Wrapper function to read amn and Lowdin orthonormalize it. The U matrix needs to be unitary or semi-unitary, so in most cases this function should be used instead of WannierIO.read_amn.
See also WannierIO.read_amn.
Brillouin.KPaths.KPathInterpolant — MethodKPathInterpolant(kpoints, symm_idx, symm_label, recip_lattice)Generate a KPathInterpolant from kpoints in seedname_band.dat/kpt/labelinfo.
Arguments
- kpoints: fractional coordinate, each column is a kpoint.
Wannier.get_symm_idx_label — Methodget_symm_idx_label(kpi::KPathInterpolant)Return the symmetry indexes and labels.
Wannier.read_w90_band — Methodread_w90_band(seedname::AbstractString, recip_lattice::AbstractMatrix)Arguments
recip_lattice: each column is a reciprocal lattice vector in Cartesian coordinates. If given, return a tuple of(KPathInterpolant, E). This is a more user-friendly version ofread_w90_band(seedname::AbstractString).
See also read_w90_band(seedname::AbstractString).
Wannier.write_w90_band — Methodwrite_w90_band(seedname::AbstractString, kpi::KPathInterpolant, E::AbstractMatrix)Write SEEDNAME_band.dat, SEEDNAME_band.kpt, SEEDNAME_band.labelinfo.dat.
This is a more user-friendly version.
See also write_w90_band(seedname, kpoints, E, x, symm_idx, symm_label).
Wannier.write_w90_kpt_label — Methodwrite_w90_kpt_label(seedname::AbstractString, kpi::KPathInterpolant)Write SEEDNAME_band.kpt and SEEDNAME_band.labelinfo.dat.
This allows generating the high-symmetry kpoints and labels from crystal structure, and use the generated kpoints in pw.x bands calculation or in the win input file for Wannier90.
Example
win = read_win("si2.win")
kp = get_kpath(win.unit_cell_cart, win.atoms_frac, win.atom_labels)
kpi = Wannier.interpolate_w90(kp, 100)
Wannier.write_w90_kpt_label("si2", kpi)Wannier.Model — MethodModel(chk::Chk)Construct a model from a WannierIO.Chk struct.
Arguments
chk: aWannierIO.Chkstruct
Keyword Arguments
E: an_wann * n_kptsarray for the eigenvalues of the Hamiltonian. If not provided, it will be set to zero.U: an_wann * n_wann * n_kptsarray for the gauge transformation.
The Chk struct does not contain eigenvalues, thus if E is not provided, it will be set to zero.
Moreover, the M matrix in Chk is already rotated by the gauge transformation, thus by default, the U matrix is set to identity. Note that although maximal localization, or disentanglement (after frozen states are chosen), do not require eigenvalues (so the user can still Wannierize the Model), it is required when writing the Model to a .chk file, in write_chk.
Additionally, be careful that the M matrix is rotated, and this rotation needs to make sure that the rotated Hamiltonian is diagonal so that E stores the diagonal eigenvalues of the Hamiltonian.
Wannier.write_chk — Methodwrite_chk(filename, model, U; exclude_bands=nothing, binary=false)Write the model to a Wannier90 .chk file, using the gauge U.
Arguments
filename: filename of the.chkfilemodel: aModelstructU:n_bands * n_wann * n_kptsarray for the gauge transformation
Keyword Arguments
exclude_bands: a list of band indices to exclude. This is irrelevant to themodel, butchkfile has this entry.binary: whether to write the.chkfile in binary format.
Wannier.write_chk — Methodwrite_chk(filename, model; exclude_bands=nothing, binary=false)Write the model to a Wannier90 .chk file.
Arguments
filename: filename of the.chkfilemodel: aModelstruct
Keyword Arguments
exclude_bands: a list of band indices to exclude. This is irrelevant to themodel, butchkfile has this entry.binary: whether to write the.chkfile in binary format.
Wannier.read_w90 — Methodread_w90(seedname::AbstractString; amn=true, orthonorm_amn=true, mmn=true, eig=true)Read win, and optionally amn, mmn, eig.
Keyword arguments
- amn: if
true, readamnfile - orthonorm_amn: Lowdin orthonormalization after reading
amnmatrices. Should betruefor most cases, since usually the inputamnmatrices are projections onto atomic orbitals, and are not unitary or semi-unitary. - mmn: if
true, readmmnfile - eig: if
true, readeigfile
Wannier.read_w90_interp — Methodread_w90_interp(seedname::AbstractString; chk=true, amn=nothing, mdrs=nothing)Return an InterpModel for Wannier interpolation.
Keyword arguments
- chk: if
true, readchkorchk.fmtfile to get the unitary matrices, otherwise readamnfile for unitary matrices. - amn: filename for
amn, read thisamnfile for unitary matrices. Only used if not readingchkandamnis given. - mdrs: use MDRS interpolation, else Wigner-Seitz interpolation. If is
nothing, detect fromwinfile; and if nouse_ws_distanceinwinfile, default totrue.
Wannier.write_w90 — Methodwrite_w90(seedname::AbstractString, model::Model)Write Model into eig, mmn, amn files.
Keyword arguments
binary: writeeig,mmn, andamnin Fortran binary format
Wannier.read_nnkp — Methodread_nnkp(filename::AbstractString)Read the nnkp file.
This is a wrapper of WannierIO.read_nnkp. It returns a BVectors instead of NamedTuple.
Wannier.write_nnkp — Functionwrite_nnkp(filename::AbstractString, bvectors::BVectors, n_wann::Integer)Write nnkp that can be used by pw2wannier90.
This is a wrapper of WannierIO.write_nnkp.
Wannier._read_w90_tb — Method_read_w90_tb(seedname::AbstractString)Read seedname_tb.dat and seedname_wsvec.dat.
Return
- R vectors
- Hamiltonian
- position operator
Wannier.read_w90_tb — Methodread_w90_tb(seedname; kpoints, atom_positions, atom_labels)Read seedname_tb.dat and seedname_wsvec.dat, return an InterpModel.
Arguments
seedname: the seedname ofseedname_tb.datandseedname_wsvec.dat
Keyword Arguments
kpoints: the kpoints used in Wannierization, each column is a fractional coordinate.atom_positions: columns are the fractional coordinates of atomsatom_labels: labels of atoms
Return
- An
InterpModel.
Usually, if you only need to interpolate operators (e.g., the band structure), that is only inverse Fourier transform invfourier is needed, then you don't need to provide kpoints. If in some cases, you want to generate Wannier gauge operators, then passing the kpoints allows you to construct a complete KRVectors, and you can subsequently call fourier on Bloch gauge operators to get Wannier gauge operators.
Since no atomic positions and labels are provided, the auto-generated Brillouin.KPath is probably wrong. It is recommended that you pass the atom_positions and atom_labels arguments so that this function can auto generate the correct Brillouin.KPath.
File manipulation
Truncate Wannier90 matrices
Here are some functions to remove some bands from mmn, eig, or UNK files, so as to skip rerunning NSCF calculations and pw2wannier90.x.
Wannier.truncate — Methodtruncate(model::Model, keep_bands::Vector{Int}, keep_wfs::Vector{Int}=nothing;
orthonorm_U::Bool=true)Truncate U, M, E matrices in model.
Arguments
model: theModelto be truncated.keep_bands: Band indexes to be kept, start from 1.keep_wfs: WF indexes to be kept, start from 1. Ifnothing, keep all.
Keyword arguments
orthonorm_U: If true, Lowdin orthonormalizeUafter truncation. TheUneeds to be (semi-)unitary, so it should always be true.
Wannier.truncate_mmn_eig — Functiontruncate_mmn_eig(seedname, keep_bands::Vector{Int}, outdir="truncate")Truncate number of bands of mmn and eig files.
Arguments
keep_bands: a vector of band indices to keep, starting from 1.outdir: the folder for writingmmnandeigfiles.
This is useful for generating valence only mmn, eig files from a valence+conduction NSCF calculation, so that no need to recompute NSCF with lower number of bands again.
Wannier.truncate_unk — Functiontruncate_unk(dir, keep_bands::Vector{Int}, outdir="truncate"; binary=true)Truncate UNK files for specified bands.
Arguments
dir: folder ofUNKfiles.keep_bands: the band indexes to keep. Start from 1.outdir: folder to write outputUNKfiles.binary: whether to write binaryUNKfiles.
Wannier.truncate_w90 — Functiontruncate_w90(seedname, keep_bands::Vector{Int}, outdir="truncate", unk=false)Truncate mmn, eig, and optionally UNK files.
Arguments
- seedname: seedname for input
mmnandeigfiles. - keep_bands: Band indexes to be kept, start from 1.
- unk: If true also truncate
UNKfiles. - outdir: folder for output files.
3D visualization files
Wannier.read_xsf — Methodread_xsf(filename::AbstractString)Read xsf file.
Return
primvec:3 * 3, Å, each column is a primitive lattice vectorconvvec:3 * 3, Å, each column is a conventional lattice vectoratoms:n_atomsString, atomic symbols or numbersatom_positions:3 * n_atoms, Å, cartesian coordinatesrgrid:RGrid, grid on whichWis definedW:nx * ny * nz, volumetric data
Only support reading 1 datagrid in BLOCK_DATAGRID_3D.
Wannier.write_xsf — Methodwrite_xsf(filename, lattice, atom_positions, atom_numbers, rgrid, W)Write xsf file.
Arguments
lattice:3 * 3, Å, each column is a lattice vectoratom_positions:3 * n_atoms, fractional coordinatesatom_numbers:n_atoms, atomic numbersrgrid:RGridW:nx * ny * nz, volumetric data
This is a more user-friendly version. The rgrid contains the information of the grid origin and spanning vectors.
See also WannierIO.write_xsf
Wannier.read_bxsf — Methodread_bxsf(filename::AbstractString)Read bxsf file.
Return
rgrid:RGrid, grid on whichEis definedfermi_energy: in eVE:n_bands * n_kx * n_ky * n_kz, energy eigenvalues
Wannier.write_bxsf — Methodwrite_bxsf(filename, lattice, atom_positions, atom_numbers, rgrid, W)Write bxsf file.
Arguments
rgrid:RGridfermi_energy: in eVE:n_bands * n_kx * n_ky * n_kz, energy eigenvalues
This is a more user-friendly version. The rgrid contains the information of the grid origin and spanning vectors.
See also WannierIO.write_bxsf
Wannier.read_cube — Methodread_cube(filename::AbstractString)Read cube file.
By default, cube use Bohr unit, here all returns are in Cartesian coordinates, Å unit.
Wannier.write_cube — Methodwrite_cube(filename, lattice, atom_positions, atom_numbers, wf_center, rgrid, W; radius=4.0)Write cube file for WF.
Arguments
lattice: each column is a lattice vector, Åatom_positions:3 * n_atoms, fractional coordinatesatom_numbers:n_atoms, atomic numberswf_centers:3, fractional coordinates w.r.t. latticergrid:RGridW:nx * ny * nz, volumetric data
Keyword arguments
radius: Å. Periodic replica of atoms are written only for the region withinradiusÅ fromwf_center.
See also write_cube(filename, filename, atom_positions, atom_numbers, origin, span_vectors, W)
Interface to DFT codes
Wannier.compute_mud — Methodcompute_mud(dir_up, dir_dn)Compute the overlap matrix between spin up and down from UNK files.
\[M^{\uparrow\downarrow}_{m n \bm{k}} = \langle u^{\uparrow}_{m \bm{k}} | u^{\downarrow}_{n \bm{k}} \rangle\]
This function compute mud matrix in real space (thus much slower), checked with QE pw2wannier90.x function which computes mud in reciprocal space, the difference is on the order of 1e-13.
This only works for norm-conserving pseudopotentials since in that case the overlap operator is just identity; for ultrasoft or PAW, a simple dot product is not enough. Also I assume the UNK files are written with a normalization factor of $N$ (total points of the FFT grid) over the unit cell, i.e., the UNK generated by QE pw2wannier90.x. If the UNK files are normalized to 1, the result should be multiplied by $N$.
Arguments
dir_up: directory of spin up UNK filesdir_dn: directory of spin down UNK files