Wannier90 files

Warning

Some of the functions, e.g. read_amn, write_amn, read_chk, write_chk, etc., support reading/writing Fortran unformatted files. However, the Fortran unformatted files are machine and compiler dependent. Therefore, it is not guaranteed that these functions work for all the cases. Currently, the functions are tested on the following platforms:

  • Linux + gfortran 11.2

Public API

WannierIO.read_amnFunction
read_amn(filename)
read_amn(filename, ::FortranText)
read_amn(filename, ::FortranBinaryStream)

Read wannier90 amn file.

Return

  • A: length-n_kpts vector, each element is a n_bands * n_wann matrix.
  • header: first line of the file

Note there are three versions of this function: the 1st one is a wrapper function that automatically detect the format (text or binary) of the file, and does some additional pretty printing to give user a quick hint of the dimensions of the A matrix; it internally calls the 2nd or the 3rd version for actual reading.

Wannier90 only has Fortran text format for amn, however I wrote a custom version of QE pw2wannier90.x that can output Fortran binary format (using Fortran stream IO) to save some disk space. The 1st function auto detect the file format so it is transparent to the user.

source
WannierIO.write_amnFunction
write_amn(filename, A; header=default_header(), binary=false)
write_amn(filename, A, ::FortranText; header=default_header())
write_amn(filename, A, ::FortranBinaryStream; header=default_header())

Write wannier90 amn file.

Arguments

  • filename: output filename
  • A: a length-n_kpts vector, each element is a n_bands * n_wann matrix

Keyword arguments

  • header: 1st line of the file
  • binary: write as Fortran unformatted file, which is the Wannier90 default. Here the binary kwargs is provided for convenience.

Same as read_amn there are three versions of this function, the 1st one is a wrapper function, it calls the 2nd or the 3rd version depending on the binary kwargs.

source
WannierIO.read_w90_bandMethod
read_w90_band(prefix)

Read prefix_band.dat, prefix_band.kpt, prefix_band.labelinfo.dat.

Arguments

  • prefix: prefix of the filenames (or called seedname in wannier90), NOT the full filename.

Return

  • x: n_kpts, x axis value of kpath, in cartesian length
  • eigenvalues: length-n_kpts vector, each element is a length-n_bands vector of band energies
  • kpoints: a vector of length n_kpts, fractional coordinates
  • kweights: a vector of length n_kpts, weights of kpoints
  • symm_point_indices: index of high-symmetry points in prefix_band.dat
  • symm_point_labels: name of high-symmetry points
source
WannierIO.read_w90_band_datMethod
read_w90_band_dat(filename)

Read prefix_band.dat file generated by wannier90.x, or prefix-band.dat file generated by postw90.x.

Return

  • x: n_kpts, x axis value of kpath, in cartesian length
  • eigenvalues: length-n_kpts vector, each elemnt is a length-n_bands vector of band energies
  • extras: optional (the postw90.x might generate a file with a third column), same size as eigenvalues, often the color of each eigenvalue, e.g., spin projection
source
WannierIO.read_w90_band_kptMethod
read_w90_band_kpt(filename)

Read a prefix_band.kpt file.

Return

  • kpoints: a vector of length n_kpts, fractional coordinates
  • kweights: a vector of length n_kpts, weights of kpoints
source
WannierIO.read_w90_band_labelinfoMethod
read_w90_band_labelinfo(filename)

Read prefix_band.labelinfo file.

Return

  • symm_point_indices: index of high-symmetry points in prefix_band.dat
  • symm_point_labels: name of high-symmetry points
source
WannierIO.write_w90_bandMethod
write_w90_band(
    prefix;
    x,
    eigenvalues,
    kpoints,
    kweights,
    symm_point_indices,
    symm_point_labels
)

Write prefix_band.dat, prefix_band.kpt, prefix_band.labelinfo.dat.

Arguments

  • prefix: prefix of prefix_band.dat, prefix_band.kpt, prefix_band.labelinfo.dat

Keyword Arguments

  • x: n_kpts, x axis value, in cartesian length
  • eigenvalues: length-n_kpts vector, each element is a length-n_bands vector of band energies
  • kpoints: length-n_kpts vector, fractional coordinates
  • kweights: a vector of length n_kpts, weights of kpoints
  • symm_point_indices: index of high-symmetry points in prefix_band.dat
  • symm_point_labels: name of high-symmetry points
source
WannierIO.write_w90_band_datMethod
write_w90_band_dat(filename; x, eigenvalues, extras)

Write prefix_band.dat file.

Arguments

  • filename: filename of prefix_band.dat

Keyword Arguments

  • x: n_kpts, x axis value, in cartesian length
  • eigenvalues: length-n_kpts vector, each element is a length-n_bands vector of band energies
  • extras: optional, same size as eigenvalues, will be written as the third column of prefix_band.dat. The prefix-band.dat file generated by postw90.x sometimes has a third column for e.g. the color of the eigenvalues
source
WannierIO.write_w90_band_kptMethod
write_w90_band_kpt(filename; kpoints, kweights)

Write prefix_band.kpt file.

Arguments

  • filename: filename of prefix_band.kpt

Keyword Arguments

  • kpoints: length-n_kpts vector, fractional coordinates
  • kweights: n_kpts, optional, weights of kpoints, default to 1.0.
source
WannierIO.write_w90_band_labelinfoMethod
write_w90_band_labelinfo(
    filename;
    x,
    kpoints,
    symm_point_indices,
    symm_point_labels
)

Write prefix_band.labelinfo file.

Arguments

  • filename: filename of prefix_band.labelinfo

Keyword Arguments

  • x: n_kpts-vector, x axis value, in cartesian length
  • kpoints: length-n_kpts vector, fractional coordinates
  • symm_point_indices: index of high-symmetry points in prefix_band.dat
  • symm_point_labels: name of high-symmetry points
source
WannierIO.read_chkMethod
read_chk(filename)
read_chk(filename, ::FortranText)
read_chk(filename, ::FortranBinary)

Read wannier90 chk checkpoint file.

Similar to read_amn, the 1st version auto detect chk file format (binary or text) and read it.

source
WannierIO.write_chkFunction
write_chk(filename, chk::Chk; binary=false)
write_chk(filename, chk::Chk, ::FortranText)
write_chk(filename, chk::Chk, ::FortranBinary)

Write wannier90 chk file.

Similar to write_amn, the 1st version is a convenience wrapper.

Keyword arguments

  • binary: write as Fortran binary file or not. Although wannier90 default is Fortran binary format, here the default is false since Fortran binary depends on compiler and platform, so it is not guaranteed to always work.
source
WannierIO.read_eigFunction
read_eig(filename)
read_eig(filename, ::FortranText)
read_eig(filename, ::FortranBinaryStream)

Read the wannier90 eig file.

Return

  • eigenvalues: a lenth-n_kpts vector, each element is a length-n_bands vector

The 1st version is a convenience wrapper function for the 2nd and 3rd versions.

source
WannierIO.write_eigFunction
write_eig(filename, eigenvalues; binary=false)
write_eig(filename, eigenvalues, ::FortranText)
write_eig(filename, eigenvalues, ::FortranBinaryStream)

Write eig file.

Arguments

  • eigenvalues: a length-n_kpts vector, each element is a length-n_bands vector

Keyword arguments

  • binary: if true write in Fortran binary format.
source
WannierIO.read_w90_hrdatMethod
read_w90_hrdat(filename)

Read prefix_hr.dat.

Return

  • Rvectors: $\mathbf{R}$-vectors on which operators are defined
  • Rdegens: degeneracies of each $\mathbf{R}$-vector
  • H: Hamiltonian $\mathbf{H}(\mathbf{R})$
  • header: the first line of the file
source
WannierIO.read_mmnFunction
read_mmn(filename)
read_mmn(filename, ::FortranText)
read_mmn(filename, ::FortranBinaryStream)

Read wannier90 mmn file.

Return

  • M: length-n_kpts vector, each element is a length-n_bvecs vector, then each element is a n_bands * n_bands matrix
  • kpb_k: length-n_kpts vector, each element is a length-n_bvecs vector of integers for the indices of the neighboring kpoints
  • kpb_G: length-n_kpts vector, each element is a lenght-n_bvecs vector of of Vec3{Int}, which are the translation vectors
  • header: 1st line of the file

The translation vector G is defined as b = kpoints[kpb_k[ik][ib]] + kpb_G[ik][ib] - kpoints[ik], where b is the ib-th bvector of the ik-th kpoint.

The 1st version is a convenience wrapper for the other two.

source
WannierIO.write_mmnFunction
write_mmn(filename; M, kpb_k, kpb_G; header=default_header(), binary=false)
write_mmn(filename; M, kpb_k, kpb_G, ::FortranText; header=default_header(), binary=false)
write_mmn(filename; M, kpb_k, kpb_G, ::FortranBinaryStream; header=default_header(), binary=false)

Write wannier90 mmn file.

Arguments

  • filename: output file name
  • M: length-n_kpts vector of n_bands * n_bands * n_bvecs arrays
  • kpb_k: length-n_kpts vector of length-n_bvecs vector of integers
  • kpb_G: length-n_kpts vector of length-n_bvecs vector of Vec3{Int} for bvectors

Keyword arguments

  • header: header string
  • binary: if true write in Fortran binary format

The 1st version is a convenience wrapper for the other two.

source
WannierIO.read_nnkpFunction
read_nnkp(filename)
read_nnkp(filename, ::Wannier90Text)
read_nnkp(filename, ::Wannier90Toml)

Read wannier90 nnkp file.

Return

  • lattice: each column is a lattice vector
  • recip_lattice: each column is a reciprocal lattice vector
  • kpoints: length-n_kpts vector, each element is Vec3, in fractional coordinates
  • kpb_k: length-n_kpts vector, each element is a length-n_bvecs vector of integers, index of kpoints
  • kpb_G: length-n_kpts vector, each element is a length-n_bvecs vector, then each element is Vec3 for translations, fractional w.r.t recip_lattice

Wannier90 nnkp file is a plain text format, the 2nd version reads nnkp file in Wannier90 format. The thrid version read a TOML-format nnkp file, which is defined by this package, see write_nnkp. The 1st version auto detects the format and parse it.

source
WannierIO.write_nnkpFunction
write_nnkp(filename; toml=false, kwargs...)
write_nnkp(filename, ::Wannier90Text; kwargs...)
write_nnkp(filename, ::Wannier90Toml; kwargs...)

Write a nnkp file that can be used by DFT codes, e.g., QE pw2wannier90.

Keyword Arguments

  • toml: write to a TOML file, otherwise write to a Wannier90 text file format
  • recip_lattice: each column is a reciprocal lattice vector
  • kpoints: length-n_kpts vector of Vec3, in fractional coordinates
  • kpb_k: length-n_kpts vector, each element is a length-n_bvecs vector of integers, index of kpoints
  • kpb_G: length-n_kpts vector, each element is a length-n_bvecs vector, then each element is a Vec3 for translation vector, fractional w.r.t. recip_lattice
  • n_wann: if given, write an auto_projections block
  • exclude_bands: if given, write the specified band indices in the exclude_bands block
  • header: first line of the file
Note

Only use auto_projections, the projections block is not supported.

source
WannierIO.read_w90_rdatMethod
read_w90_rdat(filename)

Read prefix_r.dat.

Return

  • Rvectors: $\mathbf{R}$-vectors on which operators are defined
  • r_x: $x$-component of position operator
  • r_y: $y$-component of position operator
  • r_z: $z$-component of position operator
  • header: the first line of the file
source
WannierIO.read_spnFunction
read_spn(filename)
read_spn(filename, ::FortranText)
read_spn(filename, ::FortranBinary)

Read the wannier90 spn file.

Return

  • Sx: spin x, a length-n_kpts vector, each element is a n_bands-by-n_bands matrix
  • Sy: spin y, a length-n_kpts vector, each element is a n_bands-by-n_bands matrix
  • Sz: spin z, a length-n_kpts vector, each element is a n_bands-by-n_bands matrix
  • header: 1st line of the file
source
WannierIO.write_spnFunction
write_spn(filename, Sx, Sy, Sz; binary=false, header)
write_spn(filename, Sx, Sy, Sz, ::FortranText; header)
write_spn(filename, Sx, Sy, Sz, ::FortranBinary; header)

Write the spn file.

source
WannierIO.read_w90_tbdatMethod
read_w90_tbdat(filename)

Read prefix_tb.dat.

Return

  • lattice: each column is a lattice vector in Å
  • Rvectors: $\mathbf{R}$-vectors on which operators are defined
  • Rdegens: degeneracies of each $\mathbf{R}$-vector
  • H: Hamiltonian $\mathbf{H}(\mathbf{R})$
  • r_x: $x$-component of position operator
  • r_x: $y$-component of position operator
  • r_x: $z$-component of position operator
  • header: the first line of the file
source
WannierIO.read_uHuFunction
read_uHu(filename)
read_uHu(filename, ::FortranText; transpose_band_indices=true)
read_uHu(filename, ::FortranBinary; transpose_band_indices=true)

Read the wannier90 uHu file.

Keyword Arguments

  • transpose_band_indices: QE pw2wannier90.x writes the matrix in a strange transposed manner; if reading a QE-generated uHu file, this flag should be true to restore the band indices order, so that the returned matrix has the correct order, i.e., uHu[ik][ib1, ib2][m, n] is $\langle u_{m, k + b_1}| H | u_{n, k + b_2} \rangle$

Return

  • uHu: a length-n_kpts vector, each element is a n_bvecs * n_bvecs matrix, then each element is a n_bands * n_bands matrix
  • header: 1st line of the file
source
WannierIO.write_uHuFunction
write_uHu(filename, uHu; binary=false, header)
write_uHu(filename, uHu, ::FortranText; header)
write_uHu(filename, uHu, ::FortranBinary; header)

Write the uHu file.

Keyword Arguments

source
WannierIO.read_unkFunction
read_unk(filename)
read_unk(filename, ::FortranText)
read_unk(filename, ::FortranBinary)

Read wannier90 UNK file for the periodic part of Bloch wavefunctions.

Return

  • ik: k-point index, start from 1
  • Ψ: periodic part of Bloch wavefunctions in real space, size = (n_gx, n_gy, n_gz, n_bands, n_spin)
source
WannierIO.write_unkFunction
write_unk(filename, ik, Ψ; binary=false)
write_unk(filename, ik, Ψ, ::FortranText)
write_unk(filename, ik, Ψ, ::FortranBinary)

Write UNK file for the periodic part of Bloch wavefunctions.

Arguments

  • ik: at which kpoint? start from 1
  • Ψ: Bloch wavefunctions, size(Ψ) = (n_gx, n_gy, n_gz, n_bands, n_spin)

Keyword arguments

  • binary: write as Fortran unformatted file
source
WannierIO.read_winFunction
read_win(filename; fix_inputs=true)
read_win(filename, ::Wannier90Text; fix_inputs=true)
read_win(filename, ::Wannier90Toml; fix_inputs=true)

Read wannier90 input win file.

Arguments

  • filename: The name of the input file.

Keyword Arguments

  • fix_inputs: sanity check and fix the input parameters, e.g., set num_bands = num_wann if num_bands is not specified, convert atoms_cart always to atoms_frac, etc. See also fix_win!.
source
WannierIO.write_winFunction
write_win(filename; toml=false, kwargs...)
write_win(filename, ::Wannier90Text; kwargs...)
write_win(filename, ::Wannier90Toml; kwargs...)

Write input parameters into a wannier90 win file.

The input parameters are keyword arguments, with key names same as that of wannier90.

Examples

using WannierIO

write_win(
    "silicon.win";
    num_wann=4,
    num_bands=4,
    # unit_cell_cart is a matrix, its columns are the lattice vectors in angstrom
    unit_cell_cart=[
        0.0      2.71527  2.71527
        2.71527  0.0      2.71527
        2.71527  2.71527  0.0
    ],
    # atoms_frac is a vector of pairs of atom_label and fractional coordinates
    atoms_frac=[
        :Si => [0.0, 0.0, 0.0],
        :Si => [0.25, 0.25, 0.25],
        # both `:Si` and `"Si"` are allowed
        # "Si" => [0.25, 0.25, 0.25],
    ],
    # each element in projections will be written as a line in the win file
    projections=[
        "random",
    ]
    kpoint_path=[
        [:G => [0.0, 0.0, 0.0], :X => [0.5, 0.0, 0.5]],
        [:X => [0.5, 0.0, 0.5], :U => [0.625, 0.25, 0.625]],
    ],
    mp_grid=[2, 2, 2],
    # kpoints is a matrix, its columns are the fractional coordinates
    kpoints=[
        [0.0, 0.0, 0.0],
        [0.0, 0.0, 0.5],
        [0.0, 0.5, 0.0],
        [0.0, 0.5, 0.5],
        [0.5, 0.0, 0.0],
        [0.5, 0.0, 0.5],
        [0.5, 0.5, 0.0],
        [0.5, 0.5, 0.5],
    ],
    # additional parameters can be passed as keyword arguments, e.g.,
    num_iter=500,
)
source
WannierIO.read_woutMethod
read_wout(filename)

Parse wannire90 wout file.

Return

  • lattice: each column is a lattice vector in Å
  • recip_lattice: each column is a reciprocal lattice vector in Å⁻¹
  • atom_labels: atomic symbols
  • atom_positions: in fractional coordinates
  • centers: final each WF centers in Å
  • spreads: final each WF spreads in Ų
  • ΩI, ΩD, ΩOD, Ωtotal: final spread (components) in Ų
source
WannierIO.read_w90_wsvecMethod
read_w90_wsvec(filename)

Read prefix_wsvec.dat.

Return

  • mdrs: whether use MDRS interpolation, i.e. the use_ws_distance in the header
  • Rvectors: the $\mathbf{R}$-vectors
  • Tvectors: the $\mathbf{T}_{m n \mathbf{R}}$-vectors. Returned only mdrs = true.
  • Tdegens: the degeneracies of $\mathbf{T}_{m n \mathbf{R}}$-vectors. Returned only mdrs = true.
  • n_wann: number of WFs
  • header: the first line of the file
source
WannierIO.write_w90_wsvecMethod
write_w90_wsvec(
    filename;
    Rvectors,
    n_wann,
    Tvectors,
    Tdegens,
    header
)

Write prefix_wsvec.dat.

Keyword Arguments

  • n_wann: for Wigner-Seitz Rvectors, needs to provide a n_wann for number of Wannier functions; for MDRS Rvectors, the n_wann is optional and can be automatically determined from the Tvectors
  • Tvectors and Tdegens: if provided, write in MDRS format; otherwise, write in Wigner-Seitz format

Also see the return values of read_w90_wsvec.

source

Private API

These are some lower-level types/functions that are (probably) less used, thus not exported. Of course, you can still use them by prefixing WannierIO., e.g., WannierIO.read_w90_band_dat(filename).

WannierIO.ChkType

Struct for storing matrices in prefix.chk file.

struct Chk{T<:Real}

One-to-one mapping to the wannier90 chk file, but renaming the variable names so that they are consistent with the rest of the code.

Fields

  • header: The header line, usually contains date and time

  • n_bands: number of bands, can be auto set in constructor according to dimensions of other variables

  • n_exclude_bands: number of excluded bands, can be auto set in constructor

  • exclude_bands: Indices of excluded bands, starts from 1. Vector of integers, size: n_exclude_bands

  • lattice: Matrix of size 3 x 3, each column is a lattice vector in Å unit

  • recip_lattice: Matrix of size 3 x 3, each column is a reciprocal lattice vector in Å⁻¹ unit

  • n_kpts: number of kpoints, can be auto set in constructor

  • kgrid: dimensions of kpoint grid, 3 integers

  • kpoints: kpoint coordinates, fractional, length-n_kpts vector

  • n_bvecs: number of b-vectors, can be auto set in constructor

  • n_wann: number of Wannier functions, can be auto set in constructor

  • checkpoint: a string to indicate the current step (after disentanglement, after maximal localization, ...) in wannier90

  • have_disentangled: Have finished disentanglement or not

  • ΩI: Omega invariant part of MV spreads, in Ų unit

  • dis_bands: Indices of bands taking part in disentanglement, not frozen bands! length-n_kpts vector, each element is a length-n_bands vector of bool.

    This is needed since W90 puts all the disentanglement bands in the first several rows of Udis, (and the first few columns of Udis are the frozen bands) so directly multiplying eigenvalues e.g. (Udis * U)' * diag(eigenvalues) * (Udis * U) is wrong!

  • n_dis: number of bands taking part in disentanglement at each kpoint. can be auto set in constructor from dis_bands

  • Udis: Semi-unitary matrix for disentanglement, length-n_kpts vector, each elment has size: n_bands x n_wann, i.e., the u_matrix_opt in wannier90

  • Uml: Unitary matrix for maximal localization, length-n_kpts vector, each element has size: n_wann x n_wann, i.e., the u_matrix in wannier90. The abbreviation ml stands for maximal localization, so as to differentiate from the (combined) unitary matrix U = Udis * Uml.

  • M: Wannier-gauge overlap matrix, length-n_kpts vector of length-n_bvecs vector, each element is a matrix of size n_wann x n_wann, i.e., the m_matrix in wannier90

  • r: Wannier function centers, length-n_wann vector, Cartesian coordinates in Å unit, i.e., the wannier_centres variable in wannier90

  • ω: Wannier function spreads, length-n_wann vector, Ų unit, i.e., the wannier_spreads variable in wannier90

source
WannierIO.ChkMethod
Chk(
    header,
    exclude_bands,
    lattice,
    recip_lattice,
    kgrid,
    kpoints,
    checkpoint,
    have_disentangled,
    ΩI,
    dis_bands,
    Udis,
    Uml,
    M,
    r,
    ω
)

Convenience constructor of Chk struct that auto set some fields.

source
WannierIO._check_eig_orderMethod
_check_eig_order(eigenvalues; digits)

Check that eigenvalues are in order.

Some times there are small noises, use digits to set the number of digits for comparisons.

source
WannierIO._reshape_eigMethod
_reshape_eig(idx_b, idx_k, eig)

Reshape a vector of eigenvalues into a matrix of eigenvalues.

Auto detect the number of bands and kpoints.

source
WannierIO.write_HH_RMethod
write_HH_R(filename, H, R; N, header)

Write the real space Hamiltonian to a prefix_HH_R.dat file.

Arguments

  • filename: usually prefix_HH_R.dat
  • H: a n_wann * n_wann * n_rvecs array of Hamiltonian
  • R: a n_rvecs * 3 array of integers

Keyword arguments

  • N: a n_rvecs vector of integers, the degeneracy of each R vector
  • header: a string, the header of the file
Note

Wanier90 postw90.x has a hidden input parameter effective_model, setting it to true and postw90.x will read this HH_R.dat to fill the real space Hamiltonian, and do subsequent Wannier interpolation, e.g., in BoltzWann. However, the vanila postw90.x code does not take into account the degeneracy of R vectors, and also does not use MDRS interpolation. I have modified the postw90.x code to use MDRS, and also changed a bit the number of digits for the Hamiltonian in HH_R.dat, so that it is the same as the prefix_tb.dat file, i.e., from Fortran F12.6 to E15.8.

source
WannierIO.read_u_matMethod
read_u_mat(filename)

Read wannier90 prefix_u.mat or prefix_u_dis.mat file.

Arguments

  • filename: the input file name

Return

  • U: Udis (for disentanglement) or U (for maximal localization) matrices
  • kpoints: fractional kpoint coordinates
  • header: 1st line of the file
Warning

The wannier90 output prefix_u_dis.mat internally sorts the band indices according to the disnentanglement window, therefore it can be different from the original Bloch states, see the code and comments in get_Udis.

source
WannierIO.write_u_matMethod
write_u_mat(filename, U, kpoints; header=default_header())

Write wannier90 prefix_u.mat or prefix_u_dis.mat file.

Arguments

  • filename: the input file name
  • U: Udis (for disentanglement) or U (for maximal localization) matrices
  • kpoints: fractional kpoint coordinates

Keyword arguments

  • header: 1st line of the file, optional
Warning

The wannier90 output prefix_u_dis.mat internally sorts the band indices according to the disnentanglement window, therefore it can be different from the original Bloch states, see the code and comments in get_Udis. This function just writes whatever is inside the input U matrix, without consider the order of disentanglement window.

source