Core of TightBindingApproximation
TightBindingApproximation.basicoptions
— Constantconst basicoptions = Dict(
: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.
QuantumLattices.DegreesOfFreedom.Metric
— MethodMetric(::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.
QuantumLattices.DegreesOfFreedom.Table
— MethodTable(hilbert::Hilbert{Phonon{:}}, by::OperatorIndexToTuple{(kind, :site, :direction)})
Construct a index-sequence table for a phononic system.
TightBindingApproximation.BerryCurvature
— TypeBerryCurvature{B<:ReciprocalSpace, M<:BerryCurvatureMethod, O} <: Action
Berry curvature of energy bands.
To obtain a rotation-symmetric Berry curvature, the :rcoordinate
gauge should be used. Otherwise, artificial slight rotation symmetry breaking will occur.
TightBindingApproximation.BerryCurvature
— MethodBerryCurvature(reciprocalspace::ReciprocalSpace, method::BerryCurvatureMethod; options...)
BerryCurvature(reciprocalspace::ReciprocalSpace, μ::Real, d::Real=0.1, kx::T=nothing, ky::T=nothing; gauge=:rcoordinate, options...) where {T<:Union{Nothing, Vector{Float64}}}
BerryCurvature(reciprocalspace::BrillouinZone, bands::AbstractVector{Int}, abelian::Bool=true; gauge=:icoordinate, options...)
BerryCurvature(reciprocalspace::ReciprocalZone, bands::AbstractVector{Int}, abelian::Bool=true; gauge=:rcoordinate, options...)
Construct a BerryCurvature
.
TightBindingApproximation.BerryCurvatureData
— TypeBerryCurvatureData{R<:ReciprocalSpace, C<:Array{Float64}, N<:Union{Vector{Float64}, Float64, Nothing}} <: Data
Data of Berry curvature, including:
reciprocalspace::R
: reciprocal space on which the Berry curvature is computed.values::C
: Berry curvature for a certain set of bands.chernnumber::N
: if notnothing
, Chern number of each energy band when it is a vector or total Chern number of all bands when it is a number.
TightBindingApproximation.BerryCurvatureMethod
— Typeabstract type BerryCurvatureMethod end
Abstract type for calculation of Berry curvature.
TightBindingApproximation.Bosonic
— TypeBosonic{K} <: TBAKind{K}
Bosonic quantum lattice system.
TightBindingApproximation.CompositeTBA
— TypeCompositeTBA{
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.
TightBindingApproximation.DensityOfStates
— TypeDensityOfStates(brillouinzone::BrillouinZone, bands::Union{Colon, AbstractVector{Int}}=:, orbitals::Union{Colon, AbstractVector{Int}}...=:; options...)
Construct a DensityOfStates
.
TightBindingApproximation.DensityOfStates
— TypeDensityOfStates{B<:BrillouinZone, A<:Union{Colon, AbstractVector{Int}}, L<:Tuple{Vararg{Union{Colon, AbstractVector{Int}}}}, O} <: Action
Density of states of a tight-binding system.
TightBindingApproximation.DensityOfStatesData
— TypeDensityOfStatesData <: Data
Data of density of states, including:
energies::Vector{Float64}
: energy sample points.values::Matrix{Float64}
: density of states projected onto several certain sets of orbitals.
TightBindingApproximation.EnergyBands
— TypeEnergyBands(reciprocalspace::ReciprocalSpace, bands::Union{Colon, AbstractVector{Int}}=:, orbitals::AbstractVector{Int}...; options...)
Construct an EnergyBands
.
TightBindingApproximation.EnergyBands
— TypeEnergyBands{L<:Tuple{Vararg{AbstractVector{Int}}}, R<:ReciprocalSpace, A<:Union{Colon, AbstractVector{Int}}, O} <: Action
Energy bands by tight-binding-approximation for quantum lattice systems.
TightBindingApproximation.EnergyBandsData
— TypeEnergyBandsData{R<:ReciprocalSpace, W<:Union{Vector{Matrix{Float64}}, Nothing}} <: Data
Data of energy bands, including:
reciprocalspace::R
: reciprocal space on which the energy bands are computed.values::Matrix{Float64}
: eigen energies of bands with each column storing the values on the same reciprocal point.weights::W
: if notnothing
, weights of several certain sets of orbitals projected onto the energy bands.
TightBindingApproximation.FermiSurface
— TypeFermiSurface(reciprocalspace::Union{BrillouinZone, ReciprocalZone}, μ::Real=0.0, bands::Union{Colon, AbstractVector{Int}}=:, orbitals::AbstractVector{Int}...; options...)
Construct a FermiSurface
.
TightBindingApproximation.FermiSurface
— TypeFermiSurface{L<:Tuple{Vararg{AbstractVector{Int}}}, B<:Union{BrillouinZone, ReciprocalZone}, A<:Union{Colon, AbstractVector{Int}}, O} <: Action
Fermi surface of a free fermionic system.
TightBindingApproximation.FermiSurfaceData
— TypeFermiSurfaceData{S<:ReciprocalScatter} <: Data
Data of Fermi surface, including:
values::S
: points of the Fermi surface.weights::Vector{Vector{Float64}}
: weights of several certain sets of orbitals projected onto the Fermi surface.
TightBindingApproximation.Fermionic
— TypeFermionic{K} <: TBAKind{K}
Fermionic quantum lattice system.
TightBindingApproximation.Fukui
— TypeFukui <: 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.\]
TightBindingApproximation.InelasticNeutronScatteringSpectra
— TypeInelasticNeutronScatteringSpectra{R<:ReciprocalSpace, O} <: Action
Inelastic neutron scattering spectra.
TightBindingApproximation.InelasticNeutronScatteringSpectra
— MethodInelasticNeutronScatteringSpectra(reciprocalspace::ReciprocalSpace, energies::AbstractVector{<:Number}; options...)
Construct an InelasticNeutronScatteringSpectra
.
TightBindingApproximation.InelasticNeutronScatteringSpectraData
— TypeInelasticNeutronScatteringSpectraData{R<:ReciprocalSpace} <: Data
Data of inelastic neutron scattering spectra, including:
reciprocalspace::R
: reciprocal space on which the inelastic neutron scattering spectra are computedenergies::Vector{Float64}
: energy sample points.values::Matrix{Float64}
: rescaled intensity of inelastic neutron scattering spectra.
TightBindingApproximation.Kubo
— TypeKubo{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}}\]
TightBindingApproximation.Phononic
— TypePhononic <: TBAKind{:BdG}
Phononic quantum lattice system.
TightBindingApproximation.Quadratic
— TypeQuadratic{V<:Number, C<:AbstractVector} <: OperatorPack{V, Tuple{Tuple{Int, Int}, C, C}}
The unified quadratic form in tight-binding approximation.
TightBindingApproximation.Quadraticization
— TypeQuadraticization{K<:TBAKind, T<:Table} <: LinearTransformation
The linear transformation that converts a rank-2 operator to its unified quadratic form.
TightBindingApproximation.SimpleTBA
— TypeSimpleTBA{
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.
TightBindingApproximation.TBA
— TypeTBA{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.
TightBindingApproximation.TBA
— MethodTBA{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.
TightBindingApproximation.TBAKind
— TypeTBAKind{K}
The kind of a free quantum lattice system using the tight-binding approximation.
TightBindingApproximation.TBAKind
— MethodTBAKind(T::Type{<:Term}, I::Type{<:Internal})
Depending on the kind of a Term
type and an Internal
type, get the corresponding TBA kind.
TightBindingApproximation.TBAMatrix
— TypeTBAMatrix{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.
TightBindingApproximation.TBAMatrixization
— TypeTBAMatrixization{D<:Number, V<:Union{AbstractVector{<:Number}, Nothing}} <: Matrixization
Matrixization of the Hamiltonian of a tight-binding system.
LinearAlgebra.eigen
— Functioneigen(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; kwargs...) -> Eigen
eigen(tba::Union{TBA, Algorithm{<:TBA}}, reciprocalspace::ReciprocalSpace; kwargs...) -> Tuple{Vector{<:Vector{<:Number}}, Vector{<:Matrix{<:Number}}}
Get the eigen values and eigen vectors of a free quantum lattice system.
LinearAlgebra.eigen
— Methodeigen(m::TBAMatrix) -> Eigen
Solve the eigen problem of a free quantum lattice system.
LinearAlgebra.eigvals
— Functioneigvals(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; kwargs...) -> Vector{<:Number}
eigvals(tba::Union{TBA, Algorithm{<:TBA}}, reciprocalspace::ReciprocalSpace; kwargs...) -> Vector{<:Vector{<:Number}}
Get the eigen values of a free quantum lattice system.
LinearAlgebra.eigvals
— Methodeigvals(m::TBAMatrix) -> Vector
Get the eigen values of a free quantum lattice system.
LinearAlgebra.eigvecs
— Functioneigvecs(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; kwargs...) -> Matrix{<:Number}
eigvecs(tba::Union{TBA, Algorithm{<:TBA}}, reciprocalspace::ReciprocalSpace; kwargs...) -> Vector{<:Matrix{<:Number}}
Get the eigen vectors of a free quantum lattice system.
LinearAlgebra.eigvecs
— Methodeigvecs(m::TBAMatrix) -> Matrix
Get the eigen vectors of a free quantum lattice system.
QuantumLattices.QuantumOperators.matrix
— Functionmatrix(tba::Union{TBA, Algorithm{<:TBA}}, k::Union{AbstractVector{<:Number}, Nothing}=nothing; gauge=:icoordinate, infinitesimal=infinitesimal(kind(tba.frontend)), kwargs...) -> TBAMatrix
Get the matrix representation of a free quantum lattice system.
QuantumLattices.add!
— Methodadd!(dest::AbstractMatrix, mr::TBAMatrixization, m::Quadratic; infinitesimal=0, kwargs...) -> typeof(dest)
Matrixize a quadratic form and add it to dest
.
QuantumLattices.add!
— Methodadd!(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
.
QuantumLattices.dimension
— Methoddimension(tba::Union{TBA, Algorithm{<:TBA}}) -> Int
Get the dimension of the matrix representation of a free quantum lattice system.
TightBindingApproximation.commutator
— Methodcommutator(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.
TightBindingApproximation.infinitesimal
— Methodinfinitesimal(::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.
TightBindingApproximation.Fitting
TightBindingApproximation.Fitting.SampleNode
— TypeSampleNode(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.
TightBindingApproximation.Fitting.deviation
— Methoddeviation(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.
TightBindingApproximation.Fitting.optimize!
— Functionoptimize!(
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.
TightBindingApproximation.Wannier90
TightBindingApproximation.Wannier90.W90
— TypeW90(path::AbstractString, prefix::AbstractString="wannier90"; name::Symbol=Symbol(prefix), projection::Bool=true)
Construct a quantum lattice system based on the files obtained from Wannier90.
TightBindingApproximation.Wannier90.W90
— TypeW90 <: TBA{Fermionic{:TBA}, OperatorSum{W90Hoppings, NTuple{3, Int}}, Nothing}
A quantum lattice system based on the information obtained from Wannier90.
TightBindingApproximation.Wannier90.W90
— MethodW90(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.
TightBindingApproximation.Wannier90.W90Hoppings
— TypeW90Hoppings <: OperatorPack{Matrix{Float64}, NTuple{3, Int}}
Hopping amplitudes among the Wannier orbitals in two unitcells with a fixed relative displacement.
Here,
- the
:value
attribute represents the hopping amplitude matrix of Wannier orbitals, - the
:id
attribute represents the relative displacement of the two unitcells in basis of the lattice translation vectors.
TightBindingApproximation.Wannier90.W90Matrixization
— TypeW90Matrixization{V<:AbstractVector{<:Real}} <: Matrixization
Matrixization of the Hamiltonian obtained from Wannier90.
QuantumLattices.QuantumOperators.matrix
— Functionmatrix(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.
QuantumLattices.add!
— Methodadd!(dest::AbstractMatrix, mr::W90Matrixization, m::W90Hoppings; kwargs...) -> typeof(dest)
Matrixize a group of Wannier hoppings and add it to dest
.
TightBindingApproximation.Wannier90.findblock
— Methodfindblock(name::String, content::String) -> Union{Nothing, Vector{SubString{String}}}
Find a named block in the content of Wannier90 configuration files.
TightBindingApproximation.Wannier90.readbands
— Functionreadbands(path::AbstractString, prefix::AbstractString="wannier90") -> Tuple{Vector{Float64}, Matrix{Float64}}
Read the band structure from Wannier90 files by setting "bands_plot = .true.".
TightBindingApproximation.Wannier90.readcenters
— Functionreadcenters(path::AbstractString, prefix::AbstractString="wannier90") -> Matrix{Float64}
Read the centers of Wannier functions obtained by Wannier90.
TightBindingApproximation.Wannier90.readhamiltonian
— Functionreadhamiltonian(path::AbstractString, prefix::AbstractString="wannier90") -> OperatorSum{W90Hoppings, NTuple{3, Int}}
Read the hamiltonian from the "_hr.dat" file generated by Wannier90.
TightBindingApproximation.Wannier90.readlattice
— Functionreadlattice(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.