API
QuantumHamiltonian.AbstractHilbertSpaceRepresentation
— TypeAbstractHilbertSpaceRepresentation{BR, S}
QuantumHamiltonian.AbstractOperator
— TypeAbstractOperator{S<:Number}
Represent an abstract operator in Hilbert space.
QuantumHamiltonian.AbstractOperatorRepresentation
— TypeAbstractOperatorRepresentation{S}
QuantumHamiltonian.DecomposedHilbertSpaceRepresentation
— TypeDecomposedHilbertSpaceRepresentation
Represents the decomposed Hilbert space representation. i.e. subspace of the total Hilbert space, with basis states grouped by tags.
QuantumHamiltonian.DictIndexedVector
— TypeQuantumHamiltonian.GlobalBitFlip
— TypeGlobalBitFlip
Global bit flip operation. For spin half systems, this amounts to ℤ₂ spin flip.
QuantumHamiltonian.HilbertSpace
— TypeHilbertSpace{QN, TS}
Abstract Hilbert space with quantum number type QN
.
Examples
julia> using QuantumHamiltonian
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site]);
QuantumHamiltonian.HilbertSpaceRepresentation
— TypeHilbertSpaceRepresentation{BR, HS, BasisType}
QuantumHamiltonian.HilbertSpaceSector
— TypeHilbertSpaceSector{QN}
Hilbert space sector.
QuantumHamiltonian.IntegerModulo
— TypeIntegerModulo{N}
Implement Zₙ.
QuantumHamiltonian.NullOperator
— TypeNullOperator
A null operator, i.e. 0.
QuantumHamiltonian.OperatorRepresentation
— TypeOperatorRepresentation{HSR, S, O}
Operator representation of given operator of type O
.
QuantumHamiltonian.PureOperator
— TypePureOperator{Scalar, BR}
Represents an operator $α (P₁ ⊗ P₂ ⊗ … ⊗ Pₙ)$ where $Pᵢ$ is either identity (when bitmask is set to zero), or projection $|rᵢ⟩⟨cᵢ|$ (when bitmask is set to one).
See also: pure_operator
Fields
bitmask :: BR
bitrow :: BR
bitcol :: BR
amplitude :: Scalar
QuantumHamiltonian.ReducedHilbertSpaceRepresentation
— TypeReducedHilbertSpaceRepresentation{HSR, BR, C}
Representation of the symmetry-reduced hilbert space. Currently only supports Translation group (i.e. Abelian group). ```
QuantumHamiltonian.ReducedOperatorRepresentation
— TypeReducedOperatorRepresentation{RHSR, O, S, BR}
Representation of an operator of type O
in the symmetry-reduced hilbert space representation of type RHSR
.
QuantumHamiltonian.Site
— TypeSite{QN}
A site with quantum number type QN
.
Examples
julia> using QuantumHamiltonian
julia> site = Site([State("Up", 1), State("Dn", -1)]);
QuantumHamiltonian.SortedIndexedVector
— TypeQuantumHamiltonian.SparseState
— Typestruct SparseState{Scalar<:Number, BR}
Represents a vector in unrestricted Hilbert space.
QuantumHamiltonian.SparseState
— MethodSparseState(hsrep, state_rep, tol=√eps(Float64))
Make a SparseState
from a representation
QuantumHamiltonian.State
— TypeState{QN}
State with quantum number type QN
.
Examples
julia> using QuantumHamiltonian
julia> up = State("Up", 1)
State{Tuple{Int64}}("Up", (1,))
julia> State("Dn", (-1, 1))
State{Tuple{Int64, Int64}}("Dn", (-1, 1))
QuantumHamiltonian.SumOperator
— TypeSumOperator{Scalar, BR}
Represents a sum of pure operators.
Fields
terms::Vector{PureOperator{Scalar,BR}}
Base.valtype
— Methodvaltype(lhs::Type{<:AbstractOperator{S}})
Returns the valtype
(scalar type) of the given AbstractOperator.
LatticeTools.dimension
— Methoddimension
Dimension of the Concrete Hilbert space, i.e. number of basis vectors.
LatticeTools.dimension
— Methoddimension(arg::ReducedHilbertSpaceRepresentation{HSR, BR, C}) -> Int
Dimension of the given reduced hilbert space representation, i.e. number of basis elements.
LatticeTools.dimension
— Methoddimension(site)
Hilbert space dimension of a given site (= number of states).
LatticeTools.numsites
— Methodnumsites(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call numsites
with basespace of hsr
.
LatticeTools.numsites
— Methodnumsites(hss::HilbertSpaceSector, args...;kwargs...)
Call numsites
with basespace of hss
.
QuantumHamiltonian.apply!
— Methodapply!(out, opr, state)
Perform out += opr * state
. Apply the operator representation opr
to the row vector state
and add it to the row vector out
. Return sum of errors and sum of error-squared. Call apply_serial!
if Threads.nthreads() == 1
, and apply_parallel!
otherwise.
QuantumHamiltonian.apply!
— Methodapply!(out, opr, state)
Perform out += opr * state
. Apply the operator representation opr
to the column vector state
and add it to the column vector out
. Return sum of errors and sum of error-squared. Call apply_serial!
if Threads.nthreads() == 1
, and apply_parallel!
otherwise.
QuantumHamiltonian.apply!
— Methodapply!(out, nullop, psi)
Apply operator to psi
and add it to out
.
QuantumHamiltonian.apply_parallel!
— Methodapply_parallel!(out, state, opr)
Perform out += state * opr
. Apply the operator representation opr
to the matrix state
, whose rows are vectors, and add it to the rows of the matrix out
. Multi-threaded version.
QuantumHamiltonian.apply_parallel!
— Methodapply_parallel!(out, opr, state)
Perform out += opr * state
. Apply the operator representation opr
to the matrix state
, whose columns are vectors, and add it to the columns of the matrix out
. Multi-threaded version.
QuantumHamiltonian.apply_parallel!
— Methodapply_parallel!(out, opr, state)
Perform out += opr * state
. Apply the operator representation opr
to the column vector state
and add it to the column vector out
. Multi-threaded version.
QuantumHamiltonian.apply_parallel!
— Methodapply_parallel!(out, state, opr)
Perform out += state * opr
. Apply the operator representation opr
to the row vector state
and add it to the row vector out
. Return sum of errors and sum of error-squared. Multi-threaded version.
QuantumHamiltonian.apply_serial!
— Methodapply_serial!(out, state, opr)
Perform out += state * opr
. Apply the operator representation opr
to the matrix state
, whose rows are vectors, and add it to the rows of the matrix out
. Single-threaded version.
QuantumHamiltonian.apply_serial!
— Methodapply_serial!(out, opr, state)
Perform out += opr * state
. Apply the operator representation opr
to the matrix state
, whose columns are vectors, and add it to the columns of the matrix out
. Single-threaded version.
QuantumHamiltonian.apply_serial!
— Methodapply_serial!(out, opr, state)
Perform out += opr * state
. Apply the operator representation opr
to the column vector state
and add it to the column vector out
. Single-threaded version.
QuantumHamiltonian.apply_serial!
— Methodapply_serial!(out, state, opr)
Perform out += state * opr
. Apply the operator representation opr
to the row vector state
and add it to the row vector out
. Return sum of errors and sum of error-squared. Single-threaded version.
QuantumHamiltonian.basespace
— Functionbasespace(x::AbstractHilbertSpaceRepresentation)
basespace(x::Type{T}) where {T<:AbstractHilbertSpaceRepresentation}
If the argument is an instance of AbstractHilbertSpaceRepresentation
, return the underlying Hilbert space of the Hilbert space representation. If the argument is a subtype of AbstractHilbertSpaceRepresentation
, return the type of the underlying Hilbert space. Subtypes of AbstractHilbertSpace must implement this method.
QuantumHamiltonian.basespace
— Methodbasespace(hss::HilbertSpaceSector, args...;kwargs...)
Call basespace
with basespace of hss
.
QuantumHamiltonian.basespace
— Methodbasespace(hs)
Get the base space of the HilbertSpace hs
, which is itself.
QuantumHamiltonian.bitoffset
— Methodbitoffset(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call bitoffset
with basespace of hsr
.
QuantumHamiltonian.bitoffset
— Methodbitoffset(hss::HilbertSpaceSector, args...;kwargs...)
Call bitoffset
with basespace of hss
.
QuantumHamiltonian.bitwidth
— Methodbitwidth(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call bitwidth
with basespace of hsr
.
QuantumHamiltonian.bitwidth
— Methodbitwidth(hss::HilbertSpaceSector, args...;kwargs...)
Call bitwidth
with basespace of hss
.
QuantumHamiltonian.bitwidth
— Methodbitwidth(hs, [isite])
Total number of bits
julia> using QuantumHamiltonian
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site, spin_site,]);
julia> bitwidth(hs)
3
QuantumHamiltonian.bitwidth
— Methodbitwidth(site)
Number of bits necessary to represent the states of the given site.
QuantumHamiltonian.compress
— Methodcompress(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call compress
with basespace of hsr
.
QuantumHamiltonian.compress
— Methodcompress(hss::HilbertSpaceSector, args...;kwargs...)
Call compress
with basespace of hss
.
QuantumHamiltonian.compress
— Methodcompress(bitwidths, data, BR)
Compress data array into a binary integer of type BR.
QuantumHamiltonian.compress
— Methodcompress(hs, indexarray::CartesianIndex, binary_type=UInt)
Convert a cartesian index (a of state) to its binary representation
Examples
julia> using QuantumHamiltonian
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site]);
julia> compress(hs, CartesianIndex(2,2))
0x0000000000000003
QuantumHamiltonian.compress
— Methodcompress(site, state_index, binary_type=UInt) -> binary_type
Get binary representation of the state specified by state_index
. Check bounds 1 <= state_index <= dimension(site)
, and returns binary representation of state_index-1
.
QuantumHamiltonian.extract
— Methodextract(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call extract
with basespace of hsr
.
QuantumHamiltonian.extract
— Methodextract(hss::HilbertSpaceSector, args...;kwargs...)
Call extract
with basespace of hss
.
QuantumHamiltonian.extract
— Methodextract(hs, binrep)
Convert binary representation to an array of indices (of states)
Examples
julia> using QuantumHamiltonian
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site]);
julia> extract(hs, 0x03)
CartesianIndex(2, 2)
QuantumHamiltonian.findindex
— Functionfindindex(a::AbstractIndexedVector{E}, key::E)
Return the index of the element key
if it exists in a
. Otherwise return -1.
QuantumHamiltonian.get_basis_index_amplitude
— Functionget_basis_index_amplitude(hsr::AbstractHilbertSpaceRepresentation, bin::Unsigned)
Return a tuple (index=index, amplitude=amplitude)
that corresponds to the binary bin
, i.e., the index of the basis state that overlaps with bin
and the value of the overlap ⟨b|ϕᵢ⟩. Return (index=-1, amplitude=0)
if bin
is not in the representation. Subtypes of AbstractHilbertSpaceRepresentation must implement this method.
QuantumHamiltonian.get_basis_index_amplitude
— Methodget_basis_index_amplitude(hsr, bvec)
Get the index of the basis state that overlaps with bvec, and the value of the overlap. Currentiy it is guaranteed to be at most one. Returns (i, ⟨b|ϕᵢ⟩). For the unsymmetrized HilbertSpaceRepresentation
, the amplitude is 1 of Int type. If no such basis vector exists, return (-1, 0).
QuantumHamiltonian.get_basis_iterator
— Functionget_basis_iterator(hsr::AbstractHilbertSpaceRepresentation)
Return an iterator of the list of basis binaries. Subtypes of AbstractHilbertSpaceRepresentation must implement this method.
QuantumHamiltonian.get_basis_list
— Functionget_basis_list(hsr::AbstractHilbertSpaceRepresentation)
Return a Vector of the list of basis binaries. Subtypes of AbstractHilbertSpaceRepresentation must implement this method.
QuantumHamiltonian.get_basis_state
— Functionget_basis_state(hsr::AbstractHilbertSpaceRepresentation, index::Integer)
Return the state at index index
. Subtypes of AbstractHilbertSpaceRepresentation must implement this method.
QuantumHamiltonian.get_basis_state
— Methodget_basis_state(hsr, index)
Get the basis state representation at index.
QuantumHamiltonian.get_bitmask
— Methodget_bitmask(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call get_bitmask
with basespace of hsr
.
QuantumHamiltonian.get_column_iterator
— Methodget_column_iterator(op, bcol)
Returns an iterator over the elements of the column corresponding to bit representation bc
.
QuantumHamiltonian.get_column_iterator
— Methodget_column_iterator(opr, icol)
Returns an iterator which generates a list of elements of the column icol
. Each element is represented as (irow, amplitude). May contain duplicates and invalid elements. Invalid elements are represented as (-1, amplitude).
QuantumHamiltonian.get_element
— Methodget_element(opr, irow, icol)
QuantumHamiltonian.get_operator
— Functionget_operator(x::AbstractOperatorRepresentation)
Return the operator of the operator representation. Subclass of AbstractOperatorRepresentation must define this method.
QuantumHamiltonian.get_quantum_number
— Methodget_quantum_number(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call get_quantum_number
with basespace of hsr
.
QuantumHamiltonian.get_quantum_number
— Methodget_quantum_number(hs, rep)
Get the quantum number of rep
, which is either a binary representation, or a CartesianIndex
.
QuantumHamiltonian.get_quantum_number
— Methodget_quantum_number(hss::HilbertSpaceSector, args...;kwargs...)
Call get_quantum_number
with basespace of hss
.
QuantumHamiltonian.get_quantum_number
— Methodget_quantum_number(site, state_index)
Gets the quantum number of state specified by state_index.
QuantumHamiltonian.get_quantum_numbers
— Methodget_quantum_numbers(hs)
Return a sorted list of quantum numbers of the hilbert space hs
.
QuantumHamiltonian.get_row_iterator
— Methodget_row_iterator(op, br)
Returns an iterator over the elements of the row corresponding to bit representation br
.
QuantumHamiltonian.get_row_iterator
— Methodget_row_iterator(opr, irow)
Returns an iterator which generates a list of elements of the row irow
. Each element is represented as (icol, amplitude). May contain duplicates and invalid elements. Invalid elements are represented as (-1, amplitude).
QuantumHamiltonian.get_row_iterator
— Methodget_row_iterator(ropr::ROR, irow_r::Integer)
Get the row iterator for the reduced operator representation
QuantumHamiltonian.get_site
— Methodget_site(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call get_site
with basespace of hsr
.
QuantumHamiltonian.get_site
— Methodget_site(hss::HilbertSpaceSector, args...;kwargs...)
Call get_site
with basespace of hss
.
QuantumHamiltonian.get_space
— Functionget_space(x::AbstractOperatorRepresentation)
Return the Hilbert Space representation on which the operator representation is defined. Subclass of AbstractOperatorRepresentation must define this method.
QuantumHamiltonian.get_state
— Methodget_state(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call get_state
with basespace of hsr
.
QuantumHamiltonian.get_state
— Methodget_state(hs, binrep, isite)
Get the local state at site isite
for the basis state represented by binrep
. Returns an object of type State
QuantumHamiltonian.get_state
— Methodget_state(site, binrep) where {QN, BR<:Unsigned}
Returns the state of site
represented by the bits binrep
.
QuantumHamiltonian.get_state_index
— Methodget_state_index(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call get_state_index
with basespace of hsr
.
QuantumHamiltonian.get_state_index
— Methodget_state_index(hs, binrep, isite)
Get the index of the local state at site isite
for the basis state represented by binrep
.
QuantumHamiltonian.get_state_index
— Methodget_state_index(site, binrep)
Gets the state index of the binary representation. Returns Int(binrep+1)
.
QuantumHamiltonian.get_tag
— Methodget_tag(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call get_tag
with basespace of hsr
.
QuantumHamiltonian.get_tag
— Methodget_tag(hss::HilbertSpaceSector, args...;kwargs...)
Call get_tag
with basespace of hss
.
QuantumHamiltonian.hs_get_basis_list
— Methodhs_get_basis_list(hs, allowed_quantum_numbers, binary_type=UInt)
Generate a basis for the HilbertSpaceSector
.
QuantumHamiltonian.hs_get_basis_list
— Methodhs_get_basis_list(hss, binary_type=UInt)
Generate a basis for the HilbertSpaceSector
.
QuantumHamiltonian.operatortype
— Functionoperatortype(x::AbstractOperatorRepresentation)
Return the type of the operator of the operator representation. Subclass of AbstractOperatorRepresentation must define this method.
QuantumHamiltonian.pure_operator
— Methodpure_operator(hilbert_space, isite, istate_row, istate_col, amplitude=1, [binary_type=UInt])
Creates a pure operator where projection is at one of the sites.
QuantumHamiltonian.qntype
— Methodqntype(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call qntype
with basespace of hsr
.
QuantumHamiltonian.qntype
— Methodqntype(arg::Type{HilbertSpaceSector{QN}})
Returns the quantum number type of the given hilbert space sector type.
QuantumHamiltonian.qntype
— Methodqntype(::Type{State{QN}})
Returns the quantum number type of the given state type.
QuantumHamiltonian.quantum_number_sectors
— Methodquantum_number_sectors(site) -> Vector{QN}
Gets a list of possible quantum numbers as a sorted vector of QN.
QuantumHamiltonian.represent
— Methodrepresent(hs, binary_type=UInt, basis_type=SortedIndexedVector)
Make a HilbertSpaceRepresentation with all the basis vectors of the specified HilbertSpace. This function defaults to represent_array
.
QuantumHamiltonian.represent
— Methodrepresent(hs, basis_list, basis_type=SortedIndexedVector)
Make a HilbertSpaceRepresentation with the provided list of basis vectors. This defaults to represent_array
.
QuantumHamiltonian.represent
— Methodrepresent(hilbert_space_representation, operator)
Create an OperatorRepresentation
of the operator
in the hilbert_space_representation
.
QuantumHamiltonian.represent_array
— Methodrepresent_array(hs, binary_type=UInt)
Make a HilbertSpaceRepresentation with all the basis vectors of the specified HilbertSpace using FrozenSortedArrayIndex
.
QuantumHamiltonian.represent_array
— Methodrepresent_array(hs, basis_list)
Make a HilbertSpaceRepresentation with the provided list of basis vectors using FrozenSortedArrayIndex
.
QuantumHamiltonian.represent_dict
— Methodrepresent_dict(hs, binary_type=UInt)
Make a HilbertSpaceRepresentation with the provided list of basis vectors using Dict{binary_type, Int}
.
QuantumHamiltonian.represent_dict
— Methodrepresent_dict(hs, basis_list)
Make a HilbertSpaceRepresentation with the provided list of basis vectors using Dict
.
QuantumHamiltonian.scalartype
— Methodscalartype([state or type of state])
Return the scalar type of the state.
QuantumHamiltonian.scalartype
— Methodscalartype(lhs::Type{<:AbstractOperator{S}})
Returns the scalar type of the given AbstractOperator.
QuantumHamiltonian.sectorslice
— MethodPredicate on the tags
QuantumHamiltonian.sectorslice
— MethodPredicate on the tags
QuantumHamiltonian.simplify
— Methodsimplify
Simplify the given operator.
QuantumHamiltonian.spacetype
— Functionspacetype(x::AbstractOperatorRepresentation)
spacetype(x::Type{T}) where {T<:AbstractOperatorRepresentation}
Return the type of the Hilbert space representation on which the operator representation is defined. Subclass of AbstractOperatorRepresentation must define this method.
QuantumHamiltonian.splitblock
— Methodsplitblock
Split n into b blocks.
QuantumHamiltonian.splitrange
— Methodsplitrange
QuantumHamiltonian.symmetry_reduce!
— Methodsymmetry_reduce!(out, rhsr, largevector)
Adds and not overwrites.
QuantumHamiltonian.symmetry_reduce
— Methodsymmetry_reduce(rhsr, large_vector)
Reduce a large vector into the reduced hilbert space representation. Simply throw away components that don't fit.
QuantumHamiltonian.symmetry_reduce
— Methodsymmetry_reduce(hsr, symops, amplitudes, bvec; tol=√ϵ)
Returns bᵢ => ⟨B|ϕᵢ⟩, i.e., the basis state (represented by bᵢ, and the amount of that basis state that overlaps with the input. Returns the same amplitude as the get_basis_index_amplitude
of the reduced Hilbert space representation
Basis states: |ϕᵢ⟩ with a representative bᵢ input : |B⟩
QuantumHamiltonian.tagtype
— Methodtagtype(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call tagtype
with basespace of hsr
.
QuantumHamiltonian.uncompress
— Methoduncompress(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call uncompress
with basespace of hsr
.
QuantumHamiltonian.uncompress
— Methoduncompress(hs, binrep)
Convert binary representation to an array of indices (of states)
Examples
julia> using QuantumHamiltonian
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site]);
julia> extract(hs, 0x03)
CartesianIndex(2, 2)
QuantumHamiltonian.update
— Methodupdate(hsr::AbstractHilbertSpaceRepresentation, args...;kwargs...)
Call update
with basespace of hsr
.
QuantumHamiltonian.update
— Methodupdate(hs, binrep, isite, new_state_index)
Update the binary representation of a basis state by changing the state at site isite
to a new local state specified by new_state_index
.