Core of TightBindingApproximation

TightBindingApproximation.basicoptionsConstant
const basicoptions = (
    gauge = "gauge used to perform the Fourier transformation",
    infinitesimal = "infinitesimal added to the diagonal of the matrix representation of the tight-binding Hamiltonian"
)

Basic options of tight-binding actions.

source
QuantumLattices.DegreesOfFreedom.MetricMethod
Metric(::Fermionic, hilbert::Hilbert{<:Fock{:f}}) -> OperatorIndexToTuple
Metric(::Bosonic, hilbert::Hilbert{<:Fock{:b}}) -> OperatorIndexToTuple
Metric(::Phononic, hilbert::Hilbert{<:Phonon}) -> OperatorIndexToTuple

Get the index-to-tuple metric for a free fermionic/bosonic/phononic system.

source
TightBindingApproximation.BerryCurvatureType
BerryCurvature{B<:ReciprocalSpace, M<:BerryCurvatureMethod} <: Action

Berry curvature of energy bands.

Note

To obtain a rotation-symmetric Berry curvature, the :rcoordinate gauge should be used. Otherwise, artificial slight rotation symmetry breaking will occur.

source
TightBindingApproximation.BerryCurvatureMethod
BerryCurvature(reciprocalspace::ReciprocalSpace, μ::Real, d::Real=0.1, kx::T=nothing, ky::T=nothing) where {T<:Union{Nothing, Vector{Float64}}}

Construct a BerryCurvature using the Kubo method.

source
TightBindingApproximation.BerryCurvatureDataType
BerryCurvatureData{R<:ReciprocalSpace, C<:Array{Float64}, N<:Union{Vector{Float64}, Float64, Nothing}} <: Data

Data of Berry curvature, including:

  1. reciprocalspace::R: reciprocal space on which the Berry curvature is computed.
  2. values::C: Berry curvature for a certain set of bands.
  3. chernnumber::N: if not nothing, Chern number of each energy band when it is a vector or total Chern number of all bands when it is a number.
source
TightBindingApproximation.CompositeTBAType
CompositeTBA{
    K<:TBAKind,
    L<:Union{AbstractLattice, Nothing},
    S<:Union{OperatorSet{<:Operator}, Generator{<:OperatorSet{<:Operator}}},
    Q<:Quadraticization,
    H<:Union{OperatorSet{<:Quadratic}, Generator{<:OperatorSet{<:Quadratic}}},
    C<:Union{AbstractMatrix, Nothing}
} <:TBA{K, H, C}

Composite tight-binding approximation for quantum lattice systems.

source
TightBindingApproximation.DensityOfStatesType
DensityOfStates(brillouinzone::BrillouinZone, bands::Union{Colon, AbstractVector{Int}}=:, orbitals::Union{Colon, AbstractVector{Int}}...=:)
DensityOfStates(brillouinzone::BrillouinZone, bands::Union{Colon, AbstractVector{Int}}, orbitals::Tuple{Vararg{Union{Colon, AbstractVector{Int}}}})

Construct a DensityOfStates.

source
TightBindingApproximation.DensityOfStatesDataType
DensityOfStatesData <: Data

Data of density of states, including:

  1. energies::Vector{Float64}: energy sample points.
  2. values::Matrix{Float64}: density of states projected onto several certain sets of orbitals.
source
TightBindingApproximation.EnergyBandsType
EnergyBands(reciprocalspace::ReciprocalSpace, bands::Union{Colon, AbstractVector{Int}}=:, orbitals::AbstractVector{Int}...)
EnergyBands(reciprocalspace::ReciprocalSpace, bands::Union{Colon, AbstractVector{Int}}, orbitals::Tuple{Vararg{AbstractVector{Int}}})

Construct an EnergyBands.

source
TightBindingApproximation.EnergyBandsType
EnergyBands{L<:Tuple{Vararg{AbstractVector{Int}}}, R<:ReciprocalSpace, A<:Union{Colon, AbstractVector{Int}}} <: Action

Energy bands by tight-binding-approximation for quantum lattice systems.

source
TightBindingApproximation.EnergyBandsDataType
EnergyBandsData{R<:ReciprocalSpace, W<:Union{Array{Float64, 3}, Nothing}} <: Data

Data of energy bands, including:

  1. reciprocalspace::R: reciprocal space on which the energy bands are computed.
  2. values::Matrix{Float64}: eigen energies of bands with each column storing the values on the same reciprocal point.
  3. weights::W: if not nothing, weights of several certain sets of orbitals projected onto the energy bands.
source
TightBindingApproximation.FermiSurfaceType
FermiSurface{L<:Tuple{Vararg{AbstractVector{Int}}}, B<:Union{BrillouinZone, ReciprocalZone}, A<:Union{Colon, AbstractVector{Int}}} <: Action

Fermi surface of a free fermionic system.

source
TightBindingApproximation.FermiSurfaceType
FermiSurface(reciprocalspace::Union{BrillouinZone, ReciprocalZone}, μ::Real=0.0, bands::Union{Colon, AbstractVector{Int}}=:, orbitals::AbstractVector{Int}...)
FermiSurface(reciprocalspace::Union{BrillouinZone, ReciprocalZone}, μ::Real, bands::Union{Colon, AbstractVector{Int}}, orbitals::Tuple{Vararg{AbstractVector{Int}}})

Construct a FermiSurface.

source
TightBindingApproximation.FermiSurfaceDataType
FermiSurfaceData{S<:ReciprocalScatter} <: Data

Data of Fermi surface, including:

  1. values::S: points of the Fermi surface.
  2. weights::Matrix{Float64}: weights of several certain sets of orbitals projected onto the Fermi surface.
source
TightBindingApproximation.FukuiType
Fukui <: BerryCurvatureMethod

Fukui method to calculate Berry curvature of energy bands. see Fukui et al, JPSJ 74, 1674 (2005). Hall conductivity (single band) is given by

\[\sigma_{xy} = -{e^2\over h}\sum_n c_n, \nonumber \\ c_n = {1\over 2\pi}\int_{BZ}{dk_x dk_y Ω_{xy}}, \Omega_{xy}=(\nabla\times {\bm A})_z, A_{x}=i\langle u_n|\partial_x|u_n\rangle.\]

source
TightBindingApproximation.InelasticNeutronScatteringSpectraDataType
InelasticNeutronScatteringSpectraData{R<:ReciprocalSpace} <: Data

Data of inelastic neutron scattering spectra, including:

  1. reciprocalspace::R: reciprocal space on which the inelastic neutron scattering spectra are computed
  2. energies::Vector{Float64}: energy sample points.
  3. values::Matrix{Float64}: rescaled intensity of inelastic neutron scattering spectra.
source
TightBindingApproximation.KuboType
Kubo{K<:Union{Nothing, Vector{Float64}}} <: BerryCurvatureMethod

Kubo method to calculate the total Berry curvature of occupied energy bands. The Kubo formula is given by

\[\Omega_{ij}(\bm k)=\epsilon_{ijl}\sum_{n}f(\epsilon_n(\bm k))b_n^l(\bm k)=-2{\rm Im}\sum_v\sum_c{V_{vc,i}(\bm k)V_{cv,j}(\bm k)\over [\omega_c(\bm k)-\omega_v(\bm k)]^2},\]

where

\[ V_{cv,j}={\langle u_{c\bm k}|{\partial H\over \partial {\bm k}_j}|u_{v\bm k}\rangle}\]

v and c subscripts denote valence (occupied) and conduction (unoccupied) bands, respectively. Hall conductivity in 2D space is given by

\[\sigma_{xy}=-{e^2\over h}\int_{BZ}{dk_x dk_y\over 2\pi}{\Omega_{xy}}\]

source
TightBindingApproximation.SimpleTBAType
SimpleTBA{
    K<:TBAKind,
    L<:Union{AbstractLattice, Nothing},
    H<:Union{Formula, OperatorSet{<:Quadratic}, Generator{<:OperatorSet{<:Quadratic}}},
    C<:Union{AbstractMatrix, Nothing}
} <:TBA{K, H, C}

Simple tight-binding approximation for quantum lattice systems.

source
TightBindingApproximation.TBAType
TBA{K<:TBAKind, H<:Union{Formula, OperatorSet, Generator, Frontend}, C<:Union{Nothing, AbstractMatrix}} <: Frontend

Abstract type for free quantum lattice systems using the tight-binding approximation.

source
TightBindingApproximation.TBAMethod
TBA{K}(H::Union{Formula, OperatorSet{<:Quadratic}, Generator{<:OperatorSet{<:Quadratic}}}, commutator::Union{AbstractMatrix, Nothing}=nothing) where {K<:TBAKind}
TBA{K}(lattice::Union{AbstractLattice, Nothing}, H::Union{Formula, OperatorSet{<:Quadratic}, Generator{<:OperatorSet{<:Quadratic}}}, commutator::Union{AbstractMatrix, Nothing}=nothing) where {K<:TBAKind}

TBA{K}(H::Union{OperatorSet{<:Operator}, Generator{<:OperatorSet{<:Operator}}}, q::Quadraticization, commutator::Union{AbstractMatrix, Nothing}=nothing) where {K<:TBAKind}
TBA{K}(lattice::Union{AbstractLattice, Nothing}, H::Union{OperatorSet{<:Operator}, Generator{<:OperatorSet{<:Operator}}}, q::Quadraticization, commutator::Union{AbstractMatrix, Nothing}=nothing) where {K<:TBAKind}

TBA(lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, boundary::Boundary=plain; neighbors::Union{Int, Neighbors}=nneighbor(terms))
TBA{K}(lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, boundary::Boundary=plain; neighbors::Union{Int, Neighbors}=nneighbor(terms)) where {K<:TBAKind}

Construct a tight-binding quantum lattice system.

source
TightBindingApproximation.TBAMatrixType
TBAMatrix{C<:Union{AbstractMatrix, Nothing}, T, H<:Hermitian{T, <:AbstractMatrix{T}}} <: AbstractMatrix{T}

Matrix representation of a free quantum lattice system using the tight-binding approximation.

source
LinearAlgebra.eigenFunction
eigen(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; options...) -> Eigen
eigen(tba::Union{TBA, Algorithm{<:TBA}}, reciprocalspace::ReciprocalSpace; options...) -> Tuple{Vector{<:Vector{<:Number}}, Vector{<:Matrix{<:Number}}}

Get the eigen values and eigen vectors of a free quantum lattice system.

source
LinearAlgebra.eigvalsFunction
eigvals(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; options...) -> Vector{<:Number}
eigvals(tba::Union{TBA, Algorithm{<:TBA}}, reciprocalspace::ReciprocalSpace; options...) -> Vector{<:Vector{<:Number}}

Get the eigen values of a free quantum lattice system.

source
LinearAlgebra.eigvecsFunction
eigvecs(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; options...) -> Matrix{<:Number}
eigvecs(tba::Union{TBA, Algorithm{<:TBA}}, reciprocalspace::ReciprocalSpace; options...) -> Vector{<:Matrix{<:Number}}

Get the eigen vectors of a free quantum lattice system.

source
QuantumLattices.QuantumOperators.matrixFunction
matrix(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; gauge=:icoordinate, infinitesimal=infinitesimal(kind(tba))) -> TBAMatrix

Get the matrix representation of a free quantum lattice system.

source
QuantumLattices.QuantumOperators.scalartypeMethod
scalartype(::Type{<:TBA{<:TBAKind, H}}) where {H<:Union{Formula, OperatorSet, Generator, Frontend}}
scalartype(::Type{<:Algorithm{F}}) where {F<:TBA}

Get the scalartype of a tight-binding system.

source
QuantumLattices.add!Method
add!(dest::AbstractMatrix, mr::TBAMatrixization, m::Quadratic; infinitesimal=0, kwargs...) -> typeof(dest)

Matrixize a quadratic form and add it to dest.

source
QuantumLattices.add!Method
add!(dest::OperatorSum, q::Quadraticization{<:TBAKind{:TBA}}, m::Operator{<:Number, <:ID{CoordinatedIndex{<:Index, 2}}; kwargs...) -> typeof(dest)
add!(dest::OperatorSum, q::Quadraticization{<:TBAKind{:BdG}}, m::Operator{<:Number, <:ID{CoordinatedIndex{<:Index{<:FockIndex}}, 2}}; kwargs...) -> typeof(dest)
add!(dest::OperatorSum, q::Quadraticization{<:TBAKind{:BdG}}, m::Operator{<:Number, <:ID{CoordinatedIndex{<:Index{<:PhononIndex}}, 2}}; kwargs...) -> typeof(dest)

Get the unified quadratic form of a rank-2 operator and add it to dest.

source
QuantumLattices.dimensionMethod
dimension(tba::Union{TBA, Algorithm{<:TBA}}) -> Int

Get the dimension of the matrix representation of a free quantum lattice system.

source
QuantumLattices.kindMethod
kind(tba::Union{TBA, Algorithm{<:TBA}}) -> TBAKind
kind(::Type{<:TBA{K}}) where K -> TBAKind
kind(::Type{<:Algorithm{F}}) where {F<:TBA} -> TBAKind

Get the kind of a tight-binding system.

source
TightBindingApproximation.berrycurvatureMethod
berrycurvature(vectors::Matrix{<:Matrix{<:Number}}, commutator::Union{Nothing, AbstractMatrix{<:Number}}, dS::Number, ::Val{true}) -> Array{Float64, 3}
berrycurvature(vectors::Matrix{<:Matrix{<:Number}}, commutator::Union{Nothing, AbstractMatrix{<:Number}}, dS::Number, ::Val{false}) -> Matrix{Float64}

Based on the eigen vectors over a discrete reciprocal space, get the 1) Berry curvature for individual bands and 2) total Berry curvature by the Fukui method.

Here, commutator stands for the commutation relation matrix of the basis operators, dS stands for the area element of the discrete reciprocal space.

source
TightBindingApproximation.berrycurvatureMethod
berrycurvature(tba::TBA, bc::BerryCurvature{<:Union{BrillouinZone, ReciprocalZone}, <:Fukui{true}}; options...) -> Array{Float64, 3}
berrycurvature(tba::TBA, bc::BerryCurvature{<:Union{BrillouinZone, ReciprocalZone}, <:Fukui{false}}; options...) -> Matrix{Float64}

Get the Berry curvature by the Fukui method.

source
TightBindingApproximation.commutatorMethod
commutator(k::TBAKind, hilbert::Hilbert{<:Internal}) -> Union{AbstractMatrix, Nothing}

Get the commutation relation of the single-particle operators of a free quantum lattice system using the tight-binding approximation.

source
TightBindingApproximation.infinitesimalMethod
infinitesimal(::TBAKind{:TBA}) -> 0
infinitesimal(::TBAKind{:BdG}) -> atol/5
infinitesimal(::Fermionic{:BdG}) -> 0

Infinitesimal used in the matrixization of a tight-binding approximation to be able to capture the Goldstone mode.

source

TightBindingApproximation.Fitting

TightBindingApproximation.Fitting.SampleNodeType
SampleNode(reciprocals::AbstractVector{<:AbstractVector{<:Number}}, position::Vector{<:Number}, bands::AbstractVector{Int}, values::Vector{<:Number}, ratio::Number)
SampleNode(reciprocals::AbstractVector{<:AbstractVector{<:Number}}, position::Vector{<:Number}, bands::AbstractVector{Int}, values::Vector{<:Number}, ratios::Vector{<:Number}=ones(length(bands)))

A sample node of a momentum-eigenvalues pair.

source
TightBindingApproximation.Fitting.deviationMethod
deviation(tba::Union{TBA, Algorithm{<:TBA}}, samplenode::SampleNode) -> Float64
deviation(tba::Union{TBA, Algorithm{<:TBA}}, samplesets::Vector{SampleNode}) -> Float64

Get the deviation of the eigenvalues between the sample points and model points.

source
TightBindingApproximation.Fitting.optimize!Function
optimize!(
    tba::Union{TBA, Algorithm{<:TBA}}, samplesets::Vector{SampleNode}, variables=keys(Parameters(tba));
    verbose=false, method=LBFGS(), options=Options()
) -> Tuple{typeof(tba), Optim.MultivariateOptimizationResults}

Optimize the parameters of a tight binding system whose names are specified by variables so that the total deviations of the eigenvalues between the model points and sample points minimize.

source

TightBindingApproximation.Wannier90

TightBindingApproximation.Wannier90.W90Type
W90(path::AbstractString, prefix::AbstractString="wannier90"; name::Symbol=Symbol(prefix), projection::Bool=true)

Construct a quantum lattice system based on the files obtained from Wannier90.

source
TightBindingApproximation.Wannier90.W90Method
W90(lattice::Lattice, centers::Matrix{<:Real}, H::OperatorSum{W90Hoppings})
W90(lattice::Lattice, hilbert::Hilbert, H::OperatorSum{W90Hoppings})

Construct a quantum lattice system based on the information obtained from Wannier90.

In general, the Wannier centers could deviate from their corresponding atom positions. When centers::Matrix{<:Real} is used, the the Wannier centers are assigned directly. When hilbert::Hilbert is used, the Wannier centers will be approximated by their corresponding atom positions.

source
TightBindingApproximation.Wannier90.W90HoppingsType
W90Hoppings <: OperatorPack{Matrix{Float64}, NTuple{3, Int}}

Hopping amplitudes among the Wannier orbitals in two unitcells with a fixed relative displacement.

Here,

  1. the :value attribute represents the hopping amplitude matrix of Wannier orbitals,
  2. the :id attribute represents the relative displacement of the two unitcells in basis of the lattice translation vectors.
source
QuantumLattices.QuantumOperators.matrixFunction
matrix(wan::W90, k::AbstractVector{<:Number}=SVector(0, 0, 0); gauge=:icoordinate, kwargs...) -> TBAMatrix

Get the matrix representation of a quantum lattice system based on the information obtained from Wannier90.

source
QuantumLattices.add!Method
add!(dest::AbstractMatrix, mr::W90Matrixization, m::W90Hoppings; kwargs...) -> typeof(dest)

Matrixize a group of Wannier hoppings and add it to dest.

source
TightBindingApproximation.Wannier90.readlatticeFunction
readlattice(path::AbstractString, prefix::AbstractString="wannier90"; name::Symbol=Symbol(prefix), projection::Bool=true) -> Lattice{3, Float64, 3}

Read the lattice from the ".win" configuration file of Wannier90 with a given name.

Besides, projection specifies whether only the sites with initial projections in the ".win" file are included in the constructed lattice.

source