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. KspaceStencil, 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
Wannier.ModelWannier._raw_read_w90_tbWannier.compute_mudWannier.generate_w90_kpoint_pathWannier.get_symm_point_indices_labelsWannier.has_cubic_neighborsWannier.read_amn_orthoWannier.read_bxsfWannier.read_cubeWannier.read_nnkp_compute_bweightsWannier.read_w90Wannier.read_w90_tbWannier.read_w90_with_chkWannier.read_xsfWannier.truncateWannier.truncate_mmn_eigWannier.truncate_unkWannier.truncate_w90Wannier.write_bxsfWannier.write_cubeWannier.write_nnkpWannier.write_nnkp_cubicWannier.write_w90Wannier.write_w90_kpt_labelWannier.write_w90_tbWannier.write_xsfWannierIO.read_nnkpWannierIO.read_w90_bandWannierIO.write_chkWannierIO.write_w90_band
Wannier90 files
Wannier.read_amn_ortho — Methodread_amn_ortho(filename)
Read amn file and Lowdin orthonormalize the unitary matrices.
The U matrix for Wannier functions must be unitary or semi-unitary. Thus, in most cases, this function should be used instead of WannierIO.read_amn, where the latter one just parse the amn file and return whatever is in it.
Wannier.generate_w90_kpoint_path — Methodgenerate_w90_kpoint_path(
recip_lattice,
kpoints,
symm_point_indices,
symm_point_labels
)
Generate a KPathInterpolant from the kpoint coordinates and high-symmetry kpoints.
The WannierIO.read_w90_band(prefix) function returns the required arguments of this function.
Arguments
recip_lattice: each column is a reciprocal lattice vector in Åkpoints: list of kpoint coordinates along a kpath, fractional coordinatessymm_point_indices: indices of the high-symmetry kpointssymm_point_labels: labels of the high-symmetry kpoints
Wannier.get_symm_point_indices_labels — Methodget_symm_point_indices_labels(kpi)
Return the symmetry indices and labels.
Wannier.write_w90_kpt_label — Methodwrite_w90_kpt_label(prefix, kpi)
Write kpoints into wannier90 formats: prefix_band.kpt, prefix_band.labelinfo.dat.
Arguments
prefix: the prefix of the output filesprefix_band.kptandprefix_band.labelinfo.datkpi: aKPathInterpolantobject
This allows auto generating the high-symmetry kpoints and labels from crystal structure using generate_kpath, writing them into files. Then other codes can use the kpoints for band structure calculations, e.g., QE pw.x bands calculation, or in the win input file for Wannier90.
Example
win = read_win("si2.win")
kp = generate_kpath(win.unit_cell_cart, win.atoms_frac, win.atom_labels)
kpi = Wannier.generate_w90_kpoint_path(kp, 100)
Wannier.write_w90_kpt_label("si2", kpi)WannierIO.read_w90_band — Methodread_w90_band(prefix, recip_lattice)
Arguments
recip_lattice: each column is a reciprocal lattice vector in Å.
Return
kpi: aKPathInterpolantobjecteigenvalues: the eigenvalues of the band structure
WannierIO.write_w90_band — Methodwrite_w90_band(prefix, kpi, eigenvalues)
Write prefix_band.dat, prefix_band.kpt, prefix_band.labelinfo.dat.
This is a more user-friendly version that works with KPathInterpolant; the WannierIO.write_w90_band(prefix; kwargs...) is the low-level version.
Wannier.Model — MethodModel(chk; kmesh_tol)
Construct a Model from a WannierIO.Chk struct.
Arguments
chk: aWannierIO.Chkstruct
The Chk struct does not contain eigenvalues, and the Model.E will be set to zeros.
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.
WannierIO.write_chk — Functionwrite_chk(filename, model; ...)
write_chk(
filename,
model,
gauges;
exclude_bands,
binary,
header
)
Write a Model to a wannier90 chk file.
Arguments
filename: filename of the.chkfilemodel: aModelstructU: gauge transformation matrices, default ismodel.U.
Keyword Arguments
exclude_bands: a list of band indices to be excluded.binary: write thechkfile in binary formatheader: header of thechkfile
The exclude_bands is totally not used throughout the code. However, if one want to use wannier90 to restart from the written chk file, the exclude_bands must be consistent with that in wannier90 win file. The wannier90 exclude_bands input parameter in only used in the pw2wannier90 step to remove some bands (usually semicore states) when computing amn/mmn/eig files, therefore this exclude_bands is totally irrelevant to our Model struct.
Wannier.read_w90 — Methodread_w90(prefix; ortho_amn, use_mmn_bvecs, kstencil_algo)
Read win and mmn files, and read amn/eig files if they exist.
Keyword arguments
- ortho_amn: Lowdin orthonormalization after reading
amnmatrices. Should betruefor most cases, since usually the inputamnmatrices are not guaranteed to be unitary or semi-unitary. - usemmnbvecs: use the b-vectors in
mmnfile instead of regenerating them. - kstencilalgo: algorithm to generate
KspaceStencilif `usemmnbvecsisfalse. Default isgeneratekspace_stencil`.
Wannier.read_w90_with_chk — Functionread_w90_with_chk(prefix)
read_w90_with_chk(prefix, chk)
Return a Model with U matrix filled by that from a chk file.
Arguments
- chk: path of chk file to get the unitary matrices.
Wannier.write_w90 — Methodwrite_w90(prefix, model; binary)
Write Model into eig, mmn, amn files.
Keyword arguments
binary: writeeig,mmn, andamnin Fortran binary format
Wannier.has_cubic_neighbors — Methodhas_cubic_neighbors(kpoints, kpb_k, kpb_G; kgrid_size, atol)
Check if the kpoint stencil (bvectors) contain 6 cubic neighbors.
parallel_transport requires the 6 nearest neighbors.
Arguments
kpoints: the kpoints in fractional coordinateskpb_k: indices of kpoint (inkpoints) that is translational equivalent to k+bkpb_G: translation vector that maps k to k+b
Wannier.read_nnkp_compute_bweights — Methodread_nnkp_compute_bweights(filename)
Read the nnkp file.
This function calls WannierIO.read_nnkp to parse the file, compute the bweights of b-vectors, and returns a KspaceStencil (while WannierIO.read_nnkp only returns a NamedTuple).
Wannier.write_nnkp — Methodwrite_nnkp(filename, kstencil; kwargs...)
Write nnkp that can be used by pw2wannier90.
Arguments
filename: the filename to write tokstencil: aKspaceStencilobject
Some important tags in nnkp file (can be passed as keyword arguments):
n_wann: the number of WFs, needed bypw2wannier90exclude_bands: the bands (often semicore states) to exclude, needed bypw2wannier90
For other keyword arguments, see WannierIO.write_nnkp.
Wannier.write_nnkp_cubic — Methodwrite_nnkp_cubic(filename, win)
Write a nnkp file with 6 cubic neighbors. Useful for parallel_transport.
nnkp file also contains the exclude_bands parameters, which will affect the mmn file computed by pw2wannier90.x.
Arguments
filename: the filename of the newnnkpfilewin: the Wannier90 parameters, e.g. returned byread_win
WannierIO.read_nnkp — Methodread_nnkp(filename, weights)
Arguments
filename: the filename of thennkpfileweights: the weights of the b-vectors. Ifnothing, the weights are set to 0.0. This is useful forparallel_transportusing 6 cubic neighbors, where the the $b$-vectors are not complete, thus not possible to compute the weights.
Wannier._raw_read_w90_tb — MethodOnly read tb files, without further processing
Wannier.read_w90_tb — Methodread_w90_tb(prefix)
Read prefix_tb.dat and prefix_wsvec.dat and construct tight-binding models.
Arguments
prefix: the prefix ofprefix_tb.datandprefix_wsvec.dat
Return
This will call simplify to absorb the R-vector degeneracies and T-vectors into the operator, leading to faster interpolations.
Wannier.write_w90_tb — Methodwrite_w90_tb(prefix, hamiltonian, position)
Write a tight-binding model of Hamiltonian and position operator into prefix_tb.dat and prefix_wsvec.dat files.
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.
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
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.
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