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.

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

Wannier.read_amn_orthoMethod
read_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.

source
Wannier.generate_w90_kpoint_pathMethod
generate_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 coordinates
  • symm_point_indices: indices of the high-symmetry kpoints
  • symm_point_labels: labels of the high-symmetry kpoints
source
Wannier.write_w90_kpt_labelMethod
write_w90_kpt_label(prefix, kpi)

Write kpoints into wannier90 formats: prefix_band.kpt, prefix_band.labelinfo.dat.

Arguments

  • prefix: the prefix of the output files prefix_band.kpt and prefix_band.labelinfo.dat
  • kpi: a KPathInterpolant object
Tip

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)
source
WannierIO.read_w90_bandMethod
read_w90_band(prefix, recip_lattice)

Arguments

  • recip_lattice: each column is a reciprocal lattice vector in Å.

Return

  • kpi: a KPathInterpolant object
  • eigenvalues: the eigenvalues of the band structure
Note

The WannierIO.read_w90_band(prefix) function returns a NamedTuple containing basis variables such as kpoints, symm_point_indices, etc. Here, once we know the recip_lattice, we can generate a KPathInterpolant which can be used directly in plotting functions.

source
WannierIO.write_w90_bandMethod
write_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.

source
Wannier.ModelMethod
Model(chk; kmesh_tol)

Construct a Model from a WannierIO.Chk struct.

Arguments

  • chk: a WannierIO.Chk struct
Warning

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.

source
WannierIO.write_chkFunction
write_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 .chk file
  • model: a Model struct
  • U: gauge transformation matrices, default is model.U.

Keyword Arguments

  • exclude_bands: a list of band indices to be excluded.
  • binary: write the chk file in binary format
  • header: header of the chk file
Note

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.

source
Wannier.read_w90Method
read_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 amn matrices. Should be true for most cases, since usually the input amn matrices are not guaranteed to be unitary or semi-unitary.
  • usemmnbvecs: use the b-vectors in mmn file instead of regenerating them.
  • kstencilalgo: algorithm to generate KspaceStencil if `usemmnbvecsisfalse. Default isgeneratekspace_stencil`.
source
Wannier.read_w90_with_chkFunction
read_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.
source
Wannier.write_w90Method
write_w90(prefix, model; binary)

Write Model into eig, mmn, amn files.

Keyword arguments

  • binary: write eig, mmn, and amn in Fortran binary format
source
Wannier.read_nnkp_compute_bweightsMethod
read_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).

source
Wannier.write_nnkpMethod
write_nnkp(filename, kstencil; kwargs...)

Write nnkp that can be used by pw2wannier90.

Arguments

  • filename: the filename to write to
  • kstencil: a KspaceStencil object
Tip

Some important tags in nnkp file (can be passed as keyword arguments):

  • n_wann: the number of WFs, needed by pw2wannier90
  • exclude_bands: the bands (often semicore states) to exclude, needed by pw2wannier90

For other keyword arguments, see WannierIO.write_nnkp.

source
Wannier.read_w90_tbMethod
read_w90_tb(prefix)

Read prefix_tb.dat and prefix_wsvec.dat and construct tight-binding models.

Arguments

  • prefix: the prefix of prefix_tb.dat and prefix_wsvec.dat

Return

Note

This will call simplify to absorb the R-vector degeneracies and T-vectors into the operator, leading to faster interpolations.

source
Wannier.write_w90_tbMethod
write_w90_tb(prefix, hamiltonian, position)

Write a tight-binding model of Hamiltonian and position operator into prefix_tb.dat and prefix_wsvec.dat files.

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