Core of ExactDiagonalization

ExactDiagonalization.basicoptionsConstant
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.

source
ExactDiagonalization.AbelianBasesType
AbelianBases{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:

  1. 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
  2. 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.

source
ExactDiagonalization.AbelianBasesMethod
AbelianBases(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.

source
ExactDiagonalization.AbelianBasesMethod
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.

source
ExactDiagonalization.BinaryBasesMethod
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.

source
ExactDiagonalization.BinaryBasesMethod
BinaryBases(spindws, spinups, sz::𝕊ᶻ)

Construct a set of binary bases that preserves the spin z-component but not the particle number conservation.

source
ExactDiagonalization.BinaryBasesMethod
BinaryBases(states, particle::ℕ)
BinaryBases(nstate::Integer, particle::ℕ)

Construct a set of binary bases that preserves the particle number conservation.

source
ExactDiagonalization.BinaryBasisType
BinaryBasis{I<:Unsigned}

Binary basis represented by an unsigned integer.

Here, we adopt the following common rules:

  1. In a binary basis, a bit of an unsigned integer represents a single-particle state that can be occupied (1) or unoccupied (0).
  2. 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.
  3. 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$.

source
ExactDiagonalization.BinaryBasisMethod
BinaryBasis(states; filter=index->true)
BinaryBasis{I}(states; filter=index->true) where {I<:Unsigned}

Construct a binary basis with the given occupied states.

source
ExactDiagonalization.EDType
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.

source
ExactDiagonalization.EDType
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.

source
ExactDiagonalization.EDType
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.

source
ExactDiagonalization.EDType
ED{K<:EDKind, L<:Union{AbstractLattice, Nothing}, S<:Generator{<:Operators}, M<:EDMatrixization, H<:CategorizedGenerator{<:OperatorSum{<:EDMatrix}}} <: Frontend

Exact diagonalization method of a quantum lattice system.

source
ExactDiagonalization.EDEigenDataType
EDEigenData{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.

source
ExactDiagonalization.EDMatrixType
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.

source
ExactDiagonalization.EDMatrixMethod
EDMatrix(m::SparseMatrixCSC, sector::Sector)
EDMatrix(m::SparseMatrixCSC, braket::NTuple{2, Sector})
EDMatrix(m::SparseMatrixCSC, bra::Sector, ket::Sector)

Construct a matrix representation when

  1. the bra and ket Hilbert spaces share the same bases;
  2. the bra and ket Hilbert spaces may be different;
  3. the bra and ket Hilbert spaces may or may not be the same.
source
ExactDiagonalization.EDMatrixizationMethod
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.

source
ExactDiagonalization.SectorType
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))
) -> AbelianBases

Construct the Abelian bases of a spin Hilbert space with the specified quantum number.

source
ExactDiagonalization.SectorType
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))) -> BinaryBases

Construct the binary bases of a Hilbert space with the specified quantum number.

source
ExactDiagonalization.SectorType
Sector

A sector of the Hilbert space which forms the bases of an irreducible representation of the Hamiltonian of a quantum lattice system.

source
ExactDiagonalization.SpinCoherentStateType
(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.

source
ExactDiagonalization.SpinCoherentStateType
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.

source
ExactDiagonalization.SpinCoherentStateMethod
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.

source
ExactDiagonalization.SpinCoherentStateProjectionMethod
SpinCoherentStateProjection(configuration::SpinCoherentState, polars::AbstractVector{<:Real}, azimuths::AbstractVector{<:Real})
SpinCoherentStateProjection(configuration::SpinCoherentState, np::Integer, na::Integer)

Construct a SpinCoherentStateProjection.

source
ExactDiagonalization.SpinCoherentStateProjectionDataType
SpinCoherentStateProjectionData <: Data

Data of spin coherent state projection, including:

  1. polars::Vector{Float64}: global polar angles of the spin coherent states.
  2. azimuths::Vector{Float64}: global azimuth angles of the spin coherent states.
  3. values::Matrix{Float64}: projection of the state obtained by exact diagonalization method onto the spin coherent states.
source
ExactDiagonalization.StaticTwoPointCorrelatorDataType
StaticTwoPointCorrelatorData{R<:ReciprocalSpace, V<:Array{Float64}} <: Data

Data of static two-point correlation function, including:

  1. reciprocalspace::R: reciprocal space to compute the static two-point correlation function.
  2. values::V: values of the static two-point correlation function.
source
QuantumLattices.QuantumNumbers.GradedMethod
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.

source
Base.Broadcast.broadcastMethod
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.

source
Base.Broadcast.broadcastMethod
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.

source
Base.countMethod
count(basis::BinaryBasis) -> Int
count(basis::BinaryBasis, start::Integer, stop::Integer) -> Int

Count the number of occupied single-particle states for a binary basis.

source
Base.countMethod
count(data::EDEigenData) -> Int

Count the number of eigen value-vector-sector groups contained in an EDEigenData.

source
Base.isoneMethod
isone(basis::BinaryBasis, state::Integer) -> Bool

Judge whether the specified single-particle state is occupied for a binary basis.

source
Base.iszeroMethod
iszero(basis::BinaryBasis, state::Integer) -> Bool

Judge whether the specified single-particle state is unoccupied for a binary basis.

source
Base.iterateFunction
iterate(basis::BinaryBasis)
iterate(basis::BinaryBasis, state)

Iterate over the numbers of the occupied single-particle states.

source
Base.matchMethod
match(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.

source
Base.oneMethod
one(basis::BinaryBasis, state::Integer) -> BinaryBasis

Get a new binary basis with the specified single-particle state occupied.

source
Base.rangeMethod
range(bs::AbelianBases) -> AbstractVector{Int}

Get the range of the target sector of an AbelianBases in the direct producted bases.

source
Base.zeroMethod
zero(basis::BinaryBasis, state::Integer) -> BinaryBasis

Get a new binary basis with the specified single-particle state unoccupied.

source
ExactDiagonalization.basistypeMethod
basistype(i::Integer)
basistype(::Type{I}) where {I<:Integer}

Get the binary basis type corresponding to an integer or a type of an integer.

source
ExactDiagonalization.prepare!Method
prepare!(ed::ED; timer::TimerOutput=edtimer) -> ED
prepare!(ed::Algorithm{<:ED}) -> Algorithm{<:ED}

Prepare the matrix representation.

source
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.

source
LinearAlgebra.eigenMethod
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...) -> EDEigenData

Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.

source
LinearAlgebra.eigenMethod
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) -> EDEigenData

Solve the eigen problem by the restarted Lanczos method provided by the KrylovKit package.

source
LinearAlgebra.eigenMethod
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
) -> EDEigenData

Solve the eigen problem by the restarted Lanczos method provided by the Arpack package.

source
QuantumLattices.:⊗Method
⊗(bs₁::BinaryBases, bs₂::BinaryBases) -> BinaryBases

Get the direct product of two sets of binary bases.

source
QuantumLattices.:⊗Method
⊗(basis₁::BinaryBasis, basis₂::BinaryBasis) -> BinaryBasis

Get the direct product of two binary bases.

source
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.

source
QuantumLattices.QuantumOperators.matrixFunction
matrix(index::OperatorIndex, graded::Graded, dtype::Type{<:Number}=ComplexF64) -> Matrix{dtype}

Get the matrix representation of an OperatorIndex on an Abelian graded space.

source
QuantumLattices.QuantumOperators.matrixFunction
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.

source
QuantumLattices.QuantumOperators.matrixMethod
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.

source
QuantumLattices.QuantumOperators.matrixMethod
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.

source
QuantumLattices.QuantumOperators.matrixMethod
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.

source
QuantumLattices.QuantumOperators.matrixMethod
matrix(op::Operator{V, Tuple{}}, braket::NTuple{2, Sector}, table::AbstractDict, dtype=V) where V -> SparseMatrixCSC{V, Int}

Get the CSC-formed sparse matrix representation of a scalar operator.

source