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
struct Model

A struct containing the parameters and matrices of the crystal structure.

Fields

  • lattice: columns are the lattice vectors
  • atom_positions: columns are the fractional coordinates of atoms
  • atom_labels: labels of atoms
  • kgrid: number of kpoints along 3 lattice vectors
  • kpoints: columns are the fractional coordinates of kpoints
  • bvectors: bvectors satisfying the B1 condition
  • frozen_bands: indicates which bands are frozen
  • M: 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 vectors
  • n_atoms: number of atoms
  • n_bands: number of bands
  • n_wann: number of wannier functions
  • n_kpts: number of kpoints
  • n_bvecs: number of bvectors
Note

This only cotains the necessary information for maximal localization. For Wannier interpolation, see InterpModel.

source
Wannier.ModelMethod
Model(lattice, atom_positions, atom_labels, kgrid, kpoints, bvectors, frozen_bands, M, U, E)

Construct a Model struct.

Arguments

  • lattice: columns are the lattice vectors
  • atom_positions: columns are the fractional coordinates of atoms
  • atom_labels: labels of atoms
  • kgrid: number of kpoints along 3 lattice vectors
  • kpoints: columns are the fractional coordinates of kpoints
  • bvectors: bvectors satisfying the B1 condition
  • frozen_bands: indicates which bands are frozen
  • M: 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}}$
Tip

This is more user-friendly constructor, only necessary information is required. Remaining fields are generated automatically.

source
Wannier.rotate_gaugeMethod
rotate_gauge(model::Model, U::Array{T,3}; diag_H=false)

Rotate the gauge of a Model.

Arguments

  • model: a Model struct
  • U: 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 to model.E, and the inverse of the eigenvectors to model.U, so that the model is still in the input gauge U. Otherwise, if the rotated Hamiltonian is not diagonal, raise error.
Note

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.

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.berry_connectionMethod
berry_connection(bvectors, M, U)

Compute Berry connection at each kpoint.

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.omegaMethod
omega(bvectors, M, U)

Compute WF spread.

Arguments

  • bvectors: bvecoters
  • M: n_bands * n_bands * * n_bvecs * n_kpts overlap array
  • U: n_wann * n_wann * n_kpts array
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_gradMethod
omega_grad(bvectors, M, U)

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
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