API
ExactDiagonalization.AbstractHilbertSpaceRepresentation — TypeAbstractHilbertSpaceRepresentation{S}ExactDiagonalization.AbstractOperator — TypeAbstractOperator{S<:Number}Represent an abstract operator in Hilbert space.
ExactDiagonalization.AbstractOperatorRepresentation — TypeAbstractOperatorRepresentation{S}ExactDiagonalization.FrozenSortedArrayIndex — TypeFrozenSortedArrayIndex(items)An immutable sorted array, with index lookup using binary search.
ExactDiagonalization.HilbertSpace — TypeHilbertSpace{QN}Abstract Hilbert space with quantum number type QN.
Examples
julia> using ExactDiagonalization
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site])
HilbertSpace{Tuple{Int64}}(Site{Tuple{Int64}}[Site{Tuple{Int64}}(State{Tuple{Int64}}[State{Tuple{Int64}}("Up", (1,)), State{Tuple{Int64}}("Dn", (-1,))]), Site{Tuple{Int64}}(State{Tuple{Int64}}[State{Tuple{Int64}}("Up", (1,)), State{Tuple{Int64}}("Dn", (-1,))])], [1, 1], [0, 1, 2])ExactDiagonalization.HilbertSpaceRepresentation — TypeHilbertSpaceRepresentation{HS, BR, DictType}Fields
hilbert_space :: HS
basis_list :: Vector{BR}
basis_lookup :: DictTypeExactDiagonalization.HilbertSpaceSector — TypeHilbertSpaceSector{QN}Hilbert space sector.
ExactDiagonalization.IntegerModulo — TypeIntegerModulo{N}Implement Zₙ.
ExactDiagonalization.NullOperator — TypeNullOperatorA null operator, i.e. 0.
ExactDiagonalization.OperatorRepresentation — TypeOperatorRepresentation{HSR, S, O}Operator representation of given operator of type O.
ExactDiagonalization.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 :: ScalarExactDiagonalization.ReducedHilbertSpaceRepresentation — TypeReducedHilbertSpaceRepresentation{HSR, BR, C}Representation of the symmetry-reduced hilbert space. Currently only supports Translation group (i.e. Abelian group). ```
ExactDiagonalization.ReducedOperatorRepresentation — TypeReducedOperatorRepresentation{RHSR, O, S, BR}Representation of an operator of type O in the symmetry-reduced hilbert space representation of type RHSR.
ExactDiagonalization.Site — TypeSite{QN}A site with quantum number type QN.
Examples
julia> using ExactDiagonalization
julia> site = Site([State("Up", 1), State("Dn", -1)]);ExactDiagonalization.SparseState — Typestruct SparseState{Scalar<:Number, BR}Represents a vector in unrestricted Hilbert space.
ExactDiagonalization.SparseState — MethodSparseState(hsrep, state_rep, tol=√eps(Float64))Make a SparseState from a representation
ExactDiagonalization.State — TypeState{QN}State with quantum number type QN.
Examples
julia> using ExactDiagonalization
julia> up = State("Up", 1)
State{Tuple{Int64}}("Up", (1,))
julia> State("Dn", (-1, 1))
State{Tuple{Int64,Int64}}("Dn", (-1, 1))ExactDiagonalization.SumOperator — TypeSumOperator{Scalar, BR}Represents a sum of pure operators.
Fields
terms::Vector{PureOperator{Scalar,BR}}
Base.valtype — Methodvaltype(::Type{HilbertSpace{QN}})Returns the valtype (scalar type) of the given hilbert space type.
Base.valtype — Methodvaltype(arg::Type{HilbertSpaceSector{HS, QN}})Returns the valtype (scalar type) of the given hilbert space sector type. For HilbertSpaceSector{QN}, it is always Bool.
Base.valtype — Methodvaltype(lhs::Type{<:AbstractOperator{S}})Returns the valtype (scalar type) of the given AbstractOperator.
ExactDiagonalization.apply! — Methodapply!(out, nullop, psi)Apply operator to psi and add it to out.
ExactDiagonalization.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.
ExactDiagonalization.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.
ExactDiagonalization.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.
ExactDiagonalization.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. Return sum of errors and sum of error-squared. Multi-threaded version.
ExactDiagonalization.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.
ExactDiagonalization.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. Return sum of errors and sum of error-squared. Single-threaded version.
ExactDiagonalization.basespace — Methodbasespace(hs)Get the base space of the HilbertSpace hs, which is itself.
ExactDiagonalization.basespace — Methodbasespace(hss)Get the base space of the HilbertSpaceSector, which is its parent HilbertSpace (with no quantum number restriction).
ExactDiagonalization.bitwidth — Methodbitwidth(hss::HilbertSpaceSector, args...;kwargs...)Call bitwidth with basespace of hss.
ExactDiagonalization.bitwidth — MethodTotal number of bits
julia> using ExactDiagonalization
julia> spin_site = Site([State("Up", +1), State("Dn", -1)])
Site{Tuple{Int64}}(State{Tuple{Int64}}[State{Tuple{Int64}}("Up", (1,)), State{Tuple{Int64}}("Dn", (-1,))])
julia> hs = HilbertSpace([spin_site, spin_site, spin_site,])
HilbertSpace{Tuple{Int64}}(Site{Tuple{Int64}}[Site{Tuple{Int64}}(State{Tuple{Int64}}[State{Tuple{Int64}}("Up", (1,)), State{Tuple{Int64}}("Dn", (-1,))]), Site{Tuple{Int64}}(State{Tuple{Int64}}[State{Tuple{Int64}}("Up", (1,)), State{Tuple{Int64}}("Dn", (-1,))]), Site{Tuple{Int64}}(State{Tuple{Int64}}[State{Tuple{Int64}}("Up", (1,)), State{Tuple{Int64}}("Dn", (-1,))])], [1, 1, 1], [0, 1, 2, 3])
julia> bitwidth(hs)
3ExactDiagonalization.bitwidth — Methodbitwidth(site)Number of bits necessary to represent the states of the given site.
ExactDiagonalization.compress — Methodcompress(hss::HilbertSpaceSector, args...;kwargs...)Call compress with basespace of hss.
ExactDiagonalization.compress — Methodcompress(bitwidths, data, BR)Compress data array into a binary integer of type BR.
ExactDiagonalization.compress — Methodcompress(hs, indexarray::CartesianIndex, binary_type=UInt)Convert a cartesian index (a of state) to its binary representation
Examples
julia> using ExactDiagonalization
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site]);
julia> compress(hs, CartesianIndex(2,2))
0x0000000000000003ExactDiagonalization.compress — Methodcompress(site, state_index, binary_type=UInt) -> binary_typeGet binary representation of the state specified by state_index. Check bounds 1 <= state_index <= dimension(site), and returns binary representation of state_index-1.
ExactDiagonalization.extract — Methodextract(hss::HilbertSpaceSector, args...;kwargs...)Call extract with basespace of hss.
ExactDiagonalization.extract — MethodConvert binary representation to an array of indices (of states)
Examples
julia> using ExactDiagonalization
julia> spin_site = Site([State("Up", +1), State("Dn", -1)]);
julia> hs = HilbertSpace([spin_site, spin_site]);
julia> extract(hs, 0x03)
CartesianIndex(2, 2)ExactDiagonalization.get_bitmask — Methodget_bitmask(hss::HilbertSpaceSector, args...;kwargs...)Call get_bitmask with basespace of hss.
ExactDiagonalization.get_column_iterator — Methodget_column_iterator(op, bcol)Returns an iterator over the elements of the column corresponding to bit representation bc.
ExactDiagonalization.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).
ExactDiagonalization.get_element — Methodget_element(opr, irow, icol)ExactDiagonalization.get_quantum_number — Methodget_quantum_number(hss::HilbertSpaceSector, args...;kwargs...)Call get_quantum_number with basespace of hss.
ExactDiagonalization.get_quantum_number — Methodget_quantum_numberExactDiagonalization.get_quantum_number — Methodget_quantum_number(site, state_index)Gets the quantum number of state specified by state_index.
ExactDiagonalization.get_row_iterator — Methodget_row_iterator(op, br)Returns an iterator over the elements of the row corresponding to bit representation br.
ExactDiagonalization.get_row_iterator — Methodget_row_iterator(ropr::ROR, irow_r::Integer)Get the row iterator for the reduced operator representation
ExactDiagonalization.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).
ExactDiagonalization.get_state — Methodget_state(hss::HilbertSpaceSector, args...;kwargs...)Call get_state with basespace of hss.
ExactDiagonalization.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
ExactDiagonalization.get_state — Methodget_state(site, binrep) where {QN, BR<:Unsigned}Returns the state of site represented by the bits binrep.
ExactDiagonalization.get_state_index — Methodget_state_index(hss::HilbertSpaceSector, args...;kwargs...)Call get_state_index with basespace of hss.
ExactDiagonalization.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.
ExactDiagonalization.get_state_index — Methodget_state_index(site, binrep)Gets the state index of the binary representation. Returns Int(binrep+1).
ExactDiagonalization.hs_get_basis_list — Methodhs_get_basis_list(hss, binary_type=UInt)Generate a basis for the HilbertSpaceSector.
ExactDiagonalization.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.
ExactDiagonalization.qntype — Methodqntype(arg::Type{HilbertSpaceSector{QN}})Returns the quantum number type of the given hilbert space sector type.
ExactDiagonalization.qntype — Methodqntype(::Type{HilbertSpace{QN}})Returns the quantum number type of the given hilbert space type.
ExactDiagonalization.qntype — Methodqntype(::Type{Site{QN}})Returns the quantum number type of the given site type.
ExactDiagonalization.qntype — Methodqntype(::Type{State{QN}})Returns the quantum number type of the given state type.
ExactDiagonalization.quantum_number_sectors — Methodquantum_number_sectors(hss::HilbertSpaceSector, args...;kwargs...)Call quantum_number_sectors with basespace of hss.
ExactDiagonalization.quantum_number_sectors — Methodquantum_number_sectorsExactDiagonalization.quantum_number_sectors — Methodquantum_number_sectors(site) -> Vector{QN}Gets a list of possible quantum numbers as a sorted vector of QN.
ExactDiagonalization.represent — Methodrepresent(hs, binary_type=UInt)Make a HilbertSpaceRepresentation with all the basis vectors of the specified HilbertSpace. This function defaults to represent_array.
ExactDiagonalization.represent — Methodrepresent(hs, basis_list)Make a HilbertSpaceRepresentation with the provided list of basis vectors. This defaults to represent_array.
ExactDiagonalization.represent — Methodrepresent(hilbert_space_representation, operator)Create an OperatorRepresentation of the operator in the hilbert_space_representation.
ExactDiagonalization.represent_array — Methodrepresent_array(hs, binary_type=UInt)Make a HilbertSpaceRepresentation with all the basis vectors of the specified HilbertSpace using FrozenSortedArrayIndex.
ExactDiagonalization.represent_array — Methodrepresent_array(hs, basis_list)Make a HilbertSpaceRepresentation with the provided list of basis vectors using FrozenSortedArrayIndex.
ExactDiagonalization.represent_dict — Methodrepresent_dict(hs, binary_type=UInt)Make a HilbertSpaceRepresentation with the provided list of basis vectors using Dict{binary_type, Int}.
ExactDiagonalization.represent_dict — Methodrepresent_dict(hs, basis_list)Make a HilbertSpaceRepresentation with the provided list of basis vectors using Dict.
ExactDiagonalization.scalartype — Methodscalartype(::Type{HilbertSpace{QN}})Returns the scalar type of the given hilbert space type. For HilbertSpace{QN}, it is always Bool.
ExactDiagonalization.scalartype — Methodscalartype([state or type of state])Return the scalar type of the state.
ExactDiagonalization.scalartype — Methodscalartype(arg::Type{HilbertSpaceSector{HS, QN}})Returns the scalar type of the given hilbert space sector type. For HilbertSpaceSector{QN}, it is always Bool.
ExactDiagonalization.scalartype — Methodscalartype(lhs::Type{<:AbstractOperator{S}})Returns the scalar type of the given AbstractOperator.
ExactDiagonalization.simplify — MethodsimplifySimplify the given operator.
ExactDiagonalization.splitblock — MethodsplitblockSplit n into b blocks.
ExactDiagonalization.splitrange — MethodsplitrangeExactDiagonalization.symmetry_reduce! — Methodsymmetry_reduce!(out, rhsr, largevector)Adds and not overwrites.
ExactDiagonalization.symmetry_reduce — Methodsymmetry_reduce(hsr, lattice, symmetry_irrep_component, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group.
ExactDiagonalization.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.
ExactDiagonalization.symmetry_reduce_parallel — Methodsymmetry_reduce_parallel(hsr, trans_group, frac_momentum, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (multi-threaded).
ExactDiagonalization.symmetry_reduce_parallel — Methodsymmetry_reduce_parallel(hsr, trans_group, frac_momentum, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (multi-threaded).
ExactDiagonalization.symmetry_reduce_parallel — Methodsymmetry_reduce_parallel(hsr, trans_group, frac_momentum, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (multi-threaded).
ExactDiagonalization.symmetry_reduce_parallel — Methodsymmetry_reduce_parallel(hsr, symops_and_amplitudes; tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (multi-threaded).
ExactDiagonalization.symmetry_reduce_serial — Methodsymmetry_reduce_serial(hsr, trans_group, frac_momentum, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (single threaded).
ExactDiagonalization.symmetry_reduce_serial — Methodsymmetry_reduce_serial(hsr, trans_group, frac_momentum, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (single threaded).
ExactDiagonalization.symmetry_reduce_serial — Methodsymmetry_reduce_serial(hsr, trans_group, frac_momentum, complex_type=ComplexF64, tol=√ϵ)Symmetry-reduce the HilbertSpaceRepresentation using translation group (single threaded).
ExactDiagonalization.symmetry_reduce_serial — Methodsymmetry_reduce_serial(hilbert_space_representation, symops_and_amplitudes; tol=√ϵ)The irreps have to follow certain order.
ExactDiagonalization.update — Methodupdate(hss::HilbertSpaceSector, args...;kwargs...)Call update with basespace of hss.
ExactDiagonalization.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.
LatticeTools.dimension — MethoddimensionDimension of the Concrete Hilbert space, i.e. number of basis vectors.
LatticeTools.dimension — Methoddimension(arg::ReducedHilbertSpaceRepresentation{HSR, BR, C}) -> IntDimension 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).
LinearAlgebra.mul! — MethodC(i, k) = A(i, j) * B(j, k)
C(:, k) = A(:, :) * B(:, k)