Wannier90 files
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_amn — Functionread_amn(filename)
read_amn(filename, ::FortranText)
read_amn(filename, ::FortranBinaryStream)Read wannier90 amn file.
Return
A: length-n_kptsvector, each element is an_bands * n_wannmatrix.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.
WannierIO.write_amn — Functionwrite_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 filenameA: a length-n_kptsvector, each element is an_bands * n_wannmatrix
Keyword arguments
header: 1st line of the filebinary: write as Fortran unformatted file, which is the Wannier90 default. Here thebinarykwargs 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.
WannierIO.read_w90_band — Methodread_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 lengtheigenvalues: length-n_kptsvector, each element is a length-n_bandsvector of band energieskpoints: a vector of lengthn_kpts, fractional coordinateskweights: a vector of lengthn_kpts, weights of kpointssymm_point_indices: index of high-symmetry points inprefix_band.datsymm_point_labels: name of high-symmetry points
WannierIO.read_w90_band_dat — Methodread_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 lengtheigenvalues: length-n_kptsvector, each elemnt is a length-n_bandsvector of band energiesextras: optional (thepostw90.xmight generate a file with a third column), same size aseigenvalues, often the color of each eigenvalue, e.g., spin projection
WannierIO.read_w90_band_kpt — Methodread_w90_band_kpt(filename)
Read a prefix_band.kpt file.
Return
kpoints: a vector of lengthn_kpts, fractional coordinateskweights: a vector of lengthn_kpts, weights of kpoints
WannierIO.read_w90_band_labelinfo — Methodread_w90_band_labelinfo(filename)
Read prefix_band.labelinfo file.
Return
symm_point_indices: index of high-symmetry points inprefix_band.datsymm_point_labels: name of high-symmetry points
WannierIO.write_w90_band — Methodwrite_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 ofprefix_band.dat, prefix_band.kpt, prefix_band.labelinfo.dat
Keyword Arguments
x:n_kpts, x axis value, in cartesian lengtheigenvalues: length-n_kptsvector, each element is a length-n_bandsvector of band energieskpoints: length-n_kptsvector, fractional coordinateskweights: a vector of lengthn_kpts, weights of kpointssymm_point_indices: index of high-symmetry points inprefix_band.datsymm_point_labels: name of high-symmetry points
WannierIO.write_w90_band_dat — Methodwrite_w90_band_dat(filename; x, eigenvalues, extras)
Write prefix_band.dat file.
Arguments
filename: filename ofprefix_band.dat
Keyword Arguments
x:n_kpts, x axis value, in cartesian lengtheigenvalues: length-n_kptsvector, each element is a length-n_bandsvector of band energiesextras: optional, same size aseigenvalues, will be written as the third column ofprefix_band.dat. Theprefix-band.datfile generated bypostw90.xsometimes has a third column for e.g. the color of the eigenvalues
WannierIO.write_w90_band_kpt — Methodwrite_w90_band_kpt(filename; kpoints, kweights)
Write prefix_band.kpt file.
Arguments
filename: filename ofprefix_band.kpt
Keyword Arguments
kpoints: length-n_kptsvector, fractional coordinateskweights:n_kpts, optional, weights of kpoints, default to 1.0.
WannierIO.write_w90_band_labelinfo — Methodwrite_w90_band_labelinfo(
filename;
x,
kpoints,
symm_point_indices,
symm_point_labels
)
Write prefix_band.labelinfo file.
Arguments
filename: filename ofprefix_band.labelinfo
Keyword Arguments
x:n_kpts-vector, x axis value, in cartesian lengthkpoints: length-n_kptsvector, fractional coordinatessymm_point_indices: index of high-symmetry points inprefix_band.datsymm_point_labels: name of high-symmetry points
WannierIO.get_U — Methodget_U(chk)
Extract the combined U = Udis * Uml matrices from Chk.
WannierIO.get_Udis — Methodget_Udis(chk)
Extract disentanglement Udis matrices from Chk.
WannierIO.read_chk — Methodread_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.
WannierIO.write_chk — Functionwrite_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 isfalsesince Fortran binary depends on compiler and platform, so it is not guaranteed to always work.
WannierIO.read_eig — Functionread_eig(filename)
read_eig(filename, ::FortranText)
read_eig(filename, ::FortranBinaryStream)Read the wannier90 eig file.
Return
eigenvalues: a lenth-n_kptsvector, each element is a length-n_bandsvector
The 1st version is a convenience wrapper function for the 2nd and 3rd versions.
WannierIO.write_eig — Functionwrite_eig(filename, eigenvalues; binary=false)
write_eig(filename, eigenvalues, ::FortranText)
write_eig(filename, eigenvalues, ::FortranBinaryStream)Write eig file.
Arguments
eigenvalues: a length-n_kptsvector, each element is a length-n_bandsvector
Keyword arguments
binary: if true write in Fortran binary format.
WannierIO.read_w90_hrdat — Methodread_w90_hrdat(filename)
Read prefix_hr.dat.
Return
Rvectors: $\mathbf{R}$-vectors on which operators are definedRdegens: degeneracies of each $\mathbf{R}$-vectorH: Hamiltonian $\mathbf{H}(\mathbf{R})$header: the first line of the file
WannierIO.write_w90_hrdat — Methodwrite_w90_hrdat(filename; Rvectors, Rdegens, H, header)
Write prefix_hr.dat.
Keyword arguments
See the return values of read_w90_hrdat.
WannierIO.read_mmn — Functionread_mmn(filename)
read_mmn(filename, ::FortranText)
read_mmn(filename, ::FortranBinaryStream)Read wannier90 mmn file.
Return
M: length-n_kptsvector, each element is a length-n_bvecsvector, then each element is an_bands * n_bandsmatrixkpb_k: length-n_kptsvector, each element is a length-n_bvecsvector of integers for the indices of the neighboring kpointskpb_G: length-n_kptsvector, each element is a lenght-n_bvecsvector of ofVec3{Int}, which are the translation vectorsheader: 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.
WannierIO.write_mmn — Functionwrite_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 nameM: length-n_kptsvector ofn_bands * n_bands * n_bvecsarrayskpb_k: length-n_kptsvector of length-n_bvecsvector of integerskpb_G: length-n_kptsvector of length-n_bvecsvector ofVec3{Int}for bvectors
Keyword arguments
header: header stringbinary: if true write in Fortran binary format
The 1st version is a convenience wrapper for the other two.
WannierIO.read_nnkp — Functionread_nnkp(filename)
read_nnkp(filename, ::Wannier90Text)
read_nnkp(filename, ::Wannier90Toml)Read wannier90 nnkp file.
Return
lattice: each column is a lattice vectorrecip_lattice: each column is a reciprocal lattice vectorkpoints: length-n_kptsvector, each element isVec3, in fractional coordinatesprojections: length-n_projsvector ofHydrogenOrbitalauto_projections: optional, the number of Wannier functionsn_wannfor automatic initial projectionskpb_k: length-n_kptsvector, each element is a length-n_bvecsvector of integers, index of kpointskpb_G: length-n_kptsvector, each element is a length-n_bvecsvector, then each element isVec3for translations, fractional w.r.trecip_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.
WannierIO.write_nnkp — Functionwrite_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 formatrecip_lattice: each column is a reciprocal lattice vectorkpoints: length-n_kptsvector ofVec3, in fractional coordinateskpb_k: length-n_kptsvector, each element is a length-n_bvecsvector of integers, index of kpointskpb_G: length-n_kptsvector, each element is a length-n_bvecsvector, then each element is aVec3for translation vector, fractional w.r.t.recip_latticeprojections: optional, length-n_projsvector ofHydrogenOrbitalauto_projections: optional, the number of Wannier functionsn_wannfor automatic initial projections. If given, write anauto_projectionsblockexclude_bands: if given, write the specified band indices in theexclude_bandsblockheader: first line of the file
WannierIO.read_w90_rdat — Methodread_w90_rdat(filename)
Read prefix_r.dat.
Return
Rvectors: $\mathbf{R}$-vectors on which operators are definedr_x: $x$-component of position operatorr_y: $y$-component of position operatorr_z: $z$-component of position operatorheader: the first line of the file
WannierIO.write_w90_rdat — Methodwrite_w90_rdat(filename; Rvectors, r_x, r_y, r_z, header)
Write prefix_r.dat.
Keyword arguments
See the return values of read_w90_rdat.
WannierIO.read_spn — Functionread_spn(filename)
read_spn(filename, ::FortranText)
read_spn(filename, ::FortranBinary)Read the wannier90 spn file.
Return
Sx: spin x, a length-n_kptsvector, each element is an_bands-by-n_bandsmatrixSy: spin y, a length-n_kptsvector, each element is an_bands-by-n_bandsmatrixSz: spin z, a length-n_kptsvector, each element is an_bands-by-n_bandsmatrixheader: 1st line of the file
WannierIO.write_spn — Functionwrite_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.
WannierIO.read_w90_tbdat — Methodread_w90_tbdat(filename)
Read prefix_tb.dat.
Return
lattice: each column is a lattice vector in ÅRvectors: $\mathbf{R}$-vectors on which operators are definedRdegens: degeneracies of each $\mathbf{R}$-vectorH: Hamiltonian $\mathbf{H}(\mathbf{R})$r_x: $x$-component of position operatorr_x: $y$-component of position operatorr_x: $z$-component of position operatorheader: the first line of the file
WannierIO.write_w90_tbdat — Methodwrite_w90_tbdat(
filename;
lattice,
Rvectors,
Rdegens,
H,
r_x,
r_y,
r_z,
header
)
Write prefix_tb.dat.
Keyword arguments
See the return values of read_w90_tbdat.
WannierIO.read_uHu — Functionread_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-generateduHufile, 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_kptsvector, each element is an_bvecs * n_bvecsmatrix, then each element is an_bands * n_bandsmatrixheader: 1st line of the file
WannierIO.write_uHu — Functionwrite_uHu(filename, uHu; binary=false, header)
write_uHu(filename, uHu, ::FortranText; header)
write_uHu(filename, uHu, ::FortranBinary; header)Write the uHu file.
Keyword Arguments
transpose_band_indices: seeread_uHu
WannierIO.read_uIu — Functionread_uIu(filename)
read_uIu(filename, ::FortranText; transpose_band_indices=true)
read_uIu(filename, ::FortranBinary; transpose_band_indices=true)Read the wannier90 uIu file.
Keyword Arguments
transpose_band_indices: QE pw2wannier90.x writes the matrix in a strange transposed manner; if reading a QE-generateduIufile, this flag should be true to restore the band indices order, so that the returned matrix has the correct order, i.e.,uIu[ik][ib1, ib2][m, n]is $\langle u_{m, k + b_1} | u_{n, k + b_2} \rangle$
Return
uIu: a length-n_kptsvector, each element is an_bvecs * n_bvecsmatrix, then each element is an_bands * n_bandsmatrixheader: 1st line of the file
WannierIO.write_uIu — Functionwrite_uIu(filename, uIu; binary=false, header)
write_uIu(filename, uIu, ::FortranText; header)
write_uIu(filename, uIu, ::FortranBinary; header)Write the uIu file.
Keyword Arguments
transpose_band_indices: seeread_uIu
WannierIO.read_unk — Functionread_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)
WannierIO.write_unk — Functionwrite_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
WannierIO.read_win — Functionread_win(filename; standardize=true)
read_win(filename, ::Wannier90Text; standardize=true)
read_win(filename, ::Wannier90Toml; standardize=true)Read wannier90 input win file.
Arguments
filename: The name of the input file.
Keyword Arguments
standardize: sanity check and fix the input parameters, e.g., setnum_bands = num_wannifnum_bandsis not specified, convertatoms_cartalways toatoms_frac, etc. See alsostandardize_win!.
WannierIO.write_win — Functionwrite_win(filename, params; header)
write_win(filename, params, ::Wannier90Text; header)
write_win(filename, params, ::Wannier90Toml; header)
write_win(filename; header, params...)
write_win(filename, ::Wannier90Text; header, params...)
write_win(filename, ::Wannier90Toml; header, params...)Write input parameters into a wannier90 win file.
There are two choice for passing the input parameters:
- as a
Dict(orOrderedDictto preserve ordering) to theparamsargument - as keyword arguments
params..., with argument names the same as the input parameters of wannier90
Examples
using WannierIO
# you can also use `Dict` or `OrderedDict`
params = (;
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,
)
write_win("silicon.win", params)WannierIO.read_wout — Methodread_wout(filename; iterations)
Parse wannier90 wout file.
Keyword Arguments
iterations: iftrue, parse all the iterations of disentanglement and maximal localization. Default isfalse.
Return
lattice: each column is a lattice vector in Årecip_lattice: each column is a reciprocal lattice vector in Å⁻¹atom_labels: atomic symbolsatom_positions: in fractional coordinateskgrid: kpoint grid used in Wannierizationcenters: center of each final WF, in Åspreads: spread of each final WF, in Ųsum_centers: sum of final WF centers, in Åsum_spreads: sum of final WF spreads, in ŲΩI,ΩD,ΩOD,Ωtotal: final spread (components) in Ųphase_factors: optional, global (multiplicative) phase factor for obtaining real-valued (or close to real) MLWFsim_re_ratios: optional, maximum Im/Re ratioiterations: disentanglement and max localization convergence history, only parsed when kwargiterations=true
WannierIO.read_w90_wsvec — Methodread_w90_wsvec(filename)
Read prefix_wsvec.dat.
Return
mdrs: whether use MDRS interpolation, i.e. theuse_ws_distancein the headerRvectors: the $\mathbf{R}$-vectorsTvectors: the $\mathbf{T}_{m n \mathbf{R}}$-vectors. Returned onlymdrs = true.Tdegens: the degeneracies of $\mathbf{T}_{m n \mathbf{R}}$-vectors. Returned onlymdrs = true.n_wann: number of WFsheader: the first line of the file
WannierIO.write_w90_wsvec — Methodwrite_w90_wsvec(
filename;
Rvectors,
n_wann,
Tvectors,
Tdegens,
header
)
Write prefix_wsvec.dat.
Keyword Arguments
n_wann: for Wigner-Seitz Rvectors, needs to provide an_wannfor number of Wannier functions; for MDRS Rvectors, then_wannis optional and can be automatically determined from theTvectorsTvectorsandTdegens: if provided, write in MDRS format; otherwise, write in Wigner-Seitz format
Also see the return values of read_w90_wsvec.
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.default_band_kpt_kweights — Methoddefault_band_kpt_kweights(kpoints)
Wannier90 default kweights in prefix_band.kpt is all 1.0.
WannierIO.Chk — TypeStruct 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 timen_bands: number of bands, can be auto set in constructor according to dimensions of other variablesn_exclude_bands: number of excluded bands, can be auto set in constructorexclude_bands: Indices of excluded bands, starts from 1. Vector of integers, size:n_exclude_bandslattice: Matrix of size 3 x 3, each column is a lattice vector in Å unitrecip_lattice: Matrix of size 3 x 3, each column is a reciprocal lattice vector in Å⁻¹ unitn_kpts: number of kpoints, can be auto set in constructorkgrid: dimensions of kpoint grid, 3 integerskpoints: kpoint coordinates, fractional, length-n_kptsvectorn_bvecs: number of b-vectors, can be auto set in constructorn_wann: number of Wannier functions, can be auto set in constructorcheckpoint: a string to indicate the current step (after disentanglement, after maximal localization, ...) in wannier90have_disentangled: Have finished disentanglement or notΩI: Omega invariant part of MV spreads, in Ų unitdis_bands: Indices of bands taking part in disentanglement, not frozen bands! length-n_kptsvector, each element is a length-n_bandsvector of bool.This is needed since W90 puts all the disentanglement bands in the first several rows of
Udis, (and the first few columns ofUdisare 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 fromdis_bandsUdis: Semi-unitary matrix for disentanglement, length-n_kptsvector, each elment has size:n_bandsxn_wann, i.e., theu_matrix_optin wannier90Uml: Unitary matrix for maximal localization, length-n_kptsvector, each element has size:n_wannxn_wann, i.e., theu_matrixin wannier90. The abbreviationmlstands for maximal localization, so as to differentiate from the (combined) unitary matrixU = Udis * Uml.M: Wannier-gauge overlap matrix, length-n_kptsvector of length-n_bvecsvector, each element is a matrix of sizen_wannxn_wann, i.e., them_matrixin wannier90r: Wannier function centers, length-n_wannvector, Cartesian coordinates in Å unit, i.e., thewannier_centresvariable in wannier90ω: Wannier function spreads, length-n_wannvector, Ų unit, i.e., thewannier_spreadsvariable in wannier90
WannierIO.Chk — MethodChk(
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.
Base.isapprox — Methodisapprox(a, b)
Compare two Chk objects.
WannierIO._check_eig_order — Method_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.
WannierIO._reshape_eig — Method_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.
WannierIO.write_HH_R — Methodwrite_HH_R(filename, H, R; N, header)
Write the real space Hamiltonian to a prefix_HH_R.dat file.
Arguments
filename: usuallyprefix_HH_R.datH: an_wann * n_wann * n_rvecsarray of HamiltonianR: an_rvecs * 3array of integers
Keyword arguments
N: an_rvecsvector of integers, the degeneracy of each R vectorheader: a string, the header of the file
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.
WannierIO._check_dimensions_M_kpb — Method_check_dimensions_M_kpb(M, kpb_k, kpb_G)
WannierIO._check_dimensions_kpb — MethodCheck the dimensions between the quantities are consistent.
WannierIO.HydrogenOrbital — TypeHydrogen-like analytic orbitals.
Follows the definition in Wannier90, see https://wannier90.readthedocs.io/en/latest/user_guide/wannier90/postproc/#projections-block https://wannier90.readthedocs.io/en/latest/user_guide/wannier90/projections/
struct HydrogenOrbital <: WannierIO.OrbitalFields
center: 3 real numbers of the projection center, in fractional coordinatesn: positive integer, principle quantum number $n > 0$ for the radial functionl: non-negative integer, angular momentum $l \ge 0$ of real spherical harmonics $Y_{lm}(\theta, \phi)$m: integer, magnetic quantum number $m$, $-l \leq m \leq l$α: positive real number, controlling the spread of the radial functionzaxis: 3 real numbers, the z-axis from which the polar angle $\theta$ is measured, default is[0, 0, 1]xaxis: 3 real numbers, the x-axis from which the azimuthal angle $\phi$ is measured, must be orthogonal tozaxis, default is[1, 0, 0]
WannierIO._check_dimensions_Sx_Sy_Sz — Method_check_dimensions_Sx_Sy_Sz(Sx, Sy, Sz)
WannierIO.read_u_mat — Methodread_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) orU(for maximal localization) matriceskpoints: fractional kpoint coordinatesheader: 1st line of the file
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.
WannierIO.write_u_mat — Methodwrite_u_mat(filename, U, kpoints; header=default_header())Write wannier90 prefix_u.mat or prefix_u_dis.mat file.
Arguments
filename: the input file nameU:Udis(for disentanglement) orU(for maximal localization) matriceskpoints: fractional kpoint coordinates
Keyword arguments
header: 1st line of the file, optional
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.
WannierIO._check_win_required_params — Method_check_win_required_params(kwargs)
WannierIO.standardize_win! — Methodstandardize_win!(params)
Sanity check and add missing input parameters from a win file.
See also read_win.
WannierIO._parse_wout_atoms — MethodParse block
| Site Fractional Coordinate Cartesian Coordinate (Ang) |
+----------------------------------------------------------------------------+
| Si 1 0.00000 0.00000 0.00000 | 0.00000 0.00000 0.00000 |
| Si 2 0.25000 0.25000 0.25000 | 1.35763 1.35763 1.35763 |
*----------------------------------------------------------------------------*WannierIO._parse_wout_disentangle — MethodParse block
Extraction of optimally-connected subspace
------------------------------------------
+---------------------------------------------------------------------+<-- DIS
| Iter Omega_I(i-1) Omega_I(i) Delta (frac.) Time |<-- DIS
+---------------------------------------------------------------------+<-- DIS
1 25.38943399 21.32896063 1.904E-01 0.00 <-- DIS
2 21.53095611 20.16097533 6.795E-02 0.01 <-- DIS
3 20.40788223 19.35260423 5.453E-02 0.01 <-- DIS
4 19.53883989 18.75563591 4.176E-02 0.02 <-- DIS
...
341 16.22884440 16.22884440 -1.883E-10 2.43 <-- DIS
342 16.22884440 16.22884440 -1.799E-10 2.44 <-- DIS
<<< Delta < 2.000E-10 over 3 iterations >>>
<<< Disentanglement convergence criteria satisfied >>>
Final Omega_I 16.22884440 (Ang^2)
+----------------------------------------------------------------------------+
Time to disentangle bands 2.546 (sec)WannierIO._parse_wout_lattice — MethodParse block
Lattice Vectors (Ang)
a_1 0.000000 2.715265 2.715265
a_2 2.715265 0.000000 2.715265
a_3 2.715265 2.715265 0.000000WannierIO._parse_wout_recip_lattice — MethodParse block
Reciprocal-Space Vectors (Ang^-1)
b_1 -1.157011 1.157011 1.157011
b_2 1.157011 -1.157011 1.157011
b_3 1.157011 1.157011 -1.157011WannierIO._parse_wout_wannierize — MethodParse block
*------------------------------- WANNIERISE ---------------------------------*
+--------------------------------------------------------------------+<-- CONV
| Iter Delta Spread RMS Gradient Spread (Ang^2) Time |<-- CONV
+--------------------------------------------------------------------+<-- CONV
------------------------------------------------------------------------------
Initial State
WF centre and spread 1 ( -0.000005, 0.000021, 0.000023 ) 2.56218734
WF centre and spread 2 ( 0.000013, -0.000054, 0.000016 ) 3.19493515
WF centre and spread 3 ( -0.000005, -0.000054, -0.000055 ) 3.19482997
WF centre and spread 4 ( 0.000012, 0.000015, -0.000058 ) 3.19526437
WF centre and spread 5 ( 1.357637, 1.357611, 1.357610 ) 2.56218214
WF centre and spread 6 ( 1.357619, 1.357684, 1.357617 ) 3.19532825
WF centre and spread 7 ( 1.357638, 1.357687, 1.357686 ) 3.19513205
WF centre and spread 8 ( 1.357620, 1.357617, 1.357694 ) 3.19460833
Sum of centres and spreads ( 5.430528, 5.430529, 5.430534 ) 24.29446759
0 0.243E+02 0.0000000000 24.2944680346 2.48 <-- CONV
O_D= 0.2135529 O_OD= 7.8520707 O_TOT= 24.2944680 <-- SPRD
------------------------------------------------------------------------------
Cycle: 1
WF centre and spread 1 ( -0.000005, 0.000020, 0.000022 ) 2.46316318
WF centre and spread 2 ( 0.000014, -0.000057, 0.000015 ) 3.19187236
WF centre and spread 3 ( -0.000005, -0.000057, -0.000058 ) 3.19179103
WF centre and spread 4 ( 0.000012, 0.000014, -0.000061 ) 3.19222621
WF centre and spread 5 ( 1.357637, 1.357612, 1.357611 ) 2.46315800
WF centre and spread 6 ( 1.357618, 1.357687, 1.357618 ) 3.19226713
WF centre and spread 7 ( 1.357637, 1.357690, 1.357689 ) 3.19209166
WF centre and spread 8 ( 1.357619, 1.357619, 1.357697 ) 3.19156919
Sum of centres and spreads ( 5.430528, 5.430529, 5.430534 ) 24.07813875
1 -0.216E+00 0.2558717278 24.0781391952 2.49 <-- CONV
O_D= 0.2218113 O_OD= 7.6274835 O_TOT= 24.0781392 <-- SPRD
Delta: O_D= 0.8258326E-02 O_OD= -0.2245872E+00 O_TOT= -0.2163288E+00 <-- DLTA
------------------------------------------------------------------------------
Cycle: 2
...
------------------------------------------------------------------------------
Cycle: 45
WF centre and spread 1 ( 0.000001, 0.000006, 0.000006 ) 1.95373328
WF centre and spread 2 ( 0.000016, -0.000065, 0.000019 ) 3.27910139
WF centre and spread 3 ( -0.000008, -0.000065, -0.000066 ) 3.27921479
WF centre and spread 4 ( 0.000014, 0.000016, -0.000070 ) 3.27965818
WF centre and spread 5 ( 1.357631, 1.357627, 1.357627 ) 1.95372427
WF centre and spread 6 ( 1.357616, 1.357695, 1.357615 ) 3.27949625
WF centre and spread 7 ( 1.357641, 1.357699, 1.357697 ) 3.27951005
WF centre and spread 8 ( 1.357617, 1.357616, 1.357707 ) 3.27899768
Sum of centres and spreads ( 5.430528, 5.430529, 5.430534 ) 23.58343588
45 -0.186E-09 0.0000077396 23.5834363262 2.88 <-- CONV
O_D= 0.2612087 O_OD= 7.0933833 O_TOT= 23.5834363 <-- SPRD
Delta: O_D= 0.2557594E-06 O_OD= -0.2559458E-06 O_TOT= -0.1863789E-09 <-- DLTA
------------------------------------------------------------------------------
<<< Delta < 2.000E-10 over 3 iterations >>>
<<< Wannierisation convergence criteria satisfied >>>WannierIO._parse_wout_wf_center_spread — MethodParse block
WF centre and spread 1 ( -0.000005, 0.000021, 0.000023 ) 2.56218734
WF centre and spread 2 ( 0.000013, -0.000054, 0.000016 ) 3.19493515
WF centre and spread 3 ( -0.000005, -0.000054, -0.000055 ) 3.19482997
WF centre and spread 4 ( 0.000012, 0.000015, -0.000058 ) 3.19526437
WF centre and spread 5 ( 1.357637, 1.357611, 1.357610 ) 2.56218214
WF centre and spread 6 ( 1.357619, 1.357684, 1.357617 ) 3.19532825
WF centre and spread 7 ( 1.357638, 1.357687, 1.357686 ) 3.19513205
WF centre and spread 8 ( 1.357620, 1.357617, 1.357694 ) 3.19460833
Sum of centres and spreads ( 5.430528, 5.430529, 5.430534 ) 24.29446759