Model
The Model is a Julia struct that abstracts a single Wannierization, it contains all the necessary information for maximally localization of the crystal structure.
Contents
Index
Wannier.CenterSpreadPenaltyWannier.ModelWannier.ModelWannier.SpreadWannier.SpreadCenterWannier.SpreadPenaltyBase.isapproxWannier.centerWannier.centerWannier.centerWannier.centerWannier.compute_berry_connection_kspaceWannier.isentangledWannier.isisolatedWannier.omegaWannier.omega_centerWannier.omega_gradWannier.omega_localWannier.position_opWannier.position_opWannier.position_op
Model
Wannier.Model — TypeModel(model, kstencil, overlaps)
Model(model, kstencil, overlaps, gauges)
Model(model, kstencil, overlaps, gauges, eigenvalues)
Model(
model,
kstencil,
overlaps,
gauges,
eigenvalues,
frozen_bands
)
Construct a Model from an existing Model, reuse the lattice and atomic information, but change the kstencil and overlaps.
For instance, use a cubic-6-neighbors KspaceStencil for parallel_transport.
Wannier.Model — Typestruct Model{T<:Real}A high-level data structure containing necessary parameters and matrices for constructing Wannier functions (WFs), or called, Wannierization.
In general, the problem of Wannierization is to find a set of unitary matrices $U_{\bm{k}}$ that gives a localized representation of the Bloch states $\psi_{n \bm{k}}$. Depending on the inputs, the Wannierization problem can be categorized into two classes:
- isolated manifold: when number of bands = number of Wannier functions
- entangled manifold: when number of bands > number of Wannier functions
Fields
Using these accronyms, - n_atoms: number of atoms - n_kpoints: number of kpoints - n_bvectors: number of b-vectors - n_bands: number of bands - n_wannier: number of Wannier functions the fields are defined as follows:
lattice: unit cell, 3 * 3, each column is a lattice vector in Å unitatom_positions: atomic positions, length-n_atomsvector, each element is aVec3of fractional coordinatesatom_labels: atomic labels, length-n_atomsvector of stringkstencil: stencil for finite differences on the kpoint grid, also called $\mathbf{b}$-vectors. Should satisfy completeness condition, seeKspaceStenciloverlaps: overlap matrices between neighboring wavefunctions, $M_{\bm{k},\bm{b}}$. Length-n_kpointsvector, each element is a length-n_bvectorsvector, then each element is an_bands * n_bandsmatrix. Also calledmmnmatrices in wannier90gauges: unitary or semi-unitary gauge transformation matrices, $U_{\bm{k}}$, or called the gauge matrices. Length-n_kpointsvector, each element is an_bands * n_wanniermatrixeigenvalues: energy eigenvalues, $\varepsilon_{n \bm{k}}$, length-n_kpointsvector, each element is a length-n_bandsvector, in eV unitfrozen_bands: mask for frozen bands. Length-n_kpointsvector, each element is a length-n_bandsBitVector. Iftruethe the state at that kpoint and band index is kept unchanged during the disentanglement procedure.entangled_bands: mask for bands taking part in disentanglement. Length-n_kpointsvector, each element is a length-n_bandsBitVector. Iftruethe the state at that kpoint and band index participates the disentanglement procedure.
Base.isapprox — Methodisapprox(a, b; kwargs...)
Compare two Model objects.
Wannier.isentangled — Methodisentangled(model)
Is entangled manifold?
Wannier.isisolated — Methodisisolated(model)
Is isolated manifold?
Spread
Wannier.CenterSpreadPenalty — TypePenalty for minimizing the spread as well as maximizing the "closeness" to the atoms.
Wannier.Spread — Typestruct SpreadThe Marzari-Vanderbilt (MV) spread functional.
From MV:
- $\Omega = \sum_n \langle r^2 \rangle_n - | \langle r \rangle_n |^2$
- $\langle r \rangle_n = -\frac{1}{N} \sum_{\bm{k},\bm{b}} w_{\bm{b}} \bm{b} \Im \log M_{nn}^{\bm{k},\bm{b}}$
- $\langle r^2 \rangle_n = \frac{1}{N} \sum_{\bm{k},\bm{b}} w_{\bm{b}} \bm{b} \left[ \left( 1 - | M_{nn}^{\bm{k},\bm{b}} |^2 \right) + \left( \Im \log M_{nn}^{\bm{k},\bm{b}} \right)^2 \right]$
Fields
Ω: total spread, unit ŲΩI: gauge-invarient part, unit ŲΩOD: off-diagonal part, unit ŲΩD: diagonal part, unit ŲΩ̃: Ω̃ = ΩOD + ΩD, unit Ųω: Ω of each WF, unit Ų,length(ω) = n_wannr: WF center, Cartesian coordinates, unit Å,3 * n_wann
Wannier.SpreadCenter — Typestruct SpreadCenterA struct containing both Spread and WF center penalty.
Wannier.SpreadPenalty — TypeStandard penalty for minimizing the total spread.
Wannier.center — Methodcenter(bvectors, M, U)Compute WF center in reciprocal space.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarray
Wannier.center — Methodcenter(model)Compute WF center in reciprocal space for Model.
Wannier.center — Methodcenter(bvectors, M, U)Compute WF center in reciprocal space.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarray
Wannier.center — Methodcenter(model, U)Compute WF center in reciprocal space for Model with given U gauge.
Arguments
model: theModelU:n_wann * n_wann * n_kptsarray
Wannier.compute_berry_connection_kspace — Methodcompute_berry_connection_kspace(
kstencil,
overlaps,
gauges;
imlog_diag,
force_hermiticity
)
Wannier-gauge Berry connection in kspace, WYSV Eq. 44 or MV Eq. C14
Keyword arguments
imlog_diag: use imaginary part of logrithm for diagonal elements, MV1997 Eq. 31. wannier90 default is true.
Wannier.omega — Methodomega(model, [U])
omega(bvectors, M, U)Compute WF spread for a Model, potentially for a given gauge U, or by explicitely giving bvectors and M. In case of the first bvectors = model.bvectors and M = model.M.
Wannier.omega_center — Methodomega_center(bvectors, M, U, r₀, λ)Compute WF spread with center penalty, for maximal localization.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarrayr₀:3 * n_wann, WF centers in cartesian coordinatesλ: penalty strength
Wannier.omega_grad — Methodomega_grad(bvectors, M, U, r)Compute gradient of WF spread.
Size of output dΩ/dU = n_bands * n_wann * n_kpts.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarrayr:3 * n_wann, the current WF centers in cartesian coordinates
Wannier.omega_local — Methodomega_local(bvectors, M, U)Local part of the contribution to r^2.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarray
Wannier.position_op — Methodposition_op(model)Compute WF postion operator matrix in reciprocal space for Model.
Wannier.position_op — Methodposition_op(bvectors, M, U)Compute WF postion operator matrix in reciprocal space.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarray
Wannier.position_op — Methodposition_op(model, U)Compute WF postion operator matrix in reciprocal space for Model with given U gauge.
Arguments
model: theModelU:n_wann * n_wann * n_kptsarray