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

Model

Wannier.ModelType
Model(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.

source
Wannier.ModelType
struct 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 Å unit

  • atom_positions: atomic positions, length-n_atoms vector, each element is a Vec3 of fractional coordinates

  • atom_labels: atomic labels, length-n_atoms vector of string

  • kstencil: stencil for finite differences on the kpoint grid, also called $\mathbf{b}$-vectors. Should satisfy completeness condition, see KspaceStencil

  • overlaps: overlap matrices between neighboring wavefunctions, $M_{\bm{k},\bm{b}}$. Length-n_kpoints vector, each element is a length-n_bvectors vector, then each element is a n_bands * n_bands matrix. Also called mmn matrices in wannier90

  • gauges: unitary or semi-unitary gauge transformation matrices, $U_{\bm{k}}$, or called the gauge matrices. Length-n_kpoints vector, each element is a n_bands * n_wannier matrix

  • eigenvalues: energy eigenvalues, $\varepsilon_{n \bm{k}}$, length-n_kpoints vector, each element is a length-n_bands vector, in eV unit

  • frozen_bands: mask for frozen bands. Length-n_kpoints vector, each element is a length-n_bands BitVector. If true the 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_kpoints vector, each element is a length-n_bands BitVector. If true the the state at that kpoint and band index participates the disentanglement procedure.

source

Spread

Wannier.SpreadType
struct Spread

The 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_wann
  • r: WF center, Cartesian coordinates, unit Å, 3 * n_wann
source
Wannier.centerMethod
center(bvectors, M, U)

Compute WF center in reciprocal space.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
source
Wannier.centerMethod
center(bvectors, M, U)

Compute WF center in reciprocal space.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
source
Wannier.centerMethod
center(model, U)

Compute WF center in reciprocal space for Model with given U gauge.

Arguments

  • model: the Model
  • U: n_wann * n_wann * n_kpts array
source
Wannier.compute_berry_connection_kspaceMethod
compute_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.
source
Wannier.omegaMethod
omega(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.

source
Wannier.omega_centerMethod
omega_center(bvectors, M, U, r₀, λ)

Compute WF spread with center penalty, for maximal localization.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
  • r₀: 3 * n_wann, WF centers in cartesian coordinates
  • λ: penalty strength
source
Wannier.omega_gradMethod
omega_grad(bvectors, M, U, r)

Compute gradient of WF spread.

Size of output dΩ/dU = n_bands * n_wann * n_kpts.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
  • r: 3 * n_wann, the current WF centers in cartesian coordinates
source
Wannier.omega_localMethod
omega_local(bvectors, M, U)

Local part of the contribution to r^2.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
source
Wannier.position_opMethod
position_op(bvectors, M, U)

Compute WF postion operator matrix in reciprocal space.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
source
Wannier.position_opMethod
position_op(model, U)

Compute WF postion operator matrix in reciprocal space for Model with given U gauge.

Arguments

  • model: the Model
  • U: n_wann * n_wann * n_kpts array
source