Interpolation

This page lists struct and functions for Wannier interpolation.

There are two flavors of interpolation:

  • WS: the Wigner-Seitz interpolation
  • MDRS: the Minimal-distance replica selection method

MDRS has better interpolation quality[WAN20], so it should be used in most cases. Moreover, we implement two versions of MDRS:

  • MDRSv1: the original version in Wannier90[WAN20]. This is slower than the WS interpolation.
  • MDRSv2: a slightly modified version which removes some for loops, same as what is implemented in Wannier Berri[WBERRI]. The speed is similar to the WS interpolation.

By default, we use the MDRSv2 interpolation.

Contents

Index

R vectors

RVectors

Note

The R vectors are automatically generated by calling either of these two approaches:

  • read_w90_interp to read the Wannier90 output seedname.chk.fmt file, and generate an InterpModel,
  • InterpModel a constructor to generate the InterpModel from a Wannierization Model.

So usually the user does not need to call the following functions.

Wannier.RVectorsType
struct RVectors

The R vectors for interpolation.

Fields

  • lattice: columns are lattice vectors
  • grid: number of FFT grid points in each direction, actually equal to kgrid
  • R: 3 * n_rvecs, R vectors in fractional coordinates w.r.t. lattice
  • N: n_rvecs, degeneracy of each R vector
Note

The R vectors are sorted in the same order as Wannier90.

source
Wannier.RVectorsMethod
RVectors(lattice, grid, R, N)

Constructor for RVectors.

Auto transform lattice and grid to Mat3 and Vec3, respectively.

Arguments

  • lattice: columns are lattice vectors
  • grid: number of FFT grid points in each direction, actually equal to kgrid
  • R: 3 * n_rvecs, R vectors in fractional coordinates w.r.t. lattice
  • N: n_rvecs, degeneracy of each R vector
source
Wannier.RVectorsMDRSType
struct RVectorsMDRS

The R vectors for MDRS interpolation.

Fields

  • Rvectors: Rvectors for Wigner-Seitz interpolation

For MDRSv1

  • T: n_wann * n_wann * n_rvecs, translation vectors w.r.t to lattice for MDRSv1
  • Nᵀ: n_wann * n_wann * n_rvecs, degeneracy of each T vector for MDRSv1

For MDRSv2

  • R̃vectors: RVectors containing expanded set of R + T vectors used in MDRSv2
  • R̃_RT: mapping of R̃vectors to R and T vectors
Note

The R vectors are sorted in the same order as Wannier90.

source
Wannier.RVectorsMDRSMethod
RVectorsMDRS(Rvectors, T, Nᵀ)

A friendly constructor for RVectorsMDRS.

The remaining fields R̃vectors and R̃_RT are only used in MDRSv2, the are calculated automatically based on the input arguments.

Arguments

  • Rvectors: RVectors for Wigner-Seitz interpolation
  • T: n_wann * n_wann * n_rvecs, translation vectors w.r.t to lattice for MDRSv1
  • Nᵀ: n_wann * n_wann * n_rvecs, degeneracy of each T vector for MDRSv1
source
Wannier.get_Rvectors_mdrsMethod
get_Rvectors_mdrs(lattice, rgrid, centers; atol=1e-5, max_cell=3)

Generate R vectors for MDRS interpolation (both v1 and v2).

Arguments

  • lattice: columns are lattice vectors
  • rgrid: number of FFT grid points in each direction, actually equal to kgrid
  • centers: 3 * n_wann, WF centers in fractional coordinates

Keyword arguments

  • atol: toerance for checking degeneracy, equivalent to Wannier90 input parameter ws_distance_tol
  • max_cell: number of neighboring cells to be searched, equivalent to Wannier90 input parameter ws_search_size
source
Wannier.get_Rvectors_wsMethod
get_Rvectors_ws(lattice, rgrid; atol=1e-5, max_cell=3)

Generate R vectors for Wigner-Seitz interpolation.

Arguments

  • lattice: columns are lattice vectors
  • rgrid: number of FFT grid points in each direction, actually equal to kgrid

Keyword arguments

  • atol: toerance for checking degeneracy, equivalent to Wannier90 input parameter ws_distance_tol
  • max_cell: number of neighboring cells to be searched, equivalent to Wannier90 input parameter ws_search_size
source

KRVectors

The KRVectors contains both kpoint mappings and RVectors, allowing both forward and inverse Fourier transform.

Wannier.KRVectorsType
struct KRVectors

Contains both kpoints mapping and R vectors.

Fields

  • lattice: each column is a lattice vector
  • kgrid: number of kpoints along 3 directions
  • kpoints: each column is a kpoint in fractional coordinates
  • k_xyz: kpoints mappings from ik to [ikx, iky, ikz]
  • xyz_k: kpoints mappings from [ikx, iky, ikz] to ik
  • Rvectors: R vectors, can be either RVectors or RVectorsMDRS
  • recip_lattice: each column is a reciprocal lattice vector
  • n_kpts: number of kpoints
  • n_rvecs: number of R vectors
  • n_r̃vecs: number of R̃ vectors, only for MDRS
source
Wannier.KRVectorsMethod
KRVectors(lattice, kgrid, kpoints, Rvectors)

A friendly constructor for KRVectors.

The kpoint mappings k_xyz and xyz_k are generated based on kpoints and kgrid. Remaining fields are generated automatically based on input arguments.

Arguments

  • lattice: each column is a lattice vector
  • kgrid: number of kpoints along 3 directions
  • kpoints: each column is a kpoint in fractional coordinates
  • Rvectors: R vectors, can be either RVectors or RVectorsMDRS
source
Wannier.KRVectorsMethod
KRVectors(lattice, kgrid, kpoints, k_xyz, xyz_k, Rvectors)

A friendly constructor for KRVectors.

Remaining fields are generated automatically based on input arguments.

Arguments

  • lattice: each column is a lattice vector
  • kgrid: number of kpoints along 3 directions
  • kpoints: each column is a kpoint in fractional coordinates
  • k_xyz: kpoints mappings from ik to [ikx, iky, ikz]
  • xyz_k: kpoints mappings from [ikx, iky, ikz] to ik
  • Rvectors: R vectors, can be either RVectors or RVectorsMDRS
source
Wannier.KRVectorsMethod
KRVectors(Rvectors)

Construct a KRVectors from RVectors or RVectorsMDRS.

The kpoints part are just empty.

Note

If we only need to Wannier interpolate operators, e.g., interpolate band structure from tb.dat file, then we don't need info on the kpoint grid of Wannierization, so we can just leave it empty.

source

Fourier transforms

Performs Fourier transform between k and R spaces defined by KRVectors.

Note

In general, the user only need to call the fourier(kRvectors, Oᵏ) and invfourier(kRvectors, Oᴿ, kpoints) functions directly, irrespective of the type of kRVectors. Other functions are internally used based on the type of KRVectors, thanks to multiple dispatch.

Wannier._fourier_mdrs_v1Method
_fourier_mdrs_v1(kRvectors::KRVectors{T,RVectorsMDRS{T}}, Oᵏ::Array{Complex{T},3}) where {T<:Real}

Fourier transform operator from k space to R space for MDRSv1.

The same as the forward Fourier transform of WS interpolation,

\[X_{mn}(\bm{R}) = \frac{1}{N_{\bm{k}}} \sum_{\bm{k}} \exp(-i \bm{k} \bm{R}) X_{mn}(\bm{k}),\]

where N_{\bm{k}} is the total number of kpoints.

See also _invfourier_mdrs_v1.

Arguments

  • kRvectors: for MDRS interpolation
  • Oᵏ: n_wann * n_wann * n_kpts, operator matrix in k space
source
Wannier._fourier_mdrs_v2Method
_fourier_mdrs_v2(kRvectors::KRVectors{T,RVectorsMDRS{T}}, Oᵏ::Array{Complex{T},3})
where {T<:Real}

Fourier transform operator from k space to R space for MDRSv2.

To differentiate from the ${X}_{mn}(\bm{R})$ of MDRSv1, I use $\widetilde{X}_{mn}(\widetilde{ \bm{R} })$ for MDRSv2,

\[\widetilde{X}_{mn}(\widetilde{ \bm{R} }) = \sum_{ \bm{R} } \frac{1}{ N_{ \bm{R} } \mathcal{N}_{ mn \bm{R} } } X_{mn}(\bm{R}) \sum_{j=1}^{\mathcal{N}_{mn \bm{R} }} \delta_{ \widetilde{ \bm{R} }, \bm{R} + \bm{T}_{ mn \bm{R} }^{(j)} },\]

where $N_{\bm{R}}$ is the degeneracy of R vectors (not the total number of R vectors), and $\mathcal{N}_{ mn \bm{R} }$ is the degeneracy of $\bm{T}_{ mn \bm{R} }$ vectors.

See also _invfourier_mdrs_v2.

Arguments

  • kRvectors: for MDRS interpolation
  • Oᵏ: n_wann * n_wann * n_kpts, operator matrix in k space
Note

This is a faster version of MDRS, by removing some for loops. But need to expand R vectors to R + T, so different from the R vectors in Wannier90 output seeedname_tb.dat.

source
Wannier._invfourier_mdrs_v1!Method
_invfourier_mdrs_v1!(Oᵏ,
    Rvectors::RVectorsMDRS{FT}, Oᴿ::Array{Complex{FT},3}, kpoints::AbstractMatrix{FT}
) where {FT<:Real}

In-place version of _invfourier_mdrs_v1.

Arguments

  • Oᵏ: n_wann * n_wann * n_kpts, operator matrix in k space
source
Wannier._invfourier_mdrs_v1Method
_invfourier_mdrs_v1(
    Rvectors::RVectorsMDRS{FT}, Oᴿ::Array{Complex{FT},3}, kpoints::AbstractMatrix{FT}
) where {FT<:Real}

Inverse Fourier transform operator from R space to k space for MDRSv1 interpolation.

\[X_{mn}(\bm{k}) = \sum_{\bm{R}} \frac{1}{ \mathcal{N}_{mn \bm{R}} } X_{mn}(\bm{R}) \sum_{j=1}^{\mathcal{N}_{ mn \bm{R} }} \exp\left( i \bm{k} \cdot \left( \bm{R} + \bm{T}_{ mn \bm{R} }^{(j)} \right) \right)\]

where $N_{\bm{R}}$ is the degeneracy of R vectors (not the total number of R vectors), and $\mathcal{N}_{ mn \bm{R} }$ is the degeneracy of $\bm{T}_{ mn \bm{R} }$ vectors.

See also _fourier_mdrs_v1.

Arguments

  • Oᴿ: n_wann * n_wann * n_rvecs, real space operator matrix
  • kpoints: 3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
Note

This is the MDRSv1 using plain loops. It is slower, but reproduce the Wannier90 behavior, and generate the same seedname_tb.dat file as Wannier90. The MDRSv2 has expanded set of R vectors, thus cannot be compared with tb.dat file.

source
Wannier._invfourier_mdrs_v2!Method
_invfourier_mdrs_v2!(Oᵏ,
    Rvectors::RVectorsMDRS{T}, Oᴿ̃::Array{Complex{T},3}, kpoints::AbstractMatrix{T}
) where {T<:Real}

In-place version of _invfourier_mdrs_v2.

Arguments

  • Oᵏ: n_wann * n_wann * n_kpts, operator matrix in k space
source
Wannier._invfourier_mdrs_v2Method
_invfourier_mdrs_v2(
    Rvectors::RVectorsMDRS{T}, Oᴿ̃::Array{Complex{T},3}, kpoints::AbstractMatrix{T}
) where {T<:Real}

Inverse Fourier transform operator from R space to k space for MDRSv2 interpolation.

\[X_{mn}(\bm{k}) = \sum_{ \widetilde{\bm{R}} } \exp\left( i \bm{k} \widetilde{ \bm{R} } \right) \widetilde{X}_{mn}(\widetilde{\bm{R}})\]

See also _fourier_mdrs_v2.

Arguments

  • Oᴿ̃: n_wann * n_wann * n_rvecs, real space operator matrix
  • kpoints: 3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
source
Wannier.fourierMethod
fourier(kRvectors::KRVectors{T,RVectorsMDRS{T}}, Oᵏ::Array{Complex{T},3}; version::Symbol=:v2)
where {T<:Real}

Wrapper function for Fourier transform for both MDRS v1 and v2.

Arguments

  • kRvectors: for WS interpolation
  • Oᵏ: n_wann * n_wann * n_kpts, operator matrix in k space

Keyword Arguments

  • version: :v1 or :v2, default is :v2

See also _fourier_mdrs_v1 and _fourier_mdrs_v2.

source
Wannier.fourierMethod
fourier(kRvectors::KRVectors{T,RVectors{T}}, Oᵏ::Array{Complex{T},3})
where {T<:Real}

Fourier transform operator from k space to R space for Wiger-Seitz interpolation.

\[X_{mn}(\bm{R}) = \frac{1}{N_{\bm{k}}} \sum_{\bm{k}} \exp(-i {\bm{k}} \bm{R}) X_{mn}(\bm{k}),\]

where N_{\bm{k}} is the total number of kpoints.

Arguments

  • kRvectors: for WS interpolation
  • Oᵏ: n_wann * n_wann * n_kpts, operator matrix in k space
source
Wannier.invfourier!Method
invfourier!(Oᵏ, Rvectors, Oᴿ, kpoints)

In-place version of invfourier, no allocation inside.

Arguments

  • Oᵏ: n_wann * n_wann * n_kpts, output operator matrix in k space
source
Wannier.invfourierMethod
invfourier(
    kRvectors::KRVectors{T,RVectorsMDRS{T}},
    Oᴿ::AbstractArray{Complex{T},3},
    kpoints::AbstractMatrix{T};
    version::Symbol=:v2,
) where {T<:Real}

A simple wrapper of invfourier to support kRvectors.

See also invfourier.

source
Wannier.invfourierMethod
invfourier(
    kRvectors::KRVectors{T,RVectors{T}}, Oᴿ::Array{Complex{T},3}, kpoints::AbstractMatrix{T}
) where {T<:Real}

A simple wrapper of invfourier to support kRvectors.

See also invfourier.

source
Wannier.invfourierMethod
invfourier(
    Rvectors::RVectorsMDRS{T},
    Oᴿ::AbstractArray{Complex{T},3},
    kpoints::AbstractMatrix{T};
    version::Symbol=:v2,
) where {T<:Real}

Wrapper function for inverse Fourier transform for both MDRS v1 and v2.

Arguments

  • Rvectors: for MDRS interpolation
  • Oᴿ: n_wann * n_wann * n_rvecs, operator matrix in R space
  • kpoints: 3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform

Keyword Arguments

  • version: :v1 or :v2, default is :v2

See also _invfourier_mdrs_v1 and _invfourier_mdrs_v2.

source
Wannier.invfourierMethod
invfourier(Rvectors::RVectors{T}, Oᴿ::Array{Complex{T},3}, kpoints::AbstractMatrix{T}) where {T<:Real}

Inverse Fourier transform operator from $\bm{R}$ space to $k$ space for Wigner-Seitz interpolation.

\[X_{mn}(\bm{k}) = \sum_{\bm{R}} \frac{1}{N_{\bm{R}}} \exp(i \bm{k} \bm{R}) X_{mn}(\bm{R}),\]

where $N_{\bm{R}}$ is the degeneracy of R vectors (not the total number of R vectors).

Arguments

  • Oᴿ: n_wann * n_wann * n_rvecs, real space operator matrix
  • kpoints: 3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
source
Wannier.mdrs_v1tov2Method
mdrs_v1tov2(Rvectors::RVectorsMDRS{T}, Oᴿ::Array{Complex{T},3}) where {T<:Real}

Expand an operator defined on MDRSv1 R vectors to MDRSv2 R̃ vectors.

Arguments

  • Rvectors: R vectors for MDRS interpolation
  • Oᴿ: n_wann * n_wann * n_rvecs, operator defined on R vectors

Return

  • Oᴿ̃: n_wann * n_wann * n_r̃vecs, operator defined on R̃ vectors
Note

This is useful when reading Hamiltonian from seedname_tb.dat and use MDRSv2 interpolation in Wannier.jl.

source

InterpModel

To separate the Wannier interpolation from Wannierization, we have another model, the InterpModel, which is solely for the purpose of interpolation, while Model works only for Wannierization.

Wannier.InterpModelType
struct InterpModel

The model for Wannier interpolation.

Store the real space matrices, e.g., the Hamiltonian $H(\bm{R})$.

Usually after Wannierization of a Model, we construct this InterpModel for Wannier interpolation of operators.

Fields

  • kRvectors: the kpoints and R vectors
  • kpath: the kpoint path for band structure
  • H: n_wann * n_wann * n_rvecs, the Hamiltonian in real space
  • r: n_wann * n_wann * n_rvecs, the position operator in real space
  • S: n_wann * n_wann * n_rvecs * 3, optional, the spin operator
Note

For MDRS interpolation, the n_rvecs in this docstring is actually the number of $\tilde{\bm{R}}$ vectors, i.e., KRVectors.Rvectors.n_r̃vecs. For WS interpolation, it is just KRVectors.Rvectors.n_rvecs.

source
Wannier.InterpModelMethod
InterpModel(model::Model; mdrs::Bool=true)

Construct a InterpModel from a Wannierization Model.

Arguments

  • model: the Wannierization Model

Keyword Arguments

  • mdrs: whether to use MDRS interpolation
  • kpath: if not given, use get_kpath to auto generate a kpath.
source
Wannier.InterpModelMethod
InterpModel(kRvectors, H, kpath, H, r)

A InterpModel constructor ignoring spin operator matrices.

Arguments

  • kRvectors: the kpoint and R vectors
  • kpath: the kpoint path for band structure
  • H: n_wann * n_wann * n_rvecs, the Hamiltonian in real space
  • r: n_wann * n_wann * n_rvecs * 3, the position operator in real space
source
Wannier.InterpModelMethod
InterpModel(kRvectors, H, kpath, H)

A InterpModel constructor ignoring position and spin operator matrices.

Arguments

  • kRvectors: the kpoint and R vectors
  • kpath: the kpoint path for band structure
  • H: n_wann * n_wann * n_rvecs, the Hamiltonian in real space
source

Band structure

Warning

Need some care when comparing the band interpolation between Wannier.jl and Wannier90, after running the Wannier.jl disentanglement and writing an optimized amn file for Wannier90 to interpolate band structure.

When Wannier90 read an amn file, it will (Wannier90 v3.1.0)

  1. Lowdin orthogonalize the amn matrices, in disentangle.F90:dis_project, line 1418 This should do no harm, since the optimized amn is already semi-unitary, a SVD of it should not change the optimized amn (apart from numerical noise)
  2. generate a new amn according to the frozen window, in disentangle.F90:dis_proj_froz, line 1830. This will DESTROY the optimized amn matrices, if we restart Wannier90 from the optimized amn with dis_num_iter = 0, the spreads in wout file is very different from the output of Wannier.jl, we must skip this step by commenting out ALL the dis_froz_min/max in the win file, then use Wannier90 to interpolate bands, remember also set num_iter and dis_num_ite = 0.
Wannier.cut_hamiltonianMethod
cut_hamiltonian(Rvectors::RVectorsMDRS{T}, Hᴿ::Array{Complex{T},3}, Rcut::T) where {T<:Real}

Cut real space Hamiltonian Hᴿ by the given cutoff radius.

Arguments

  • Rvectors: RVectorsMDRS
  • Hᴿ: n_wann * n_wann * n_rvecs, Hamiltonian in real space.
  • Rcut: cutoff radius in angstrom, the Hamiltonian H(R) corresponds to the R vectors with norm larger than Rcut will be set to 0.

Return

  • H: n_wann * n_wann * n_rvecs, cutted Hamiltonian in real space.

Example

R, H, pos = read_w90_tb("mos2")
# Expand the H(R) to R̃, I assume the H(R) is read from w90 tb.dat file
H1 = Wannier.mdrs_v1tov2(R, H)
H2 = cut_hamiltonian(R, H1, 7.0)
source
Wannier.diag_HkMethod
diag_Hk(Hᵏ:: AbstractArray{T, 3}) where {T<:Complex}

Diagonalize k space Hamiltonian H.

\[H = V E V^{-1}\]

Arguments

  • H: n_wann * n_wann * n_kpts, k space Hamiltonian

Return

  • E: n_wann * n_kpts, energy eigen values
  • V: n_wann * n_wann * n_kpts, V[:, i, ik] is the i-th eigen vector at ik-th kpoint
source
Wannier.get_HkMethod
get_Hk(E, U)

Construct k space Hamiltonian Hᵏ.

\[H_{\bm{k}} = U_{\bm{k}}^\dagger [\epsilon_{n \bm{k}}] U_{\bm{k}},\]

where $[\epsilon_{n \bm{k}}]$ is a diagonal matrix with $\epsilon_{n \bm{k}}$ as the diagonal elements.

Arguments

  • E: n_bands * n_kpts, energy eigenvalue
  • U: n_bands * n_wann * n_kpts, gauge matrices
source
Wannier.interpolateMethod
interpolate(model::InterpModel, kpi::KPathInterpolant)

Interpolate band structure along the given kpath.

Arguments

  • model: InterpModel
  • kpi: KPathInterpolant
source
Wannier.interpolateMethod
interpolate(model::InterpModel)

Interpolate band structure along kpath.

The model.kpath will be used.

Arguments

  • model: InterpModel
Note

The kpath has the same density as Wannier90's default, i.e., kpath_num_points = 100.

source
Wannier.interpolateMethod
interpolate(model::InterpModel{T}, kpoints::Matrix{T}) where {T<:Real}

Interpolate energy eigenvalues at kpoints.

Arguments

  • model: InterpModel
  • kpoints: 3 * n_kpts, kpoints to be interpolated, in fractional coordinates, can be nonuniform.
source
Wannier.sort_hamiltonian_by_normMethod
sort_hamiltonian_by_norm(Rvectors::RVectorsMDRS{T}, Hᴿ::Array{Complex{T},3}) where {T<:Real}

Sort the Hamiltonian $H(m{R})$ by its norm in descending order.

Arguments

  • Rvectors: RVectorsMDRS
  • Hᴿ: n_wann * n_wann * n_rvecs, Hamiltonian in real space.

Return

  • idx: the index to sort the Hamiltonian, i.e. Hᴿ[:, :, idx] is the sorted Hamiltonian.
  • normR: the norm of R vectors, normR[idx] is the sorted R vectors in descending order.
  • normH: the norm of H, normH[idx] is the sorted Hᴿ in descending order.

Example

R, H, pos = read_w90_tb("mos2")
# Expand the H(R) to R̃, I assume the H(R) is read from w90 tb.dat file
H1 = Wannier.mdrs_v1tov2(R, H)
idx, normR, normH = sort_hamiltonian_by_norm(R, H1)
source

Fermi surface

Note

The fermi_surface function will use WS or MDRS interpolation based on the type of Rvectors. However, Wannier90 only use WS interpolation when plotting Fermi surface (even if the use_ws_distance is set as true in the win file). So the fermi_surface function will output different result than Wannier90 if using MDRS.

Wannier.fermi_surfaceMethod
fermi_surface(Rvectors, H; n_k)

Interpolate Fermi surface.

Arguments

  • Rvectors: RVectors or RVectorsMDRS
  • H: n_wann * n_wann * n_r̃vecs, Hamiltonian in R space
  • n_k: integer or 3-vector, number of interpolated kpoints along three directions

Return

  • kpoints: 3 * n_kx * n_ky * n_kz, interpolated kpoints in fractional coordinates
  • E: n_wann * n_kx * n_ky * n_kz, interpolated eigenvalues
Note

The output n_kx = n_k + 1, since the bxsf format requires general grid, i.e., the last kpoint is the periodic image of the first one. This also restores the behavior of Wannier90.

Note

For MDRS interpolation, the H should be defined on the $\bm{R}$ vectors instead of the MDRSv2 $\tilde{\bm{R}}$ vectors; the expansion of $H(\bm{R})$ to $H(\tilde{\bm{R}})$ is done internally. This means that if you read the seedname_tb.dat file, then you can directly pass the H to this function.

source

Derivatives

Wannier._get_DMethod
_get_D(Rvectors, Hᴿ, kpoints; use_degen_pert=false, degen=1e-4)

Compute the matrix D in YWVS Eq. 25 (or Eq. 32 if use_degen_pert = true).

Arguments

  • Rvectors: RVectorsMDRS
  • Hᴿ: n_wann * n_wann * n_r̃vecs, the Hamiltonian matrix
  • kpoints: 3 * n_kpts, each column is a fractional kpoints coordinates

Keyword arguments

  • use_degen_pert: use perturbation treatment for degenerate eigenvalues
  • degen: degeneracy threshold in eV

Return

  • E: n_wann * n_kpts, energy eigenvalues
  • U: n_wann * n_wann * n_kpts, the unitary transformation matrix
  • Haᴴ: n_wann * n_wann * n_kpts * 3, the covariant part of derivative of Hamiltonian in Hamiltonian gauge, the $\bar{H}_{\alpha}^{(H)}$ in YWVS Eq. 26
  • D: n_wann * n_wann * n_kpts * 3, the matrix D in YWVS Eq. 25 or Eq. 32
Warning

If use_degen_pert = true, the degenerate subspace is rotated such that $\bar{H}_{x}^{(H)}$ is diagonal, note only the $x$ direction. I cannot diagonalize simultaneously all the three directions.

source
Wannier.effmass_fdMethod
effmass_fd(Rvectors, Hᴿ, kpoints; Δk=1e-3)

Compute the inverse of effective mass using finite differences of 2nd order.

Apply twice PRB 93, 205147 (2016) Eq. 80.

source
Wannier.get_d2H_dadbMethod
get_d2H_dadb(Rvectors, Hᴿ, kpoints)

Compute the second derivative of the Hamiltonian $H$ with respect to three Cartesian directions.

YWVS Eq. 28.

source
Wannier.get_dH_daMethod
get_dH_da(Rvectors, Hᴿ, kpoints)

Compute the derivative of the Hamiltonian $H$ with respect to three Cartesian directions.

YWVS Eq. 26.

source
Wannier.velocityMethod
velocity(Rvectors, Hᴿ, kpoints; use_degen_pert=false, degen=1e-4)

Compute velocity along three Cartesian directions.

YWVS Eq. 27.

Arguments

  • Rvectors: RVectorsMDRS
  • Hᴿ: n_wann * n_wann * n_r̃vecs, the Hamiltonian matrix
  • kpoints: 3 * n_kpts, each column is a fractional kpoints coordinates

Keyword arguments

  • use_degen_pert: use perturbation treatment for degenerate eigenvalues
  • degen: degeneracy threshold in eV

Return

  • E: n_wann * n_kpts, energy eigenvalues
  • v: n_wann * n_kpts * 3, velocity along three Cartesian directions, in unit hbar * m / s
Warning

Wannier90 by default set use_degen_pert = false. In 3D, and for N degenerate states, the velocity is a tensor of size N * N * 3, where 3 is for three Cartesian directions. Thus I cannot simultaneously diagonalize the tensor for all 3 directions. This means I can only use perturbation treatment for one of the directions, and only in that direction the velocity matrix is diagonal. So for degenerate states, the velocity is not well defined, and the results are meaningless, instead one should use the full velocity matrix which also include non-diagonal part, see get_dH_da.

source
Wannier.velocity_fdMethod
velocity_fd(Rvectors, Hᴿ, kpoints; Δk=1e-3)

Compute the velocity using finite differences of 2nd order.

PRB 93, 205147 (2016) Eq. 80.

source
  • WAN20Pizzi, G.; Vitale, V.; Arita, R.; Blügel, S.; Freimuth, F.; Géranton, G.; Gibertini, M.; Gresch, D.; Johnson, C.; Koretsune, T.; Ibañez-Azpiroz, J.; Lee, H.; Lihm, J.-M.; Marchand, D.; Marrazzo, A.; Mokrousov, Y.; Mustafa, J. I.; Nohara, Y.; Nomura, Y.; Paulatto, L.; Poncé, S.; Ponweiser, T.; Qiao, J.; Thöle, F.; Tsirkin, S. S.; Wierzbowska, M.; Marzari, N.; Vanderbilt, D.; Souza, I.; Mostofi, A. A. & Yates, J. R., Wannier90 as a community code: new features and applications, Journal of Physics: Condensed Matter, 2020, 32, 165902
  • WBERRITsirkin, S. S., High performance Wannier interpolation of Berry curvature and related quantities with WannierBerri code, npj Computational Materials, 2021, 7