ExactDiagonalization
Quantum numbers
Quantum numbers can be considered as the conserved labels for the bases of a Hilbert space when a quantum system hosts some symmetries. Here we only implement Abelian quantum numbers because non-Abelian ones are far more complicated yet much less used.
AbelianQuantumNumber
RepresentationSpace
Core of ExactDiagonalization
BinaryBases
AbelianBases
Manual
ExactDiagonalization.basicoptions — Constant
const basicoptions = (
nev = "number of eigenvalues to be computed",
which = "type of eigenvalues to be computed",
tol = "tolerance of the computation",
maxiter = "maximum iteration of the computation",
v₀ = "initial state",
krylovdim = "maximum dimension of the Krylov subspace that will be constructed",
verbosity = "verbosity level"
)Basic options of actions for exact diagonalization method.
ExactDiagonalization.edtimer — Constant
const edtimer = TimerOutput()Default shared timer for all exact diagonalization methods.
ExactDiagonalization.Abelian — Type
const Abelian = AbelianQuantumNumberType alias for AbelianQuantumNumber.
ExactDiagonalization.Abelian — Method
Abelian(bs::AbelianBases)Get the Abelian quantum number of a set of spin bases.
ExactDiagonalization.Abelian — Method
Abelian(bs::BinaryBases)Get the Abelian quantum number of a set of binary bases.
ExactDiagonalization.AbelianBases — Type
AbelianBases{A<:Abelian, N} <: SectorA set of Abelian bases, that is, a set of bases composed from the product of local Abelian Graded spaces.
To improve the efficiency of the product of local Abelian Graded spaces, we adopt a two-step strategy:
- partition the local spaces into several groups in each of which the local spaces are direct producted and rearranged according to the Abelian quantum numbers, and then
- glue the results obtained in the previous step so that a sector with a certain Abelian quantum number can be targeted.
In principle, a binary-tree strategy can be more efficient, but our two-step strategy is enough for a quantum system that can be solved by the exact diagonalization method.
The partition of the local Abelian Graded spaces is assigned by a NTuple{N, Vector{Int}}, with each of its element contains the sequences of the grouped local spaces specified by a table.
ExactDiagonalization.AbelianBases — Method
AbelianBases(locals::AbstractVector{Int}, partition::NTuple{N, AbstractVector{Int}}=partition(length(locals))) where NConstruct a set of spin bases that subjects to no quantum number conservation.
ExactDiagonalization.AbelianBases — Method
AbelianBases(locals::Vector{Graded{A}}, quantumnumber::A, partition::NTuple{N, AbstractVector{Int}}=partition(length(locals))) where {N, A<:Abelian}Construct a set of spin bases that preserves a certain symmetry specified by the corresponding Abelian quantum number.
ExactDiagonalization.AbelianGradedSpace — Type
AbelianGradedSpace{QN<:AbelianQuantumNumber} <: RepresentationSpace{QN}A quantum representation space of an Abelian group that has been decomposed into the direct sum of its irreducible representations.
ExactDiagonalization.AbelianGradedSpace — Type
AbelianGradedSpace(quantumnumbers::AbstractVector{<:AbelianQuantumNumber}, dimensions::AbstractVector{<:Integer}, dual::Bool=false; ordercheck::Bool=false, duplicatecheck::Bool=false, degeneracycheck::Bool=false)Construct an Abelian graded space.
Here:
quantumnumbersspecifies the Abelian quantum numbers labeling the irreducible representations of the corresponding Abelian group which must be sorted in the ascending order. Such an ordering should be manually guaranteed by the user. When and only when the keyword argumentordercheckistrue, the constructor will check whether this condition is satisfied and raise an error if it doesn't. Besides,quantumnumbersmust not contain duplicate Abelian quantum numbers, manually guaranteed by the user as well. This condition can be checked when and only when bothordercheck==trueandduplicatecheck==true. An error will be raised if this check fails.dimensionsspecifies the degenerate dimensions of the corresponding Abelian quantum numbers. Apparently, each degenerate dimension must be positive, which should also be manually guaranteed by the user. When and only when the keyword argumentdegeneracycheckistrue, the constructor will check whether this condition is satisfied and raise an error if it doesn't.dualspecifies whether the graded space is the dual representation of the corresponding Abelian group, which roughly speaking, can be viewed as the direction of the arrow of the Abelian quantum numbers labeling the irreducible representations. We assumedual==truecorresponds to the in-arrow anddual==falsecorresponds to the out-arrow.
ExactDiagonalization.AbelianGradedSpace — Method
AbelianGradedSpace(pairs; dual::Bool=false)
AbelianGradedSpace(pairs::Pair...; dual::Bool=false)
AbelianGradedSpace{QN}(pairs; dual::Bool=false) where {QN<:AbelianQuantumNumber}
AbelianGradedSpace{QN}(pairs::Pair...; dual::Bool=false) where {QN<:AbelianQuantumNumber}Construct an Abelian graded space.
In this function, the Abelian quantum numbers will be sorted automatically, therefore, their orders need not be worried. Duplicate and dimension checks of the quantum numbers are also carried out and errors will be raised if either such checks fails.
ExactDiagonalization.AbelianGradedSpaceProd — Type
AbelianGradedSpaceProd{N, QN<:AbelianQuantumNumber} <: CompositeAbelianGradedSpace{N, QN}Direct product of Abelian graded spaces.
ExactDiagonalization.AbelianGradedSpaceSum — Type
AbelianGradedSpaceSum{N, QN<:AbelianQuantumNumber} <: CompositeAbelianGradedSpace{N, QN}Direct sum of Abelian graded spaces.
ExactDiagonalization.AbelianQuantumNumber — Type
AbelianQuantumNumberAbstract type of Abelian quantum numbers.
An Abelian quantum number is the label of a irreducible representation of an Abelian group acted on a quantum representation space.
ExactDiagonalization.AbelianQuantumNumberProd — Type
AbelianQuantumNumberProd{T<:Tuple{Vararg{SimpleAbelianQuantumNumber}}} <: AbelianQuantumNumberDeligne tensor product of simple Abelian quantum numbers.
ExactDiagonalization.AbelianQuantumNumberProd — Method
AbelianQuantumNumberProd(contents::SimpleAbelianQuantumNumber...)
AbelianQuantumNumberProd(contents::Tuple{Vararg{SimpleAbelianQuantumNumber}})Construct a Deligne tensor product of simple Abelian quantum numbers.
ExactDiagonalization.AbelianQuantumNumberProd — Method
AbelianQuantumNumberProd{T}(vs::Vararg{Number, N}) where {N, T<:NTuple{N, SimpleAbelianQuantumNumber}}
AbelianQuantumNumberProd{T}(vs::NTuple{N, Number}) where {N, T<:NTuple{N, SimpleAbelianQuantumNumber}}Construct a Deligne tensor product of simple Abelian quantum numbers by their values.
ExactDiagonalization.AbstractGreenFunction — Type
AbstractGreenFunction{T<:Number} <: FunctionAbstract type for Green's functions obtained by Krylov space expansion.
ExactDiagonalization.AbstractGreenFunction — Method
(gf::AbstractGreenFunction)(ω::Number) -> Matrix{promote_type(typeof(ω), eltype(gf))}Get the values of an AbstractGreenFunction at ω.
ExactDiagonalization.BinaryBases — Type
BinaryBases{A<:Abelian, B<:BinaryBasis, T<:AbstractVector{B}} <: SectorA set of binary bases.
ExactDiagonalization.BinaryBases — Method
BinaryBases(spindws, spinups, spinfulparticle::Abelian[ℕ ⊠ 𝕊ᶻ])
BinaryBases(spindws, spinups, spinfulparticle::Abelian[𝕊ᶻ ⊠ ℕ])Construct a set of binary bases that preserves both the particle number and the spin z-component conservation.
ExactDiagonalization.BinaryBases — Method
BinaryBases(spindws, spinups, sz::𝕊ᶻ)Construct a set of binary bases that preserves the spin z-component but not the particle number conservation.
ExactDiagonalization.BinaryBases — Method
BinaryBases(states, particle::ℕ)
BinaryBases(nstate::Integer, particle::ℕ)Construct a set of binary bases that preserves the particle number conservation.
ExactDiagonalization.BinaryBases — Method
BinaryBases(states)
BinaryBases(nstate::Integer)Construct a set of binary bases that subjects to no quantum number conservation.
ExactDiagonalization.BinaryBasis — Type
BinaryBasis{I<:Unsigned}Binary basis represented by an unsigned integer.
Here, we adopt the following common rules:
- In a binary basis, a bit of an unsigned integer represents a single-particle state that can be occupied (
1) or unoccupied (0). - The position of this bit in the unsigned integer counting from the right corresponds to the sequence of the single-particle state specified by a table.
- When representing a many-body state by creation operators, they are arranged in ascending order according to their sequences.
In this way, any many-body state of canonical fermionic or hardcore bosonic systems can be represented ambiguously by the binary bases, e.g., $c^†_2c^†_3c^†_4|\text{Vacuum}\rangle$ is represented by $1110$.
ExactDiagonalization.BinaryBasis — Method
BinaryBasis(states; filter=index->true)
BinaryBasis{I}(states; filter=index->true) where {I<:Unsigned}Construct a binary basis with the given occupied states.
ExactDiagonalization.BinaryBasisRange — Type
BinaryBasisRange{I<:Unsigned} <: VectorSpace{BinaryBasis{I}}A continuous range of binary basis from 0 to 2^n-1.
ExactDiagonalization.CompositeAbelianGradedSpace — Type
CompositeAbelianGradedSpace{N, QN<:AbelianQuantumNumber} <: RepresentationSpace{QN}Abstract type of composite Abelian graded spaces.
ExactDiagonalization.ED — Type
ED(
lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, boundary::Boundary=plain, dtype::Type{<:Number}=valtype(terms);
neighbors::Union{Int, Neighbors}=nneighbor(terms)
)
ED(
lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, table::AbstractDict, sectors::OneOrMore{Sector}=Sector(hilbert; table=table), boundary::Boundary=plain, dtype::Type{<:Number}=valtype(terms);
neighbors::Union{Int, Neighbors}=nneighbor(terms)
)Construct the exact diagonalization method for a quantum lattice system.
ExactDiagonalization.ED — Type
ED(
lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, quantumnumbers::OneOrMore{Abelian}, boundary::Boundary=plain, dtype::Type{<:Number}=valtype(terms);
neighbors::Union{Int, Neighbors}=nneighbor(terms)
)
ED(
lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, table::AbstractDict, quantumnumbers::OneOrMore{Abelian}, boundary::Boundary=plain, dtype::Type{<:Number}=valtype(terms);
neighbors::Union{Int, Neighbors}=nneighbor(terms)
)Construct the exact diagonalization method for a quantum lattice system.
ExactDiagonalization.ED — Type
ED(system::Generator{<:Operators}, table::AbstractDict, sectors::OneOrMore{Sector}, dtype::Type{<:Number}=scalartype(system))
ED(lattice::Union{AbstractLattice, Nothing}, system::Generator{<:Operators}, table::AbstractDict, sectors::OneOrMore{Sector}, dtype::Type{<:Number}=scalartype(system))Construct the exact diagonalization method for a quantum lattice system.
ExactDiagonalization.ED — Type
ED{K<:EDKind, L<:Union{AbstractLattice, Nothing}, S<:Generator{<:Operators}, M<:EDMatrixization, H<:CategorizedGenerator{<:OperatorSum{<:EDMatrix}}} <: FrontendExact diagonalization method of a quantum lattice system.
ExactDiagonalization.EDEigen — Type
EDEigen{S<:Tuple{Vararg{Union{Abelian, Sector}}}} <: ActionEigen system by exact diagonalization method.
ExactDiagonalization.EDEigen — Method
EDEigen(sectors::Union{Abelian, Sector}...)
EDEigen(sectors:::Tuple{Vararg{Union{Abelian, Sector}}})Construct an EDEigen.
ExactDiagonalization.EDEigenData — Type
EDEigenData{V<:Number, T<:Number, S<:Sector} <: DataEigen decomposition in exact diagonalization method.
Compared to the usual eigen decomposition Eigen, EDEigenData contains a :sectors attribute to store the sectors of Hilbert space in which the eigen values and eigen vectors are computed. Furthermore, given that in different sectors the dimensions of the sub-Hilbert spaces can also be different, the :vectors attribute of EDEigenData is a vector of vector instead of a matrix.
ExactDiagonalization.EDKind — Type
EDKind{K}Kind of the exact diagonalization method applied to a quantum lattice system.
ExactDiagonalization.EDKind — Method
EDKind(::Type{<:FockIndex})Kind of the exact diagonalization method applied to a canonical quantum Fock lattice system.
ExactDiagonalization.EDKind — Method
EDKind(::Type{<:SpinIndex})Kind of the exact diagonalization method applied to a canonical quantum spin lattice system.
ExactDiagonalization.EDMatrix — Type
EDMatrix{M<:SparseMatrixCSC, S<:Sector} <: OperatorPack{M, Tuple{S, S}}Matrix representation of quantum operators between a ket Hilbert space and a bra Hilbert space.
ExactDiagonalization.EDMatrix — Method
EDMatrix(m::SparseMatrixCSC, sector::Sector)
EDMatrix(m::SparseMatrixCSC, braket::NTuple{2, Sector})
EDMatrix(m::SparseMatrixCSC, bra::Sector, ket::Sector)Construct a matrix representation when
- the bra and ket Hilbert spaces share the same bases;
- the bra and ket Hilbert spaces may be different;
- the bra and ket Hilbert spaces may or may not be the same.
ExactDiagonalization.EDMatrixization — Type
EDMatrixization{D<:Number, T<:AbstractDict, S<:Sector} <: MatrixizationMatrixization of a quantum lattice system on a target Hilbert space.
ExactDiagonalization.EDMatrixization — Method
EDMatrixization{D}(table::AbstractDict, sector::S, sectors::S...) where {D<:Number, S<:Sector}
EDMatrixization{D}(table::AbstractDict, brakets::Vector{Tuple{S, S}}) where {D<:Number, S<:Sector}Construct a matrixization.
ExactDiagonalization.Graded — Type
const Graded = AbelianGradedSpaceType alias for AbelianGradedSpace.
ExactDiagonalization.Graded — Method
Graded{ℤ₁}(spin::Spin)
Graded{𝕊ᶻ}(spin::Spin)Decompose a local spin space into an Abelian graded space that preserves 1) no symmetry, and 2) spin-z component symmetry.
ExactDiagonalization.GreenFunction — Type
GreenFunction(
operators::AbstractVector{<:QuantumOperator}, ed::Algorithm{<:ED}, kind::Symbol=:greater;
E₀::Union{Real, Nothing}=nothing, Ω::Union{AbstractVector{<:Number}, Nothing}=nothing, sector₀::Union{Sector, Nothing}=nothing, maxdim::Integer=200, tol::Real=1e-12, kwargs...
)
GreenFunction(
operators::AbstractVector{<:QuantumOperator}, ed::ED, kind::Symbol=:greater;
E₀::Union{Real, Nothing}=nothing, Ω::Union{AbstractVector{<:Number}, Nothing}=nothing, sector₀::Union{Sector, Nothing}=nothing, maxdim::Integer=200, tol::Real=1e-12, timer::TimerOutput=edtimer, kwargs...
)Construct a GreenFunction.
ExactDiagonalization.GreenFunction — Type
GreenFunction{T<:Number, V<:Real} <: AbstractGreenFunction{T}Green function obtained by Krylov space expansion.
ExactDiagonalization.GreenFunction — Method
(gf::GreenFunction)(dest::AbstractMatrix{<:Number}, ω::Number; sign::Bool=false) -> typeof(dest)Get the values of a GreenFunction at ω and add the result to dest (when sign is false) or substrate the result from dest (when sign is true).
ExactDiagonalization.GroundStateExpectation — Type
GroundStateExpectation{D<:Number, O<:Array{<:Union{Operator, Operators}}} <: ActionGround state expectation of operators.
ExactDiagonalization.GroundStateExpectationData — Type
GroundStateExpectationData{A<:Array{<:Number}} <: DataData of ground state expectation of operators, including:
values::A: values of the ground state expectation.
ExactDiagonalization.RepresentationSpace — Type
RepresentationSpace{QN<:AbelianQuantumNumber} <: VectorSpace{QN}Abstract type of quantum representation spaces of Abelian groups.
ExactDiagonalization.RetardedGreenFunction — Type
RetardedGreenFunction(operators::AbstractVector{<:QuantumOperator}, ed::Algorithm{<:ED}, sign::Bool=false; maxdim::Integer=200, tol::Real=1e-12, kwargs...)
RetardedGreenFunction(operators::AbstractVector{<:QuantumOperator}, ed::ED, sign::Bool=false; maxdim::Integer=200, tol::Real=1e-12, timer::TimerOutput=edtimer, kwargs...)Construct a RetardedGreenFunction.
ExactDiagonalization.RetardedGreenFunction — Type
RetardedGreenFunction{T<:Number, V<:Real} <: AbstractGreenFunction{T}Retarded Green's function obtained by Krylov space expansion.
ExactDiagonalization.RetardedGreenFunction — Method
(gf::RetardedGreenFunction)(dest::AbstractMatrix{<:Number}, ω::Number) -> typeof(dest)Get the values of a RetardedGreenFunction at ω and add the result to dest.
ExactDiagonalization.Sector — Type
Sector(hilbert::Hilbert{<:Fock}, basistype::Type{<:Unsigned}=UInt) -> BinaryBases
Sector(quantumnumber::ℕ, hilbert::Hilbert{<:Fock}, basistype::Type{<:Unsigned}=UInt; table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> BinaryBases
Sector(quantumnumber::Union{𝕊ᶻ, Abelian[ℕ ⊠ 𝕊ᶻ], Abelian[𝕊ᶻ ⊠ ℕ]}, hilbert::Hilbert{<:Fock}, basistype::Type{<:Unsigned}=UInt; table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> BinaryBasesConstruct the binary bases of a Hilbert space with the specified quantum number.
ExactDiagonalization.Sector — Type
Sector(hilbert::Hilbert{<:Spin}, partition::Tuple{Vararg{AbstractVector{Int}}}=partition(length(hilbert))) -> AbelianBases
Sector(
quantumnumber::Abelian, hilbert::Hilbert{<:Spin}, partition::Tuple{Vararg{AbstractVector{Int}}}=partition(length(hilbert));
table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))
) -> AbelianBasesConstruct the Abelian bases of a spin Hilbert space with the specified quantum number.
ExactDiagonalization.Sector — Type
SectorA sector of the Hilbert space which forms the bases of an irreducible representation of the Hamiltonian of a quantum lattice system.
ExactDiagonalization.SectorFilter — Type
SectorFilter{S} <: LinearTransformationFilter the target bra and ket Hilbert spaces.
ExactDiagonalization.SimpleAbelianQuantumNumber — Type
SimpleAbelianQuantumNumber <: AbelianQuantumNumberAbstract type of simple Abelian quantum numbers. That is, it contains only one label.
ExactDiagonalization.SpinCoherentState — Type
SpinCoherentState <: CompositeDict{Int, Tuple{Float64, Float64}}Spin coherent state on a block of lattice sites.
The structure of the spin coherent state is specified by a Dict{Int, Tuple{Float64, Float64}}, which contains the site-(θ, ϕ) pairs with site being the site index in a lattice and (θ, ϕ) denoting the polar and azimuth angles in radians of the classical magnetic moment on this site.
ExactDiagonalization.SpinCoherentState — Type
(state::SpinCoherentState)(bases::AbelianBases, table::AbstractDict, dtype=ComplexF64) -> Vector{dtype}Get the vector representation of a spin coherent state with the given Abelian bases and table.
ExactDiagonalization.SpinCoherentState — Method
SpinCoherentState(structure::AbstractDict{Int, <:AbstractVector{<:Number}})
SpinCoherentState(structure::AbstractDict{Int, <:NTuple{2, Number}}; unit::Symbol=:radian)Construct a spin coherent state on a block of lattice sites.
ExactDiagonalization.SpinCoherentStateProjection — Type
SpinCoherentStateProjection <: ActionProjection of states obtained by exact diagonalization method onto spin coherent states.
ExactDiagonalization.SpinCoherentStateProjection — Method
SpinCoherentStateProjection(configuration::SpinCoherentState, polars::AbstractVector{<:Real}, azimuths::AbstractVector{<:Real})
SpinCoherentStateProjection(configuration::SpinCoherentState, np::Integer, na::Integer)Construct a SpinCoherentStateProjection.
ExactDiagonalization.SpinCoherentStateProjectionData — Type
SpinCoherentStateProjectionData <: DataData of spin coherent state projection, including:
polars::Vector{Float64}: global polar angles of the spin coherent states.azimuths::Vector{Float64}: global azimuth angles of the spin coherent states.values::Matrix{Float64}: projection of the state obtained by exact diagonalization method onto the spin coherent states.
ExactDiagonalization.StaticTwoPointCorrelator — Type
StaticTwoPointCorrelator{O<:Union{Operator, Operators}, R<:ReciprocalSpace} <: ActionStatic two-point correlation function.
ExactDiagonalization.StaticTwoPointCorrelatorData — Type
StaticTwoPointCorrelatorData{R<:ReciprocalSpace, V<:Array{Float64}} <: DataData of static two-point correlation function, including:
reciprocalspace::R: reciprocal space to compute the static two-point correlation function.values::V: values of the static two-point correlation function.
ExactDiagonalization.fℤ₂ — Type
fℤ₂ <: ℤ{2}Fermion parity.
ExactDiagonalization.sℤ₂ — Type
sℤ₂ <: ℤ{2}Spin parity.
ExactDiagonalization.ℕ — Type
ℕ <: <: 𝕌₁Concrete Abelian quantum number of the particle number.
ExactDiagonalization.ℤ — Type
ℤ{N} <: SimpleAbelianQuantumNumberAbstract type of ℤₙ quantum numbers.
ExactDiagonalization.ℤ₁ — Type
ℤ₁ <: ℤ{1}ℤ₁, i.e., the trivial quantum number.
ExactDiagonalization.𝕊ᶻ — Type
𝕊ᶻ <: 𝕌₁Concrete Abelian quantum number of the z-component of a spin.
ExactDiagonalization.𝕌₁ — Type
𝕌₁ <: SimpleAbelianQuantumNumberAbstract type of 𝕌₁ quantum numbers.
QuantumLattices.DegreesOfFreedom.CompositeIndex — Method
(index::CompositeIndex)(quantumnumber::Abelian) -> AbelianGet the resulting Abelian quantum number after a CompositeIndex acts upon an initial Abelian quantum number.
QuantumLattices.DegreesOfFreedom.Index — Method
(index::Index)(quantumnumber::Abelian) -> AbelianGet the resulting Abelian quantum number after an Index acts upon an initial Abelian quantum number.
QuantumLattices.DegreesOfFreedom.Metric — Method
Metric(::EDKind{:Abelian}, ::Hilbert{<:Spin}) -> OperatorIndexToTupleGet the index-to-tuple metric for a canonical quantum spin lattice system.
QuantumLattices.DegreesOfFreedom.Metric — Method
Metric(::EDKind{:Binary}, ::Hilbert{<:Fock}) -> OperatorIndexToTupleGet the index-to-tuple metric for a canonical quantum Fock lattice system.
QuantumLattices.QuantumOperators.OperatorProd — Method
(op::OperatorProd)(quantumnumber::Abelian) -> AbelianGet the resulting Abelian quantum number after an OperatorProd acts upon an initial Abelian quantum number.
QuantumLattices.QuantumOperators.OperatorSet — Method
(ops::OperatorSet)(quantumnumber::Abelian) -> AbelianGet the resulting Abelian quantum number after an OperatorSet acts upon an initial Abelian quantum number.
QuantumLattices.QuantumSystems.FockIndex — Method
(index::FockIndex)(quantumnumber::ℤ₁) -> ℤ₁
(index::FockIndex)(quantumnumber::ℕ) -> ℕ
(index::FockIndex)(quantumnumber::𝕊ᶻ) -> 𝕊ᶻ
(index::FockIndex)(quantumnumber::(ℕ ⊠ 𝕊ᶻ)) -> ℕ ⊠ 𝕊ᶻ
(index::FockIndex)(quantumnumber::(𝕊ᶻ ⊠ ℕ)) -> 𝕊ᶻ ⊠ ℕGet the resulting Abelian quantum number after a FockIndex acts upon an initial Abelian quantum number.
Base.Broadcast.broadcast — Method
broadcast(::Type{Sector}, quantumnumbers::OneAtLeast{Abelian}, hilbert::Hilbert, args...; table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> NTuple{fieldcount(typeof(quantumnumbers)), Sector}Construct a set of sectors based on the quantum numbers and a Hilbert space.
Base.Broadcast.broadcast — Method
broadcast(
::Type{Sector}, quantumnumbers::OneAtLeast{Abelian}, hilbert::Hilbert{<:Spin}, partition::NTuple{N, AbstractVector{Int}}=partition(length(hilbert));
table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))
) where N -> NTuple{fieldcount(typeof(quantumnumbers)), AbelianBases}Construct a set of Abelian based based on the quantum numbers and a Hilbert space.
Base.adjoint — Method
adjoint(gs::AbelianGradedSpace) -> typeof(gs)Get the dual of an Abelian graded space.
Base.count — Method
count(basis::BinaryBasis) -> Int
count(basis::BinaryBasis, start::Integer, stop::Integer) -> IntCount the number of occupied single-particle states for a binary basis.
Base.count — Method
count(data::EDEigenData) -> IntCount the number of eigen value-vector-sector groups contained in an EDEigenData.
Base.cumsum — Method
cumsum(rs::RepresentationSpace, i::Union{Integer, CartesianIndex}) -> IntGet the accumulative degenerate dimension up to the ith Abelian quantum number contained in a representation space.
Base.cumsum — Method
cumsum(gs::AbelianGradedSpace{QN}, qn::QN) where {QN<:AbelianQuantumNumber} -> IntGet the accumulative dimension of an Abelian graded space up to a certain Abelian quantum number contained in an Abelian graded space.
Base.eltype — Method
eltype(gf::AbstractGreenFunction)
eltype(::Type{<:AbstractGreenFunction{T}}) where {T<:Number}Get the eltype of an AbstractGreenFunction.
Base.getindex — Method
getindex(gs::AbelianGradedSpace, indexes::AbstractVector{<:Integer}) -> typeof(gs)
getindex(gs::AbelianGradedSpace{QN}, quantumnumbers::AbstractVector{QN}) where {QN<:AbelianQuantumNumber} -> AbelianGradedSpace{QN}Get a subset of an Abelian graded space.
Base.getindex — Method
getindex(gs::AbelianGradedSpace, i::Union{Integer, CartesianIndex})Get the ith Abelian quantum number contained in an Abelian graded space.
Base.getindex — Method
getindex(qn::AbelianQuantumNumberProd, i::Integer) -> SimpleAbelianQuantumNumberGet the ith simple Abelian quantum number in a Deligne tensor product.
Base.getindex — Method
getindex(::Type{Abelian}, ::Type{T}) where {T<:AbelianQuantumNumber} -> Type{T}Overloaded [] for Abelian, i.e., the support of syntax Abelian[T] where T<:AbelianQuantumNumber, which is helpful for the construction of tensor producted Abelian quantum numbers.
Base.isone — Method
isone(basis::BinaryBasis, state::Integer) -> BoolJudge whether the specified single-particle state is occupied for a binary basis.
Base.iszero — Method
iszero(basis::BinaryBasis, state::Integer) -> BoolJudge whether the specified single-particle state is unoccupied for a binary basis.
Base.iterate — Function
iterate(basis::BinaryBasis)
iterate(basis::BinaryBasis, state)Iterate over the numbers of the occupied single-particle states.
Base.length — Method
length(gs::AbelianGradedSpace) -> IntGet the number of inequivalent irreducible representations (i.e., the Abelian quantum numbers) of an Abelian graded space.
Base.length — Method
length(qn::AbelianQuantumNumberProd) -> IntGet the length of a Deligne tensor product.
Base.match — Method
match(sector₁::Sector, sector₂::Sector) -> BoolJudge whether two sectors match each other, that is, whether they can be used together as the bra and ket spaces.
Base.merge — Method
merge(rs::AbelianGradedSpaceProd) -> Tuple{AbelianGradedSpace{eltype(rs)}, Dict{eltype(rs), Vector{NTuple{rank(rs), eltype(rs)}}}}Get the decomposition of the direct product of several Abelian graded spaces and its corresponding fusion processes.
For a set of Abelian graded spaces (gs₁, gs₂, ...), their direct product space can contain several equivalent irreducible representations because for different sets of Abelian quantum numbers (qn₁, qn₂, ...) where qnᵢ∈gsᵢ, the fusion, i.e., ⊗(qn₁, qn₂, ...) may give the same result qn. This function returns the decomposition of the direct product of (gs₁, gs₂, ...) as well as all the fusion processes of each quantum number contained in the decomposition.
Base.pairs — Method
pairs(rs::RepresentationSpace, ::typeof(dimension)) -> RepresentationSpacePairs
pairs(rs::RepresentationSpace, ::typeof(range)) -> RepresentationSpacePairsReturn an iterator that iterates over the pairs of the Abelian quantum numbers and their corresponding (slices of the) degenerate dimensions contained in a representation space.
Base.range — Method
range(bs::AbelianBases) -> AbstractVector{Int}Get the range of the target sector of an AbelianBases in the direct producted bases.
Base.range — Method
range(gs::AbelianGradedSpace, qn::Union{Integer, CartesianIndex}) -> UnitRange{Int}
range(gs::AbelianGradedSpace{QN}, qn::QN) where {QN<:AbelianQuantumNumber} -> UnitRange{Int}Get the slice of the degenerate dimension of an Abelian quantum number contained in an Abelian graded space.
Base.range — Method
range(rs::AbelianGradedSpaceProd, i::Union{Integer, CartesianIndex}) -> AbstractVector{Int}Get the slice of the degenerate dimension of the ith Abelian quantum number in the direct product of several Abelian graded spaces.
Base.range — Method
range(rs::AbelianGradedSpaceSum, i::Union{Integer, CartesianIndex}) -> UnitRange{Int}Get the slice of the degenerate dimension of the ith Abelian quantum number in the direct sum of several Abelian graded spaces.
Base.range — Method
range(rs::AbelianGradedSpaceProd{N, QN}, qns::NTuple{N, QN}) where {N, QN<:AbelianQuantumNumber} -> AbstractVector{Int}Get the slice of the degenerate dimension of the Abelian quantum number fused by qns in the direct product of several Abelian graded spaces.
Base.split — Method
split(target::QN, rs::AbelianGradedSpaceProd{N, QN}; nmax::Real=20) where {N, QN<:AbelianQuantumNumber} -> Set{NTuple{N, QN}}Find a set of splittings of the target Abelian quantum number with respect to the direct product of several Abelian graded spaces.
Base.values — Method
values(qn::AbelianQuantumNumberProd) -> NTuple{rank(qn), Number}Get the values of the simple Abelian quantum numbers in a Deligne tensor product.
Base.values — Method
values(qn::SimpleAbelianQuantumNumber) -> Tuple{Number}Get the value of a simple Abelian quantum number and return it as the sole element of a tuple.
ExactDiagonalization.:⊠ — Method
⊠(bs::BinaryBases, another::Abelian) -> BinaryBases
⊠(another::Abelian, bs::BinaryBases) -> BinaryBasesDeligne tensor product the quantum number of a set of binary bases with another quantum number.
ExactDiagonalization.:⊠ — Method
⊠(qn::SimpleAbelianQuantumNumber, qns::SimpleAbelianQuantumNumber...) -> AbelianQuantumNumberProd
⊠(qn₁::SimpleAbelianQuantumNumber, qn₂::AbelianQuantumNumberProd) -> AbelianQuantumNumberProd
⊠(qn₁::AbelianQuantumNumberProd, qn₂::SimpleAbelianQuantumNumber) -> AbelianQuantumNumberProd
⊠(qn₁::AbelianQuantumNumberProd, qn₂::AbelianQuantumNumberProd) -> AbelianQuantumNumberProdDeligne tensor product of Abelian quantum numbers.
ExactDiagonalization.:⊠ — Method
⊠(QN::Type{<:SimpleAbelianQuantumNumber}, QNS::Type{<:SimpleAbelianQuantumNumber}...) -> Type{AbelianQuantumNumberProd{Tuple{QN, QNS...}}}
⊠(::Type{QN}, ::Type{AbelianQuantumNumberProd{T}}) where {QN<:SimpleAbelianQuantumNumber, T<:Tuple{Vararg{SimpleAbelianQuantumNumber}}} -> Type{AbelianQuantumNumberProd{Tuple{QN, fieldtypes(T)...}}}
⊠(::Type{AbelianQuantumNumberProd{T}}, ::Type{QN}) where {T<:Tuple{Vararg{SimpleAbelianQuantumNumber}}, QN<:SimpleAbelianQuantumNumber} -> Type{AbelianQuantumNumberProd{Tuple{fieldtypes(T)...}, QN}}
⊠(::Type{AbelianQuantumNumberProd{T₁}}, ::Type{AbelianQuantumNumberProd{T₂}}) where {T₁<:Tuple{Vararg{SimpleAbelianQuantumNumber}}, T₂<:Tuple{Vararg{SimpleAbelianQuantumNumber}}} -> Type{AbelianQuantumNumberProd{Tuple{fieldtypes(T₁)..., fieldtypes(T₂)...}}}Deligne tensor product of Abelian quantum numbers.
ExactDiagonalization.basistype — Method
basistype(i::Integer)
basistype(::Type{I}) where {I<:Integer}Get the binary basis type corresponding to an integer or a type of an integer.
ExactDiagonalization.findindex — Method
findindex(position::Integer, gs::AbelianGradedSpace, guess::Integer) -> IntFind the index of an Abelian quantum number in an Abelian graded space beginning at guess whose position in the complete dimension range is position.
ExactDiagonalization.prepare! — Method
prepare!(ed::ED; timer::TimerOutput=edtimer) -> ED
prepare!(ed::Algorithm{<:ED}) -> Algorithm{<:ED}Prepare the matrix representation.
ExactDiagonalization.productable — Method
productable(sector₁::Sector, sector₂::Sector) -> BoolJudge whether two sectors could be direct producted.
ExactDiagonalization.regularize! — Method
regularize!(quantumnumbers::Vector{<:AbelianQuantumNumber}, dimensions::Vector{Int}; check::Bool=false) -> Tuple{typeof(quantumnumbers), typeof(dimensions), Vector{Int}}In place regularization of the input Abelian quantum numbers and their corresponding degenerate dimensions.
After the regularization, the Abelian quantum numbers will be sorted in the ascending order and duplicates will be merged together. The degenerate dimensions will be processed accordingly. When check is true, this function also check whether all input degenerate dimensions are positive. The regularized Abelian quantum numbers and degenerate dimensions, as well as the permutation vector that sorts the input Abelian quantum numbers, will be returned.
ExactDiagonalization.regularize — Method
regularize(quantumnumbers::AbstractVector{<:AbelianQuantumNumber}, dimension::AbstractVector{<:Integer}; check::Bool=false) -> Tuple{Vector{eltype(quantumnumbers)}, Vector{Int}, Vector{Int}}Regularize of the input Abelian quantum numbers and their corresponding degenerate dimensions.
See regularize!.
ExactDiagonalization.release! — Method
release!(ed::ED; gc::Bool=true) -> ED
release!(ed::Algorithm{<:ED}; gc::Bool=true) -> Algorithm{<:ED}Release the memory source used in preparing the matrix representation. If gc is true, call the garbage collection immediately.
ExactDiagonalization.sumable — Method
sumable(sector₁::Sector, sector₂::Sector) -> BoolJudge whether two sectors could be direct summed.
LinearAlgebra.eigen — Method
eigen(ed::ED, sectors::Union{Abelian, Sector}...; timer::TimerOutput=edtimer, release::Bool=false, kwargs...) -> EDEigenData
eigen(ed::Algorithm{<:ED}, sectors::Union{Abelian, Sector}...; release::Bool=false, kwargs...) -> EDEigenDataSolve the eigen problem by the restarted Lanczos method provided by the Arpack package.
LinearAlgebra.eigen — Method
eigen(m::EDMatrix; nev::Int=1, which::Symbol=:SR, tol::Real=1e-12, maxiter::Int=300, v₀::Union{AbstractVector{<:Number}, Int}=dimension(m.bra), krylovdim::Int=max(20, 2*nev+1), verbosity::Int=0) -> EDEigenDataSolve the eigen problem by the restarted Lanczos method provided by the KrylovKit package.
LinearAlgebra.eigen — Method
eigen(
ms::OperatorSum{<:EDMatrix};
nev::Int=1,
which::Symbol=:SR,
tol::Real=1e-12,
maxiter::Int=300,
v₀::Union{Dict{<:Abelian, <:Union{AbstractVector{<:Number}, Int}}, Dict{<:Sector, <:Union{AbstractVector{<:Number}, Int}}}=Dict(Abelian(m.ket)=>dimension(m.ket) for m in ms),
krylovdim::Int=max(20, 2*nev+1),
verbosity::Int=0,
timer::TimerOutput=edtimer
) -> EDEigenDataSolve the eigen problem by the restarted Lanczos method provided by the Arpack package.
LinearAlgebra.rank — Method
rank(qn::AbelianQuantumNumberProd) -> Int
rank(::Type{<:AbelianQuantumNumberProd}) -> IntGet the rank of a Deligne tensor product of simple Abelian quantum numbers.
LinearAlgebra.rank — Method
rank(rs::CompositeAbelianGradedSpace) -> Int
rank(::Type{<:CompositeAbelianGradedSpace{N}}) where N -> IntGet the number of Abelian graded spaces in the direct sum.
LinearAlgebra.rank — Method
rank(gf::GreenFunction) -> IntGet the rank of a GreenFunction.
LinearAlgebra.rank — Method
rank(gf::RetardedGreenFunction) -> IntGet the rank of a RetardedGreenFunction.
QuantumLattices.:⊕ — Method
⊕(gs::AbelianGradedSpace, gses::AbelianGradedSpace...) -> AbelianGradedSpaceSum
⊕(gs::AbelianGradedSpace, rs::AbelianGradedSpaceSum) -> AbelianGradedSpaceSum
⊕(rs::AbelianGradedSpaceSum, gs::AbelianGradedSpace) -> AbelianGradedSpaceSum
⊕(rs₁::AbelianGradedSpaceSum, rs₂::AbelianGradedSpaceSum) -> AbelianGradedSpaceSumGet the direct sum of some Abelian graded spaces.
QuantumLattices.:⊗ — Method
⊗(gs::AbelianGradedSpace, gses::AbelianGradedSpace...) -> AbelianGradedSpaceProd
⊗(gs::AbelianGradedSpace, rs::AbelianGradedSpaceProd) -> AbelianGradedSpaceProd
⊗(rs::AbelianGradedSpaceProd, gs::AbelianGradedSpace) -> AbelianGradedSpaceProd
⊗(rs₁::AbelianGradedSpaceProd, rs₂::AbelianGradedSpaceProd) -> AbelianGradedSpaceProdGet the direct product of some Abelian graded spaces.
QuantumLattices.:⊗ — Method
⊗(qn::AbelianQuantumNumber, qns::AbelianQuantumNumber...) -> eltype(qns)Get the direct product of some AbelianQuantumNumbers.
QuantumLattices.:⊗ — Method
⊗(bs₁::BinaryBases, bs₂::BinaryBases) -> BinaryBasesGet the direct product of two sets of binary bases.
QuantumLattices.:⊗ — Method
⊗(basis₁::BinaryBasis, basis₂::BinaryBasis) -> BinaryBasisGet the direct product of two binary bases.
QuantumLattices.DegreesOfFreedom.partition — Method
partition(n::Int) -> NTuple{2, Vector{Int}}Get the default partition of n local Hilbert spaces.
QuantumLattices.QuantumOperators.matrix — Function
matrix(index::OperatorIndex, graded::Graded, dtype::Type{<:Number}=ComplexF64) -> Matrix{dtype}Get the matrix representation of an OperatorIndex on an Abelian graded space.
QuantumLattices.QuantumOperators.matrix — Function
matrix(ops::Operators, braket::NTuple{2, Sector}, table::AbstractDict, dtype=scalartype(ops)) -> SparseMatrixCSC{dtype, Int}Get the CSC-formed sparse matrix representation of a set of operators.
Here, table specifies the order of the operator indexes.
QuantumLattices.QuantumOperators.matrix — Function
matrix(index::SpinIndex, graded::Graded, dtype::Type{<:Number}=ComplexF64) -> Matrix{dtype}Get the matrix representation of a SpinIndex on an Abelian graded space.
QuantumLattices.QuantumOperators.matrix — Function
matrix(op::OperatorIndex, braket::NTuple{2, Sector}, table::AbstractDict, dtype=Int) -> SparseMatrixCSC{dtype, Int}Get the CSC-formed sparse matrix representation of an operator index.
QuantumLattices.QuantumOperators.matrix — Method
matrix(ed::ED, sectors::Union{Abelian, Sector}...; timer::TimerOutput=edtimer, release::Bool=false) -> OperatorSum{<:EDMatrix}
matrix(ed::Algorithm{<:ED}, sectors::Union{Abelian, Sector}...; release::Bool=false) -> OperatorSum{<:EDMatrix}Get the sparse matrix representation of a quantum lattice system in the target space.
QuantumLattices.QuantumOperators.matrix — Method
matrix(op::Operator{V, <:OneAtLeast{OperatorIndex}}, braket::NTuple{2, AbelianBases}, table::AbstractDict, dtype=V) where V -> SparseMatrixCSC{dtype, Int}Get the CSC-formed sparse matrix representation of an operator.
Here, table specifies the order of the operator indexes.
QuantumLattices.QuantumOperators.matrix — Method
matrix(op::Operator{V, <:OneAtLeast{OperatorIndex}}, braket::NTuple{2, BinaryBases}, table::AbstractDict, dtype=V) where V -> SparseMatrixCSC{dtype, Int}Get the CSC-formed sparse matrix representation of an operator.
Here, table specifies the order of the operator indexes.
QuantumLattices.QuantumOperators.matrix — Method
matrix(op::Operator{V, Tuple{}}, braket::NTuple{2, Sector}, table::AbstractDict, dtype=V) where V -> SparseMatrixCSC{dtype, Int}Get the CSC-formed sparse matrix representation of a scalar operator.
QuantumLattices.Spatials.period — Method
period(qn::AbelianQuantumNumberProd, i::Integer) -> Number
period(::Type{AbelianQuantumNumberProd{T}}, i::Integer) where {T<:Tuple{Vararg{SimpleAbelianQuantumNumber}}} -> NumberGet the period of the ith simple Abelian number contained in a Deligne tensor product.
QuantumLattices.Spatials.period — Method
period(qn::SimpleAbelianQuantumNumber) -> Number
period(::Type{QN}) where {QN<:SimpleAbelianQuantumNumber} -> NumberGet the period of a simple Abelian quantum number.
QuantumLattices.Spatials.periods — Method
periods(qn::AbelianQuantumNumber) -> Tuple{Vararg{Number}}
periods(::Type{QN}) where {QN<:AbelianQuantumNumber} -> Tuple{Vararg{Number}}Get the periods of Abelian quantum numbers.
QuantumLattices.decompose — Method
decompose(rs::CompositeAbelianGradedSpace; expand::Bool=true) -> Tuple{AbelianGradedSpace{eltype(rs)}, Vector{Int}}Decompose a composite of several Abelian graded spaces to the canonical one.
When expand is false, the corresponding permutation vector that sorts the Abelian quantum numbers will be returned as well. When expand is true, the expanded dimension indexes of the permutation vector that sorts the Abelian quantum numbers will be returned as well.
QuantumLattices.dimension — Method
dimension(gs::AbelianGradedSpace, qn::Union{Integer, CartesianIndex}) -> Int
dimension(gs::AbelianGradedSpace{QN}, qn::QN) where {QN<:AbelianQuantumNumber} -> IntGet the degenerate dimension of an Abelian quantum number contained in an Abelian graded space.
QuantumLattices.dimension — Method
dimension(rs::AbelianGradedSpaceProd, i::Union{Integer, CartesianIndex}) -> IntGet the degenerate dimension of the ith Abelian quantum number in the direct product of several Abelian graded spaces.
QuantumLattices.dimension — Method
dimension(rs::AbelianGradedSpaceProd) -> IntGet the total dimension of the direct product of several Abelian graded spaces.
QuantumLattices.dimension — Method
dimension(rs::AbelianGradedSpaceSum, i::Union{Integer, CartesianIndex}) -> IntGet the degenerate dimension of the ith Abelian quantum number in the direct sum of several Abelian graded spaces.
QuantumLattices.dimension — Method
dimension(rs::AbelianGradedSpaceSum) -> IntGet the total dimension of the direct sum of several Abelian graded spaces.
QuantumLattices.dimension — Method
dimension(gs::AbelianGradedSpace) -> IntGet the total dimension of an Abelian graded space.
QuantumLattices.dimension — Method
dimension(gf::GreenFunction) -> IntGet the dimension of the Krylov space expanded to obtained a GreenFunction.
QuantumLattices.dimension — Method
dimension(rs::AbelianGradedSpaceProd{N, QN}, qns::NTuple{N, QN}) where {N, QN<:AbelianQuantumNumber} -> IntGet the degenerate dimension of the Abelian quantum number fused by qns in the direct product of several Abelian graded spaces.
QuantumLattices.reset! — Method
reset!(
gf::GreenFunction, H::AbstractMatrix{<:Number}, V::AbstractVector{<:AbstractVector{<:Number}}, E₀::Real, kind::Symbol;
ranks::AbstractVector{<:Integer}=1:rank(gf), dimensions::AbstractVector{<:Integer}=1:dimension(gf), tol::Real=1e-12
) -> typeof(gf)Reset (a block) of GreenFunction.
QuantumLattices.value — Method
value(qn::AbelianQuantumNumberProd, i::Integer) -> NumberGet the value of the ith simple Abelian quantum number in a Deligne tensor product.
QuantumLattices.value — Method
value(qn::SimpleAbelianQuantumNumber) -> NumberGet the value of a simple Abelian quantum number.