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

Contents

Index

Read/write

WannierIO.read_winMethod
read_win(filename; fix_inputs=true)

Read the input file of Wannier90.

Arguments

  • filename::AbstractString: 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_winMethod
write_win(filename::AbstractString; kwargs...)

Write input parameters into a Wannier90 win file.

The input parameters are keyword arguments, with key names same as that of Wannier90. The only exception is atoms_frac, which contains both atom labels and fractional coordinates; however, here it is split into two keywords: atom_labels which is a vector of element names, and atoms_frac which is a matrix.

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
    ],
    # both [:Si, :Si] and ["Si", "Si"] are allowed
    atom_labels = ["Si", "Si"],
    # atoms_frac is a matrix, its columns are the fractional coordinates
    atoms_frac=[
        0.0  0.25
        0.0  0.25
        0.0  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.5  0.5  0.5  0.5
        0.0  0.0  0.5  0.5  0.0  0.0  0.5  0.5
        0.0  0.5  0.0  0.5  0.0  0.5  0.0  0.5
    ],
    # additional parameters can be passed as keyword arguments, e.g.,
    num_iter=500,
)
source
WannierIO.read_woutMethod
read_wout(filename::AbstractString)

Parse wout file.

Return

  • lattice: in Å, each column is a lattice vector
  • atom_labels: atomic symbols
  • atom_positions: in fractional coordinates, each column is a coordinate vector
  • centers: WF centers in Å, each column is a WF center vector
  • spreads: WF spreads in Å^2
source
WannierIO._read_amn_binMethod
_read_amn_bin(filename::AbstractString)

Read Fortran binary (unformatted) amn file.

Using Fortran stream IO.

source
WannierIO.read_amnMethod
read_amn(filename::AbstractString)

Read the amn file.

Return

  • A: a n_bands * n_wann * n_kpts array.
source
WannierIO.write_amnMethod
write_amn(filename, A::Array{Complex,3}, header; binary=false)

Write amn file.

Arguments

  • filename: output filename
  • A: a n_bands * n_wann * n_kpts array
  • header: 1st line of the file

Keyword arguments

  • binary: write as Fortran unformatted file
source
WannierIO.write_amnMethod
write_amn(filename, A::Array{Complex,3})

Write amn file.

Default header is "Created by WannierIO.jl CURRENT_DATE".

Arguments

  • filename: output filename
  • A: a n_bands * n_wann * n_kpts array

Keyword arguments

  • binary: write as Fortran unformatted file
source
WannierIO.read_mmnMethod
read_mmn(filename::AbstractString)

Read mmn file.

Returns a n_bands * n_bands * n_bvecs * n_kpts array.

source
WannierIO.write_mmnMethod
write_mmn(filename, M::Array{ComplexF64,4}, kpb_k, kpb_b, header; binary=false)

Write mmn file.

Arguments

  • filename: output file name
  • M: n_bands * n_bands * n_bvecs * n_kpts array
  • kpb_k: n_bvecs * n_kpts array
  • kpb_b: 3 * n_bvecs * n_kpts array
  • header: header string

Keyword arguments

  • binary: if true write in Fortran binary format
source
WannierIO.write_mmnMethod
write_mmn(filename, M::Array{ComplexF64,4}, kpb_k, kpb_b; binary=false)

Write mmn file.

Default header is "Created by WannierIO.jl CURRENT_DATE".

Arguments

  • filename: output file name
  • M: n_bands * n_bands * n_bvecs * n_kpts array
  • kpb_k: n_bvecs * n_kpts array
  • kpb_b: 3 * n_bvecs * n_kpts array

Keyword arguments

  • binary: if true write in Fortran binary format
source
WannierIO.read_eigMethod
read_eig(filename::AbstractString)

Read the eig file.

Return

  • E: a n_bands * n_kpts array
source
WannierIO.write_eigMethod
write_eig(filename::AbstractString, E::AbstractArray; binary=false)

Write eig file.

Arguments

  • E: n_bands * n_kpts

Keyword arguments

  • binary: if true write in Fortran binary format
source
WannierIO.read_spnMethod
read_spn(filename::AbstractString)

Read the spn file.

Return a n_bands * n_bands * n_kpts * 3 array.

source
WannierIO.write_spnMethod
write_spn(filename::String, S::Array{Complex,4})

Write the spn file.

The header is "Created by WannierIO.jl CURRENT_DATE".

source
WannierIO.read_unkMethod
read_unk(filename::AbstractString)

Read UNK file for the periodic part of Bloch wavefunctions.

Return

  • ik: k-point index
  • Ψ: periodic part of Bloch wavefunctions in real space
source
WannierIO.write_unkMethod
write_unk(filename::AbstractString, ik::Integer, Ψ::Array{Complex,4}; binary=false)

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)

Keyword arguments

  • binary: write as Fortran unformatted file
source
WannierIO.read_nnkpMethod
read_nnkp(filename::AbstractString)

Read the nnkp file.

Return

  • recip_lattice: each column is a reciprocal lattice vector
  • kpoints: 3 * n_kpts, in fractional coordinates
  • kpb_k: n_bvecs * n_kpts, index of kpoints
  • kpb_b: 3 * n_bvecs * n_kpts, fractional w.r.t recip_lattice
source
WannierIO.write_nnkpMethod
write_nnkp(filename::AbstractString, bvectors::BVectors, n_wann::Integer)

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

Arguments

  • recip_lattice: each column is a reciprocal lattice vector
  • kpoints: 3 * n_kpts, in fractional coordinates
  • kpb_k: n_bvecs * n_kpts, index of kpoints
  • kpb_b: 3 * n_bvecs * n_kpts, fractional w.r.t recip_lattice
  • n_wann: if given, write a auto_projections block
  • exclude_bands: if given, write the specified band indexes in the exclude_bands block
Note

Only a preliminary version, use auto_projections, no exclude_bands.

source
WannierIO.write_chkMethod
write_chk(filename, chk::Chk)

Write formatted or binary chk file.

Keyword arguments

  • binary: write as Fortran unformatted file
source
WannierIO.read_w90_tbdatMethod
read_w90_tbdat(filename::AbstractString)

Read seedname_tb.dat.

Return

  • lattice: each column is a lattice vector
  • R: $\bm{R}$ vectors on which operators are defined
  • N: degeneracies of $\bm{R}$ vectors
  • H: Hamiltonian $\bm{H}(\bm{R})$
  • r: position operator
source
WannierIO.write_HH_RMethod
write_HH_R(filename, H, R; N=nothing, header=nothing)

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

Arguments

  • filename: usually seedname_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 seedname_tb.dat file, i.e., from Fortran F12.6 to E15.8.

source
WannierIO.read_w90_bandMethod
read_w90_band(seedname::AbstractString)

Read SEEDNAME_band.dat, SEEDNAME_band.kpt, SEEDNAME_band.labelinfo.dat.

Return

  • kpoints: 3 * n_kpts, fractional coordinates
  • E: n_bands * n_kpts, band energies
  • x: n_kpts, x axis value, in cartesian length
  • symm_idx: index of high-symmetry points in seedname_band.dat
  • symm_label: name of high-symmetry points
source
WannierIO.read_w90_band_datMethod
read_w90_band_dat(filename::AbstractString)

Read seedname_band.dat file.

Return

  • x: n_kpts, x axis value, in cartesian length
  • E: n_bands * n_kpts, band energies
source
WannierIO.read_w90_band_kptMethod
read_w90_band_kpt(filename::AbstractString)

Read seedname_band.kpt file.

Return

  • kpoints: 3 * n_kpts, fractional coordinates
  • weights: n_kpts, weights of kpoints
source
WannierIO.read_w90_band_labelinfoMethod
read_w90_band_labelinfo(filename::AbstractString)

Read seedname_band.labelinfo file.

Return

  • symm_idx: index of high-symmetry points in seedname_band.dat
  • symm_label: name of high-symmetry points
source
WannierIO.write_w90_bandMethod
write_w90_band(seedname, kpoints, E, x, symm_idx, symm_label)

Write SEEDNAME_band.dat, SEEDNAME_band.kpt, SEEDNAME_band.labelinfo.dat.

Arguments

  • seedname: seedname of SEEDNAME_band.dat, SEEDNAME_band.kpt, SEEDNAME_band.labelinfo.dat
  • kpoints: 3 * n_kpts, fractional coordinates
  • E: n_bands * n_kpts, band energies
  • x: n_kpts, x axis value, in cartesian length
  • symm_idx: index of high-symmetry points in seedname_band.dat
  • symm_label: name of high-symmetry points
source
WannierIO.write_w90_band_datMethod
write_w90_band_dat(filename, x, E)

Write seedname_band.dat file.

Arguments

  • filename: filename of seedname_band.dat
  • x: n_kpts, x axis value, in cartesian length
  • E: n_bands * n_kpts, band energies
source
WannierIO.write_w90_band_kptMethod
write_w90_band_kpt(filename, kpoints, weights=nothing)

Write seedname_band.kpt file.

Arguments

  • filename: filename of seedname_band.kpt
  • kpoints: 3 * n_kpts, fractional coordinates
  • weights: n_kpts, optional, weights of kpoints
source
WannierIO.write_w90_band_labelinfoMethod
write_w90_band_labelinfo(filename, symm_idx, symm_label, x, kpoints)

Write seedname_band.labelinfo file.

Arguments

  • filename: filename of seedname_band.labelinfo
  • symm_idx: index of high-symmetry points in seedname_band.dat
  • symm_label: name of high-symmetry points
  • x: n_kpts, x axis value, in cartesian length
  • kpoints: 3 * n_kpts, fractional coordinates
source