Fourier series

Fourier series represent functions as linear combinations of sinusoids whose frequencies are integer multiples of a fundamental mode. In Wannier interpolation, the periodicity of an operator, $\hat{O}$, in a crystal admits a Fourier series expansion of the Bloch operator $\hat{O}(\bm{k}) = e^{-i\bm{k}\cdot\hat{\bm{r}}} \hat{O} e^{i\bm{k}\cdot\hat{\bm{r}}}$ as linear combinations of sinusoids at integer multiples of the real-space lattice vectors, forming a Bravais lattice. When exponentially localized Wannier functions exist, it suffices to use a modest number of lattice vectors in this expansion. The goal of this page of documentation is to describe the features, interface, and conventions of Fourier series evaluation as implemented by this library.

Lattice vectors

It is conventional to construct a basis for a reciprocal lattice $\{\bm{b}_j\}$ from a basis for a real-space lattice $\{\bm{a}_i\}$ such that

\[\bm{b}_j \cdot \bm{a}_i = 2\pi\delta_{ij}\]

Then we write the momentum space variable in the reciprocal lattice vector basis

\[\bm{k} = \sum_{j=1}^{d} k_j \bm{b}_j\]

Without loss of generality, the $k_j$ can be taken in the domain $[0,1]$. Additionally, any real space lattice vector can be written as an integer linear combination of the $\{\bm{a}_i\}$. We lastly mention that when computing Brillouin zone integrals, we keep track of this change of basis with the Jacobian determinant of the transformation,

\[\int_{\text{BZ}} d \bm{k} = |\det{B}| \int_{[0,1]^3} d^3 k.\]

Here, $B$ is a matrix whose columns are the $\{\bm{b}_j\}$.

Series coefficients

Wannier-interpolated Hamiltonians are small matrices obtained by projecting the ab-initio Hamiltonian onto a subspace of bands with energies near the Fermi level, a process called downfolding. These Hamiltonians can be expressed by a Fourier series as in the sum below

\[H(\bm{k}) = \sum_{\bm{R}} e^{i\bm{k}\cdot\bm{R}} H_{\bm{R}}\]

where the coefficients $H_{\bm{R}}$ are the matrix-valued Fourier coefficients. Truncating the sum over $\bm{R}$ at a modest number of modes can be done for Wannier Hamiltonians in the maximally-localized orbital basis, for which $H(\bm{k})$ is a smooth and periodic function and thus the truncation error of its Fourier series converges super-algebraically with respect to the number of modes.

Hamiltonian recipe

In model systems, a Bloch Hamiltonian can often be written down analytically. The recipe to write it as a Fourier series has two-steps

  1. Identify the real and reciprocal Bravais lattices, $\{\bm{a}_i\}$ and $\{\bm{b}_j\}$, and rewrite all of the phase dependences of the Hamiltonian as $\bm{k}\cdot\bm{R}$ with $\bm{k}$ written in the basis of $\{\bm{b}_j\}$ and the lattice vector $\bm{R}$ written as an integer linear combination of $\{\bm{a}_i\}$.
  2. Factor the Hamiltonian into a linear combination of normal modes indexed by the distinct $\bm{R}$ vectors. If the Hamiltonian is matrix-valued, this can be done one matrix element at a time.

For a non-trivial example, see the DOS interpolation for Graphene.

Interface

See FourierSeriesEvaluators.jl for the AbstractFourierSeries interface, which allows evaluation of the series with a function-like f(x) syntax. The Wannier90 interface will automatically construct these objects from the necessary output files.

Types

The concrete types listed below all implement the AbstractFourierSeries interface and should cover most use cases. They are organized by several abstract types to provide indicate physical properties of the Wannier-interpolated operator, such as a coordinate system or a gauge.

AutoBZ.AbstractWannierInterpType
AbstractWannierInterp{N,T,iip} <: AbstractFourierSeries{N,T,iip}

Abstract supertype for all Wannier-interpolated quantities in AutoBZ

source
AutoBZ.AbstractGaugeInterpType
AbstractGaugeInterp{G,N,T,iip} <: AbstractWannierInterp{N,T,iip}

An abstract subtype of AbstractFourierSeries representing Fourier series evaluators for Wannier-interpolated quantities with a choice of gauge, G, which is typically Hamiltonian or Wannier. A gauge is a choice of basis for the function space of the operator. For details, see to_gauge!.

source
AutoBZ.AbstractHamiltonianInterpType
AbstractHamiltonianInterp{G,N,T,iip,isherm} <: AbstractGaugeInterp{G,N,T,iip}

Abstract type representing Hamiltonians, which are matrix-valued Hermitian Fourier series. They should also have period 1, but produce derivatives with wavenumber 1 (not $2\pi$), in order to be consistent with definitions of the velocity operator.

source
AutoBZ.HamiltonianInterpType
HamiltonianInterp(f::AbstractFourierSeries; gauge=:Wannier)

A wrapper for FourierSeries with an additional gauge that allows for convenient diagonalization of the result. For details see to_gauge!.

source
AutoBZ.SOCHamiltonianInterpType
SOCHamiltonianInterp(f::Freq2RadSeries, λ; gauge=Wannier())

A wrapper for a Fourier series in a given gauge that has spin-orbit coupling represented by the matrix λ. In particular, this interpolates [f(k) 0; 0 f(k)] + λ.

source
AutoBZ.BerryConnectionInterpType
BerryConnectionInterp{P}(a::ManyFourierSeries, B; coord)

Interpolate the Berry connection in basis coord. a must evaluate the components of the connection in coordinate basis P, and B is the coordinate transformation from P to coord.

source
AutoBZ.AbstractVelocityInterpType
AbstractVelocityInterp{C,B,G,N,T,iip} <:AbstractCoordInterp{B,G,N,T,iip}

An abstract subtype of AbstractCoordInterp also containing information the velocity component, C, which is typically Whole, Inter, or Intra. These choices refer to the diagonal (intra) or off-diagonal (inter) matrix elements of the velocity operator in the eigebasis of H(k). For details see to_vcomp_gauge!. Since the velocity depends on the Hamiltonian, subtypes should evaluate (H(k), (v_1(k), v_2(k), ...)).

source
AutoBZ.GradientVelocityInterpType
GradientVelocityInterp(h::AbstractHamiltonianInterp, A; gauge, coord, vcomp)

Evaluate the Hamiltonian and its gradient, which doesn't produce gauge-covariant velocities. The Hamiltonian h must be in the Wannier gauge, but this will give the result in the requested gauge. A must be the coordinate transformation from the lattice basis to the desired coord system. vcomp selects the contribution to the velocities.

source
AutoBZ.CovariantVelocityInterpType
CovariantVelocityInterp(hv::GradientVelocityInterp, a::BerryConnectionInterp)

Uses the Berry connection to return fully gauge-covariant velocities. Returns a tuple of the Hamiltonian and the three velocity matrices.

source

Gauges

AutoBZ.WannierType
Wannier <: AbstractGauge

Singleton type representing the Wannier gauge. In this gauge, the Fourier series representation of the operator is evaluated as-is, usually resulting in a dense matrix at each $\bm{k}$ point. When evaluating a Green's function, choosing this gauge requires inverting a dense matrix, which scales as $\mathcal{O}(N^3)$ where $N$ is the number of orbitals.

source
AutoBZ.HamiltonianType
Hamiltonian() <: AbstractGauge

Singleton type representing the Hamiltonian gauge. This gauge is the eigenbasis of H, and all operators will be rotated to this basis. The Hamiltonian will also be returned as an Eigen factorization. In this basis, H is a diagonal matrix, so calculating a Green's function is an $\mathcal{O}(N)$ operation where $N$ is the number of bands.

source

Coordinate systems

AutoBZ.CartesianType
Cartesian <: AbstractCoordinate

Singleton type representing Cartesian coordinates.

source
AutoBZ.LatticeType
Lattice <: AbstractCoordinate

Singleton type representing lattice coordinates. The matrix $B$ whose columns are reciprocal lattice vectors converts this basis to the Cartesian basis.

source

Velocity components

AutoBZ.WholeType
Whole <: AbstractVelocityComponent

Singleton type representing whole velocities, i.e. the total velocity operator.

source
AutoBZ.IntraType
Intra <: AbstractVelocityComponent

Singleton type representing intraband velocities, which are the off-diagonal terms of the velocity operator in the Hamiltonian gauge

source
AutoBZ.InterType
Inter <: AbstractVelocityComponent

Singleton type representing interband velocities, which are the diagonal terms of the velocity operator in the Hamiltonian gauge

source

Special matrix types

AutoBZ.SOCMatrixType
SOCMatrix(A)

Wrapper for a matrix A that should be used as a block-diagonal matrix [A 0 0 A] as when dealing with spin-orbit coupling

source
AutoBZ.SSymmetricCompactType
SSymmetricCompact{N, T, L} <: StaticMatrix{N, N, T}

A StaticArray subtype that represents a Symmetric matrix. Unlike LinearAlgebra.Symmetric, SSymmetricCompact stores only the lower triangle of the matrix (as an SVector). The lower triangle is stored in column-major order. For example, for an SSymmetricCompact{3}, the indices of the stored elements can be visualized as follows:

┌ 1 ⋅ ⋅ ┐
| 2 4 ⋅ |
└ 3 5 6 ┘

Type parameters:

  • N: matrix dimension;
  • T: element type for lower triangle;
  • L: length of the SVector storing the lower triangular elements.

Note that L is always the Nth triangular number.

An SSymmetricCompact may be constructed either:

  • from an AbstractVector containing the lower triangular elements; or
  • from a Tuple containing both upper and lower triangular elements in column major order; or
  • from another StaticMatrix.

For the latter two cases, only the lower triangular elements are used; the upper triangular elements are ignored.

source