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 inWannier90[WAN20]. This is slower than the WS interpolation.MDRSv2: a slightly modified version which removes some for loops, same as what is implemented inWannier Berri[WBERRI]. The speed is similar to the WS interpolation.
By default, we use the MDRSv2 interpolation.
Contents
Index
Wannier.InterpModelWannier.InterpModelWannier.InterpModelWannier.InterpModelWannier.KRVectorsWannier.KRVectorsWannier.KRVectorsWannier.KRVectorsWannier.RVectorsWannier.RVectorsWannier.RVectorsMDRSWannier.RVectorsMDRSWannier._fourier_mdrs_v1Wannier._fourier_mdrs_v2Wannier._get_DWannier._invfourier_mdrs_v1Wannier._invfourier_mdrs_v1!Wannier._invfourier_mdrs_v2Wannier._invfourier_mdrs_v2!Wannier.check_weightsWannier.cut_hamiltonianWannier.diag_HkWannier.effmass_fdWannier.fermi_surfaceWannier.fermi_surfaceWannier.fourierWannier.fourierWannier.get_HkWannier.get_Rvectors_mdrsWannier.get_Rvectors_wsWannier.get_d2H_dadbWannier.get_dH_daWannier.interpolateWannier.interpolateWannier.interpolateWannier.invfourierWannier.invfourierWannier.invfourierWannier.invfourierWannier.invfourier!Wannier.mdrs_v1tov2Wannier.sort_hamiltonian_by_normWannier.velocityWannier.velocity_fd
R vectors
RVectors
The R vectors are automatically generated by calling either of these two approaches:
read_w90_interpto read theWannier90outputseedname.chk.fmtfile, and generate anInterpModel,InterpModela constructor to generate theInterpModelfrom a WannierizationModel.
So usually the user does not need to call the following functions.
Wannier.RVectors — Typestruct RVectorsThe R vectors for interpolation.
Fields
lattice: columns are lattice vectorsgrid: number of FFT grid points in each direction, actually equal tokgridR:3 * n_rvecs, R vectors in fractional coordinates w.r.t. latticeN:n_rvecs, degeneracy of each R vector
The R vectors are sorted in the same order as Wannier90.
Wannier.RVectors — MethodRVectors(lattice, grid, R, N)Constructor for RVectors.
Auto transform lattice and grid to Mat3 and Vec3, respectively.
Arguments
lattice: columns are lattice vectorsgrid: number of FFT grid points in each direction, actually equal tokgridR:3 * n_rvecs, R vectors in fractional coordinates w.r.t. latticeN:n_rvecs, degeneracy of each R vector
Wannier.RVectorsMDRS — Typestruct RVectorsMDRSThe R vectors for MDRS interpolation.
Fields
Rvectors:Rvectorsfor Wigner-Seitz interpolation
For MDRSv1
T:n_wann * n_wann * n_rvecs, translation vectors w.r.t to lattice for MDRSv1Nᵀ:n_wann * n_wann * n_rvecs, degeneracy of eachTvector for MDRSv1
For MDRSv2
R̃vectors:RVectorscontaining expanded set ofR + Tvectors used in MDRSv2R̃_RT: mapping ofR̃vectorstoRandTvectors
The R vectors are sorted in the same order as Wannier90.
Wannier.RVectorsMDRS — MethodRVectorsMDRS(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:RVectorsfor Wigner-Seitz interpolationT:n_wann * n_wann * n_rvecs, translation vectors w.r.t to lattice for MDRSv1Nᵀ:n_wann * n_wann * n_rvecs, degeneracy of eachTvector for MDRSv1
Wannier.check_weights — Methodcheck_weights(R::RVectors)Sanity check for the degeneracies of R vectors.
Wannier.get_Rvectors_mdrs — Methodget_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 vectorsrgrid: number of FFT grid points in each direction, actually equal tokgridcenters:3 * n_wann, WF centers in fractional coordinates
Keyword arguments
atol: toerance for checking degeneracy, equivalent toWannier90input parameterws_distance_tolmax_cell: number of neighboring cells to be searched, equivalent toWannier90input parameterws_search_size
Wannier.get_Rvectors_ws — Methodget_Rvectors_ws(lattice, rgrid; atol=1e-5, max_cell=3)Generate R vectors for Wigner-Seitz interpolation.
Arguments
lattice: columns are lattice vectorsrgrid: number of FFT grid points in each direction, actually equal tokgrid
Keyword arguments
atol: toerance for checking degeneracy, equivalent toWannier90input parameterws_distance_tolmax_cell: number of neighboring cells to be searched, equivalent toWannier90input parameterws_search_size
KRVectors
The KRVectors contains both kpoint mappings and RVectors, allowing both forward and inverse Fourier transform.
Wannier.KRVectors — Typestruct KRVectorsContains both kpoints mapping and R vectors.
Fields
lattice: each column is a lattice vectorkgrid: number of kpoints along 3 directionskpoints: each column is a kpoint in fractional coordinatesk_xyz: kpoints mappings fromikto[ikx, iky, ikz]xyz_k: kpoints mappings from[ikx, iky, ikz]toikRvectors: R vectors, can be eitherRVectorsorRVectorsMDRSrecip_lattice: each column is a reciprocal lattice vectorn_kpts: number of kpointsn_rvecs: number of R vectorsn_r̃vecs: number of R̃ vectors, only for MDRS
Wannier.KRVectors — MethodKRVectors(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 vectorkgrid: number of kpoints along 3 directionskpoints: each column is a kpoint in fractional coordinatesRvectors: R vectors, can be eitherRVectorsorRVectorsMDRS
Wannier.KRVectors — MethodKRVectors(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 vectorkgrid: number of kpoints along 3 directionskpoints: each column is a kpoint in fractional coordinatesk_xyz: kpoints mappings fromikto[ikx, iky, ikz]xyz_k: kpoints mappings from[ikx, iky, ikz]toikRvectors: R vectors, can be eitherRVectorsorRVectorsMDRS
Wannier.KRVectors — MethodKRVectors(Rvectors)Construct a KRVectors from RVectors or RVectorsMDRS.
The kpoints part are just empty.
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.
Fourier transforms
Performs Fourier transform between k and R spaces defined by KRVectors.
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_v1 — Method_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 interpolationOᵏ:n_wann * n_wann * n_kpts, operator matrix in k space
Wannier._fourier_mdrs_v2 — Method_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 interpolationOᵏ:n_wann * n_wann * n_kpts, operator matrix in k space
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.
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
Wannier._invfourier_mdrs_v1 — Method_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 matrixkpoints:3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
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.
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
Wannier._invfourier_mdrs_v2 — Method_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 matrixkpoints:3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
Wannier.fourier — Methodfourier(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 interpolationOᵏ:n_wann * n_wann * n_kpts, operator matrix in k space
Keyword Arguments
version::v1or:v2, default is:v2
See also _fourier_mdrs_v1 and _fourier_mdrs_v2.
Wannier.fourier — Methodfourier(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 interpolationOᵏ:n_wann * n_wann * n_kpts, operator matrix in k space
Wannier.invfourier! — Methodinvfourier!(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
Wannier.invfourier — Methodinvfourier(
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.
Wannier.invfourier — Methodinvfourier(
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.
Wannier.invfourier — Methodinvfourier(
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 interpolationOᴿ:n_wann * n_wann * n_rvecs, operator matrix in R spacekpoints:3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
Keyword Arguments
version::v1or:v2, default is:v2
See also _invfourier_mdrs_v1 and _invfourier_mdrs_v2.
Wannier.invfourier — Methodinvfourier(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 matrixkpoints:3 * n_kpts, kpoints to be interpolated, fractional coordinates, can be nonuniform
Wannier.mdrs_v1tov2 — Methodmdrs_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 interpolationOᴿ: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
This is useful when reading Hamiltonian from seedname_tb.dat and use MDRSv2 interpolation in Wannier.jl.
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.InterpModel — Typestruct InterpModelThe 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 vectorskpath: the kpoint path for band structureH:n_wann * n_wann * n_rvecs, the Hamiltonian in real spacer:n_wann * n_wann * n_rvecs, the position operator in real spaceS:n_wann * n_wann * n_rvecs * 3, optional, the spin operator
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.
Wannier.InterpModel — MethodInterpModel(model::Model; mdrs::Bool=true)Construct a InterpModel from a Wannierization Model.
Arguments
model: the WannierizationModel
Keyword Arguments
mdrs: whether to use MDRS interpolationkpath: if not given, useget_kpathto auto generate a kpath.
Wannier.InterpModel — MethodInterpModel(kRvectors, H, kpath, H, r)A InterpModel constructor ignoring spin operator matrices.
Arguments
kRvectors: the kpoint and R vectorskpath: the kpoint path for band structureH:n_wann * n_wann * n_rvecs, the Hamiltonian in real spacer:n_wann * n_wann * n_rvecs * 3, the position operator in real space
Wannier.InterpModel — MethodInterpModel(kRvectors, H, kpath, H)A InterpModel constructor ignoring position and spin operator matrices.
Arguments
kRvectors: the kpoint and R vectorskpath: the kpoint path for band structureH:n_wann * n_wann * n_rvecs, the Hamiltonian in real space
Band structure
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)
- Lowdin orthogonalize the
amnmatrices, indisentangle.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) - generate a new
amnaccording to the frozen window, indisentangle.F90:dis_proj_froz, line 1830. This will DESTROY the optimizedamnmatrices, if we restartWannier90from the optimizedamnwithdis_num_iter = 0, the spreads inwoutfile is very different from the output ofWannier.jl, we must skip this step by commenting out ALL thedis_froz_min/maxin thewinfile, then useWannier90to interpolate bands, remember also setnum_iteranddis_num_ite = 0.
Wannier.cut_hamiltonian — Methodcut_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:RVectorsMDRSHᴿ:n_wann * n_wann * n_rvecs, Hamiltonian in real space.Rcut: cutoff radius in angstrom, the HamiltonianH(R)corresponds to the R vectors with norm larger thanRcutwill 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)Wannier.diag_Hk — Methoddiag_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 valuesV:n_wann * n_wann * n_kpts,V[:, i, ik]is the i-th eigen vector atik-th kpoint
Wannier.get_Hk — Methodget_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 eigenvalueU:n_bands * n_wann * n_kpts, gauge matrices
Wannier.interpolate — Methodinterpolate(model::InterpModel, kpi::KPathInterpolant)Interpolate band structure along the given kpath.
Arguments
model:InterpModelkpi:KPathInterpolant
Wannier.interpolate — Methodinterpolate(model::InterpModel)Interpolate band structure along kpath.
The model.kpath will be used.
Arguments
model:InterpModel
The kpath has the same density as Wannier90's default, i.e., kpath_num_points = 100.
Wannier.interpolate — Methodinterpolate(model::InterpModel{T}, kpoints::Matrix{T}) where {T<:Real}Interpolate energy eigenvalues at kpoints.
Arguments
model:InterpModelkpoints:3 * n_kpts, kpoints to be interpolated, in fractional coordinates, can be nonuniform.
Wannier.sort_hamiltonian_by_norm — Methodsort_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:RVectorsMDRSHᴿ: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 sortedRvectors in descending order.normH: the norm of H,normH[idx]is the sortedHᴿ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)Fermi surface
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_surface — MethodWannier.fermi_surface — Methodfermi_surface(Rvectors, H; n_k)Interpolate Fermi surface.
Arguments
Rvectors:RVectorsorRVectorsMDRSH:n_wann * n_wann * n_r̃vecs, Hamiltonian in R spacen_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 coordinatesE:n_wann * n_kx * n_ky * n_kz, interpolated eigenvalues
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.
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.
Derivatives
Wannier._get_D — Method_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:RVectorsMDRSHᴿ:n_wann * n_wann * n_r̃vecs, the Hamiltonian matrixkpoints:3 * n_kpts, each column is a fractional kpoints coordinates
Keyword arguments
use_degen_pert: use perturbation treatment for degenerate eigenvaluesdegen: degeneracy threshold in eV
Return
E:n_wann * n_kpts, energy eigenvaluesU:n_wann * n_wann * n_kpts, the unitary transformation matrixHaᴴ: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. 26D:n_wann * n_wann * n_kpts * 3, the matrix D in YWVS Eq. 25 or Eq. 32
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.
Wannier.effmass_fd — Methodeffmass_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.
Wannier.get_d2H_dadb — Methodget_d2H_dadb(Rvectors, Hᴿ, kpoints)Compute the second derivative of the Hamiltonian $H$ with respect to three Cartesian directions.
YWVS Eq. 28.
Wannier.get_dH_da — Methodget_dH_da(Rvectors, Hᴿ, kpoints)Compute the derivative of the Hamiltonian $H$ with respect to three Cartesian directions.
YWVS Eq. 26.
Wannier.velocity — Methodvelocity(Rvectors, Hᴿ, kpoints; use_degen_pert=false, degen=1e-4)Compute velocity along three Cartesian directions.
YWVS Eq. 27.
Arguments
Rvectors:RVectorsMDRSHᴿ:n_wann * n_wann * n_r̃vecs, the Hamiltonian matrixkpoints:3 * n_kpts, each column is a fractional kpoints coordinates
Keyword arguments
use_degen_pert: use perturbation treatment for degenerate eigenvaluesdegen: degeneracy threshold in eV
Return
E:n_wann * n_kpts, energy eigenvaluesv:n_wann * n_kpts * 3, velocity along three Cartesian directions, in unithbar * m / s
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.
Wannier.velocity_fd — Methodvelocity_fd(Rvectors, Hᴿ, kpoints; Δk=1e-3)Compute the velocity using finite differences of 2nd order.
PRB 93, 205147 (2016) Eq. 80.
- 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