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.

Tip

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.

Tip

The following abbreviations are used throughout the code and documentation:

  • n_bands for number of bands
  • n_wann for number of WFs
  • n_kpts for number of kpoints
  • n_bvecs for number of b-vectors
  • n_atoms for number of atoms
  • U for amn or the gauge matrices
  • M for mmn matrices
  • E for eig matrices
  • S for spn matrices
Note

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

Wannier90 files

Brillouin.KPaths.KPathInterpolantMethod
KPathInterpolant(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.
source
Wannier.write_w90_kpt_labelMethod
write_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)
source
Wannier.ModelMethod
Model(chk::Chk)

Construct a model from a WannierIO.Chk struct.

Arguments

  • chk: a WannierIO.Chk struct

Keyword Arguments

  • E: a n_wann * n_kpts array for the eigenvalues of the Hamiltonian. If not provided, it will be set to zero.
  • U: a n_wann * n_wann * n_kpts array for the gauge transformation.
Warning

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.

source
Wannier.write_chkMethod
write_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 .chk file
  • model: a Model struct
  • U: n_bands * n_wann * n_kpts array for the gauge transformation

Keyword Arguments

  • exclude_bands: a list of band indices to exclude. This is irrelevant to the model, but chk file has this entry.
  • binary: whether to write the .chk file in binary format.
source
Wannier.write_chkMethod
write_chk(filename, model; exclude_bands=nothing, binary=false)

Write the model to a Wannier90 .chk file.

Arguments

  • filename: filename of the .chk file
  • model: a Model struct

Keyword Arguments

  • exclude_bands: a list of band indices to exclude. This is irrelevant to the model, but chk file has this entry.
  • binary: whether to write the .chk file in binary format.
source
Wannier.read_w90Method
read_w90(seedname::AbstractString; amn=true, orthonorm_amn=true, mmn=true, eig=true)

Read win, and optionally amn, mmn, eig.

Keyword arguments

  • amn: if true, read amn file
  • orthonorm_amn: Lowdin orthonormalization after reading amn matrices. Should be true for most cases, since usually the input amn matrices are projections onto atomic orbitals, and are not unitary or semi-unitary.
  • mmn: if true, read mmn file
  • eig: if true, read eig file
source
Wannier.read_w90_interpMethod
read_w90_interp(seedname::AbstractString; chk=true, amn=nothing, mdrs=nothing)

Return an InterpModel for Wannier interpolation.

Keyword arguments

  • chk: if true, read chk or chk.fmt file to get the unitary matrices, otherwise read amn file for unitary matrices.
  • amn: filename for amn, read this amn file for unitary matrices. Only used if not reading chk and amn is given.
  • mdrs: use MDRS interpolation, else Wigner-Seitz interpolation. If is nothing, detect from win file; and if no use_ws_distance in win file, default to true.
source
Wannier.write_w90Method
write_w90(seedname::AbstractString, model::Model)

Write Model into eig, mmn, amn files.

Keyword arguments

  • binary: write eig, mmn, and amn in Fortran binary format
source
Wannier.read_nnkpMethod
read_nnkp(filename::AbstractString)

Read the nnkp file.

Note

This is a wrapper of WannierIO.read_nnkp. It returns a BVectors instead of NamedTuple.

source
Wannier.write_nnkpFunction
write_nnkp(filename::AbstractString, bvectors::BVectors, n_wann::Integer)

Write nnkp that can be used by pw2wannier90.

Note

This is a wrapper of WannierIO.write_nnkp.

source
Wannier._read_w90_tbMethod
_read_w90_tb(seedname::AbstractString)

Read seedname_tb.dat and seedname_wsvec.dat.

Return

  • R vectors
  • Hamiltonian
  • position operator
source
Wannier.read_w90_tbMethod
read_w90_tb(seedname; kpoints, atom_positions, atom_labels)

Read seedname_tb.dat and seedname_wsvec.dat, return an InterpModel.

Arguments

  • seedname: the seedname of seedname_tb.dat and seedname_wsvec.dat

Keyword Arguments

  • kpoints: the kpoints used in Wannierization, each column is a fractional coordinate.
  • atom_positions: columns are the fractional coordinates of atoms
  • atom_labels: labels of atoms

Return

Note

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.

Warning

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.

source

File manipulation

Truncate Wannier90 matrices

Tip

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.truncateMethod
truncate(model::Model, keep_bands::Vector{Int}, keep_wfs::Vector{Int}=nothing;
    orthonorm_U::Bool=true)

Truncate U, M, E matrices in model.

Arguments

  • model: the Model to be truncated.
  • keep_bands: Band indexes to be kept, start from 1.
  • keep_wfs: WF indexes to be kept, start from 1. If nothing, keep all.

Keyword arguments

  • orthonorm_U: If true, Lowdin orthonormalize U after truncation. The U needs to be (semi-)unitary, so it should always be true.
source
Wannier.truncate_mmn_eigFunction
truncate_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 writing mmn and eig files.
Tip

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.

source
Wannier.truncate_unkFunction
truncate_unk(dir, keep_bands::Vector{Int}, outdir="truncate"; binary=true)

Truncate UNK files for specified bands.

Arguments

  • dir: folder of UNK files.
  • keep_bands: the band indexes to keep. Start from 1.
  • outdir: folder to write output UNK files.
  • binary: whether to write binary UNK files.
source
Wannier.truncate_w90Function
truncate_w90(seedname, keep_bands::Vector{Int}, outdir="truncate", unk=false)

Truncate mmn, eig, and optionally UNK files.

Arguments

  • seedname: seedname for input mmn and eig files.
  • keep_bands: Band indexes to be kept, start from 1.
  • unk: If true also truncate UNK files.
  • outdir: folder for output files.
source

3D visualization files

Wannier.read_xsfMethod
read_xsf(filename::AbstractString)

Read xsf file.

Return

  • primvec: 3 * 3, Å, each column is a primitive lattice vector
  • convvec: 3 * 3, Å, each column is a conventional lattice vector
  • atoms: n_atoms String, atomic symbols or numbers
  • atom_positions: 3 * n_atoms, Å, cartesian coordinates
  • rgrid: RGrid, grid on which W is defined
  • W: nx * ny * nz, volumetric data
Note

Only support reading 1 datagrid in BLOCK_DATAGRID_3D.

source
Wannier.write_xsfMethod
write_xsf(filename, lattice, atom_positions, atom_numbers, rgrid, W)

Write xsf file.

Arguments

  • lattice: 3 * 3, Å, each column is a lattice vector
  • atom_positions: 3 * n_atoms, fractional coordinates
  • atom_numbers: n_atoms, atomic numbers
  • rgrid: RGrid
  • W: 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

source
Wannier.read_bxsfMethod
read_bxsf(filename::AbstractString)

Read bxsf file.

Return

  • rgrid: RGrid, grid on which E is defined
  • fermi_energy: in eV
  • E: n_bands * n_kx * n_ky * n_kz, energy eigenvalues
source
Wannier.write_bxsfMethod
write_bxsf(filename, lattice, atom_positions, atom_numbers, rgrid, W)

Write bxsf file.

Arguments

  • rgrid: RGrid
  • fermi_energy: in eV
  • E: 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

source
Wannier.read_cubeMethod
read_cube(filename::AbstractString)

Read cube file.

Note

By default, cube use Bohr unit, here all returns are in Cartesian coordinates, Å unit.

source
Wannier.write_cubeMethod
write_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 coordinates
  • atom_numbers: n_atoms, atomic numbers
  • wf_centers: 3, fractional coordinates w.r.t. lattice
  • rgrid: RGrid
  • W: nx * ny * nz, volumetric data

Keyword arguments

  • radius: Å. Periodic replica of atoms are written only for the region within radius Å from wf_center.

See also write_cube(filename, filename, atom_positions, atom_numbers, origin, span_vectors, W)

source

Interface to DFT codes

Wannier.compute_mudMethod
compute_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.

Warning

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 files
  • dir_dn: directory of spin down UNK files
source