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.ModelWannier.ModelWannier.SpreadWannier.berry_connectionWannier.centerWannier.centerWannier.centerWannier.omegaWannier.omegaWannier.omegaWannier.omega_gradWannier.omega_gradWannier.omega_localWannier.position_opWannier.position_opWannier.position_opWannier.rotate_gauge
Model
Wannier.Model — Typestruct ModelA struct containing the parameters and matrices of the crystal structure.
Fields
lattice: columns are the lattice vectorsatom_positions: columns are the fractional coordinates of atomsatom_labels: labels of atomskgrid: number of kpoints along 3 lattice vectorskpoints: columns are the fractional coordinates of kpointsbvectors: bvectors satisfying the B1 conditionfrozen_bands: indicates which bands are frozenM:n_bands * n_bands * n_bvecs * n_kpts, overlap matrix $M_{\bm{k},\bm{b}}$U:n_bands * n_wann * n_kpts, (semi-)unitary gauge rotation matrix $U_{\bm{k}}$E:n_bands * n_kpts, energy eigenvalues $\epsilon_{n \bm{k}}$recip_lattice: columns are the reciprocal lattice vectorsn_atoms: number of atomsn_bands: number of bandsn_wann: number of wannier functionsn_kpts: number of kpointsn_bvecs: number of bvectors
This only cotains the necessary information for maximal localization. For Wannier interpolation, see InterpModel.
Wannier.Model — MethodModel(lattice, atom_positions, atom_labels, kgrid, kpoints, bvectors, frozen_bands, M, U, E)Construct a Model struct.
Arguments
lattice: columns are the lattice vectorsatom_positions: columns are the fractional coordinates of atomsatom_labels: labels of atomskgrid: number of kpoints along 3 lattice vectorskpoints: columns are the fractional coordinates of kpointsbvectors: bvectors satisfying the B1 conditionfrozen_bands: indicates which bands are frozenM:n_bands * n_bands * n_bvecs * n_kpts, overlap matrix $M_{\bm{k},\bm{b}}$U:n_bands * n_wann * n_kpts, (semi-)unitary gauge rotation matrix $U_{\bm{k}}$E:n_bands * n_kpts, energy eigenvalues $\epsilon_{n \bm{k}}$
This is more user-friendly constructor, only necessary information is required. Remaining fields are generated automatically.
Wannier.rotate_gauge — Methodrotate_gauge(model::Model, U::Array{T,3}; diag_H=false)Rotate the gauge of a Model.
Arguments
model: aModelstructU:n_bands * n_wann * n_kpts, (semi-)unitary gauge rotation matrix $U_{\bm{k}}$
Keyword Arguments
diag_H: if after rotation, the Hamiltonian is not diagonal, then diagonalize it and save the eigenvalues tomodel.E, and the inverse of the eigenvectors tomodel.U, so that themodelis still in the input gaugeU. Otherwise, if the rotated Hamiltonian is not diagonal, raise error.
The original Model.U will be discarded; the M, and E matrices will be rotated by the input U. However, since E is not the Hamiltonian matrices but only the eigenvalues, if diag_H = false, this function only support rotations that keep the Hamiltonian in diagonal form.
Spread
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.berry_connection — Methodberry_connection(bvectors, M, U)Compute Berry connection at each kpoint.
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.omega — Methodomega(model)Compute WF spread for Model.
Wannier.omega — Methodomega(bvectors, M, U)Compute WF spread.
Arguments
bvectors: bvecotersM:n_bands * n_bands * * n_bvecs * n_kptsoverlap arrayU:n_wann * n_wann * n_kptsarray
Wannier.omega — Methodomega(model, U)Compute WF spread for Model using given U gauge.
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_grad — Methodomega_grad(bvectors, M, U)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_kptsarray
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