Core of ExactDiagonalization
ExactDiagonalization.edtimer
— Constantconst edtimer = TimerOutput()
Default shared timer for all exact diagonalization methods.
ExactDiagonalization.AbelianBases
— TypeAbelianBases{A<:Abelian, N} <: Sector
A 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
— MethodAbelianBases(locals::AbstractVector{Int}, partition::NTuple{N, AbstractVector{Int}}=partition(length(locals))) where N
Construct a set of spin bases that subjects to no quantum number conservation.
ExactDiagonalization.AbelianBases
— MethodAbelianBases(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.BinaryBases
— TypeBinaryBases{A<:Abelian, B<:BinaryBasis, T<:AbstractVector{B}} <: Sector
A set of binary bases.
ExactDiagonalization.BinaryBases
— MethodBinaryBases(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
— MethodBinaryBases(spindws, spinups, sz::𝕊ᶻ)
Construct a set of binary bases that preserves the spin z-component but not the particle number conservation.
ExactDiagonalization.BinaryBases
— MethodBinaryBases(states, particle::ℕ)
BinaryBases(nstate::Integer, particle::ℕ)
Construct a set of binary bases that preserves the particle number conservation.
ExactDiagonalization.BinaryBases
— MethodBinaryBases(states)
BinaryBases(nstate::Integer)
Construct a set of binary bases that subjects to no quantum number conservation.
ExactDiagonalization.BinaryBasis
— TypeBinaryBasis{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
— MethodBinaryBasis(states; filter=index->true)
BinaryBasis{I}(states; filter=index->true) where {I<:Unsigned}
Construct a binary basis with the given occupied states.
ExactDiagonalization.BinaryBasisRange
— TypeBinaryBasisRange{I<:Unsigned} <: VectorSpace{BinaryBasis{I}}
A continuous range of binary basis from 0 to 2^n-1.
ExactDiagonalization.ED
— TypeED(system::Generator{<:Operators}, targetspace::TargetSpace, dtype::Type{<:Number}=scalartype(system); timer::TimerOutput=edtimer, delay::Bool=false)
ED(lattice::Union{AbstractLattice, Nothing}, system::Generator{<:Operators}, targetspace::TargetSpace, dtype::Type{<:Number}=scalartype(system); timer::TimerOutput=edtimer, delay::Bool=false)
Construct the exact diagonalization method for a quantum lattice system.
ExactDiagonalization.ED
— TypeED(
lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, targetspace::TargetSpace=TargetSpace(hilbert), boundary::Boundary=plain, dtype::Type{<:Number}=valtype(terms);
neighbors::Union{Int, Neighbors}=nneighbor(terms), timer::TimerOutput=edtimer, delay::Bool=false
)
Construct the exact diagonalization method for a quantum lattice system.
ExactDiagonalization.ED
— TypeED(
lattice::AbstractLattice, hilbert::Hilbert, terms::OneOrMore{Term}, quantumnumbers::OneOrMore{Abelian}, boundary::Boundary=plain, dtype::Type{<:Number}=valtype(terms);
neighbors::Union{Int, Neighbors}=nneighbor(terms), timer::TimerOutput=edtimer, delay::Bool=false
)
Construct the exact diagonalization method for a quantum lattice system.
ExactDiagonalization.ED
— TypeED{K<:EDKind, L<:Union{AbstractLattice, Nothing}, S<:Generator{<:Operators}, M<:EDMatrixization, H<:CategorizedGenerator{<:OperatorSum{<:EDMatrix}}} <: Frontend
Exact diagonalization method of a quantum lattice system.
ExactDiagonalization.EDEigenData
— TypeEDEigenData{V<:Number, T<:Number, S<:Sector} <: Data
Eigen 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
— TypeEDKind{K}
Kind of the exact diagonalization method applied to a quantum lattice system.
ExactDiagonalization.EDKind
— MethodEDKind(::Type{<:FockIndex})
Kind of the exact diagonalization method applied to a canonical quantum Fock lattice system.
ExactDiagonalization.EDKind
— MethodEDKind(::Type{<:SpinIndex})
Kind of the exact diagonalization method applied to a canonical quantum spin lattice system.
ExactDiagonalization.EDMatrix
— TypeEDMatrix{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
— MethodEDMatrix(m::SparseMatrixCSC, sector::Sector)
EDMatrix(m::SparseMatrixCSC, braket::NTuple{2, Sector})
Construct a matrix representation when
- the ket and bra Hilbert spaces share the same bases;
- the ket and bra Hilbert spaces may be different.
ExactDiagonalization.EDMatrixization
— TypeEDMatrixization{D<:Number, S<:Sector, T<:AbstractDict} <: Matrixization
Matrixization of a quantum lattice system on a target Hilbert space.
ExactDiagonalization.EDMatrixization
— MethodEDMatrixization{D}(target::TargetSpace) where {D<:Number}
Construct a matrixization.
ExactDiagonalization.Sector
— TypeSector(hilbert::Hilbert{<:Spin}, partition::Tuple{Vararg{AbstractVector{Int}}}=partition(length(hilbert)); table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> AbelianBases
Sector(hilbert::Hilbert{<:Spin}, quantumnumber::Abelian, partition::Tuple{Vararg{AbstractVector{Int}}}=partition(length(hilbert)); table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> AbelianBases
Construct the Abelian bases of a spin Hilbert space with the specified quantum number.
ExactDiagonalization.Sector
— TypeSector(hilbert::Hilbert{<:Fock}, basistype::Type{<:Unsigned}=UInt) -> BinaryBases
Sector(hilbert::Hilbert{<:Fock}, quantumnumber::ℕ, basistype::Type{<:Unsigned}=UInt; table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> BinaryBases
Sector(hilbert::Hilbert{<:Fock}, quantumnumber::Union{𝕊ᶻ, Abelian[ℕ ⊠ 𝕊ᶻ], Abelian[𝕊ᶻ ⊠ ℕ]}, basistype::Type{<:Unsigned}=UInt; table::AbstractDict=Table(hilbert, Metric(EDKind(hilbert), hilbert))) -> BinaryBases
Construct the binary bases of a Hilbert space with the specified quantum number.
ExactDiagonalization.Sector
— TypeSector
A sector of the Hilbert space which forms the bases of an irreducible representation of the Hamiltonian of a quantum lattice system.
ExactDiagonalization.SectorFilter
— TypeSectorFilter{S} <: LinearTransformation
Filter the target bra and ket Hilbert spaces.
ExactDiagonalization.TargetSpace
— TypeTargetSpace{S<:Sector, T<:AbstractDict} <: VectorSpace{S}
Target Hilbert space in which the exact diagonalization method is performed, which could be the direct sum of several sectors.
ExactDiagonalization.TargetSpace
— MethodTargetSpace(hilbert::Hilbert, args...)
TargetSpace(hilbert::Hilbert, table::AbstractDict, args...)
TargetSpace(hilbert::Hilbert, quantumnumbers::OneOrMore{Abelian}, args...)
Construct a target space from the total Hilbert space and the associated quantum numbers.
ExactDiagonalization.TargetSpace
— MethodTargetSpace(hilbert::Hilbert{<:Spin}, quantumnumbers::OneOrMore{Abelian}, table::AbstractDict, partition::NTuple{N, AbstractVector{Int}}=partition(length(hilbert))) where N
Construct a target space from the total Hilbert space and the associated quantum numbers.
ExactDiagonalization.TargetSpace
— MethodTargetSpace(hilbert::Hilbert{<:Fock}, quantumnumbers::OneOrMore{Abelian}, table::AbstractDict, basistype::Type{<:Unsigned}=UInt)
Construct a target space from the total Hilbert space and the associated quantum numbers.
QuantumLattices.DegreesOfFreedom.Metric
— MethodMetric(::EDKind{:Abelian}, ::Hilbert{<:Spin}) -> OperatorIndexToTuple
Get the index-to-tuple metric for a canonical quantum spin lattice system.
QuantumLattices.DegreesOfFreedom.Metric
— MethodMetric(::EDKind{:Binary}, ::Hilbert{<:Fock}) -> OperatorIndexToTuple
Get the index-to-tuple metric for a canonical quantum Fock lattice system.
QuantumLattices.QuantumNumbers.Abelian
— MethodAbelian(bs::AbelianBases)
Get the Abelian quantum number of a set of spin bases.
QuantumLattices.QuantumNumbers.Abelian
— MethodAbelian(bs::BinaryBases)
Get the Abelian quantum number of a set of binary bases.
QuantumLattices.QuantumNumbers.Graded
— MethodGraded{ℤ₁}(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.
Base.count
— Methodcount(basis::BinaryBasis) -> Int
count(basis::BinaryBasis, start::Integer, stop::Integer) -> Int
Count the number of occupied single-particle states for a binary basis.
Base.isone
— Methodisone(basis::BinaryBasis, state::Integer) -> Bool
Judge whether the specified single-particle state is occupied for a binary basis.
Base.iszero
— Methodiszero(basis::BinaryBasis, state::Integer) -> Bool
Judge whether the specified single-particle state is unoccupied for a binary basis.
Base.iterate
— Functioniterate(basis::BinaryBasis)
iterate(basis::BinaryBasis, state)
Iterate over the numbers of the occupied single-particle states.
Base.match
— Methodmatch(sector₁::Sector, sector₂::Sector) -> Bool
Judge whether two sectors match each other, that is, whether they can be used together as the bra and ket spaces.
Base.one
— Methodone(basis::BinaryBasis, state::Integer) -> BinaryBasis
Get a new binary basis with the specified single-particle state occupied.
Base.range
— Methodrange(bs::AbelianBases) -> AbstractVector{Int}
Get the range of the target sector of an AbelianBases
in the direct producted bases.
Base.zero
— Methodzero(basis::BinaryBasis, state::Integer) -> BinaryBasis
Get a new binary basis with the specified single-particle state unoccupied.
ExactDiagonalization.basistype
— Methodbasistype(i::Integer)
basistype(::Type{I}) where {I<:Integer}
Get the binary basis type corresponding to an integer or a type of an integer.
ExactDiagonalization.prepare!
— Methodprepare!(ed::ED; timer::TimerOutput=edtimer) -> ED
prepare!(ed::Algorithm{<:ED}) -> Algorithm{<:ED}
Prepare the matrix representation.
ExactDiagonalization.productable
— Methodproductable(sector₁::Sector, sector₂::Sector) -> Bool
Judge whether two sectors could be direct producted.
ExactDiagonalization.release!
— Methodrelease!(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
— Methodsumable(sector₁::Sector, sector₂::Sector) -> Bool
Judge whether two sectors could be direct summed.
LinearAlgebra.eigen
— Methodeigen(ed::ED, sectors::Union{Abelian, Sector}...; timer::TimerOutput=edtimer, kwargs...) -> EDEigenData
eigen(ed::Algorithm{<:ED}, sectors::Union{Abelian, Sector}...; kwargs...) -> EDEigenData
Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.
LinearAlgebra.eigen
— Methodeigen(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) -> EDEigenData
Solve the eigen problem by the restarted Lanczos method provided by the KrylovKit package.
LinearAlgebra.eigen
— Methodeigen(
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
) -> EDEigenData
Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.
QuantumLattices.:⊕
— Method⊕(target::TargetSpace, sectors::Sector...) -> TargetSpace
Get the direct sum of a target space with several sectors.
QuantumLattices.:⊗
— Method⊗(bs₁::BinaryBases, bs₂::BinaryBases) -> BinaryBases
Get the direct product of two sets of binary bases.
QuantumLattices.:⊗
— Method⊗(basis₁::BinaryBasis, basis₂::BinaryBasis) -> BinaryBasis
Get the direct product of two binary bases.
QuantumLattices.:⊠
— Method⊠(bs::BinaryBases, another::Abelian) -> BinaryBases
⊠(another::Abelian, bs::BinaryBases) -> BinaryBases
Deligne tensor product the quantum number of a set of binary bases with another quantum number.
QuantumLattices.DegreesOfFreedom.partition
— Methodpartition(n::Int) -> NTuple{2, Vector{Int}}
Get the default partition of n local Hilbert spaces.
QuantumLattices.QuantumOperators.matrix
— Functionmatrix(op::Operator, braket::NTuple{2, BinaryBases}, table::AbstractDict, dtype=valtype(op); kwargs...) -> 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
— Functionmatrix(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
— Functionmatrix(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
— Functionmatrix(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
— Functionmatrix(op::Operator, braket::NTuple{2, AbelianBases}, table::AbstractDict, dtype=valtype(op); kwargs...) -> 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
— Methodmatrix(ed::ED, sectors::Union{Abelian, Sector}...; timer::TimerOutput=edtimer, kwargs...) -> OperatorSum{<:EDMatrix}
matrix(ed::Algorithm{<:ED}, sectors::Union{Abelian, Sector}...; kwargs...) -> OperatorSum{<:EDMatrix}
Get the sparse matrix representation of a quantum lattice system in the target space.