API

ExactDiagonalization.HilbertSpaceType
HilbertSpace{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])
source
ExactDiagonalization.PureOperatorType
PureOperator{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
source
ExactDiagonalization.SiteType
Site{QN}

A site with quantum number type QN.

Examples

julia> using ExactDiagonalization

julia> site = Site([State("Up", 1), State("Dn", -1)]);
source
ExactDiagonalization.StateType
State{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))
source
Base.valtypeMethod
valtype(::Type{HilbertSpace{QN}})

Returns the valtype (scalar type) of the given hilbert space type.

source
Base.valtypeMethod
valtype(arg::Type{HilbertSpaceSector{HS, QN}})

Returns the valtype (scalar type) of the given hilbert space sector type. For HilbertSpaceSector{QN}, it is always Bool.

source
Base.valtypeMethod
valtype(lhs::Type{<:AbstractOperator{S}})

Returns the valtype (scalar type) of the given AbstractOperator.

source
ExactDiagonalization.apply!Method
apply!(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.

source
ExactDiagonalization.apply!Method
apply!(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.

source
ExactDiagonalization.apply_parallel!Method
apply_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.

source
ExactDiagonalization.apply_parallel!Method
apply_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.

source
ExactDiagonalization.apply_serial!Method
apply_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.

source
ExactDiagonalization.apply_serial!Method
apply_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.

source
ExactDiagonalization.bitwidthMethod

Total 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)
3
source
ExactDiagonalization.compressMethod
compress(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))
0x0000000000000003
source
ExactDiagonalization.compressMethod
compress(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.

source
ExactDiagonalization.extractMethod

Convert 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)
source
ExactDiagonalization.get_column_iteratorMethod
get_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).

source
ExactDiagonalization.get_row_iteratorMethod
get_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).

source
ExactDiagonalization.get_stateMethod
get_state(hs, binrep, isite)

Get the local state at site isite for the basis state represented by binrep. Returns an object of type State

source
ExactDiagonalization.pure_operatorMethod
pure_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.

source
ExactDiagonalization.representMethod
represent(hilbert_space_representation, operator)

Create an OperatorRepresentation of the operator in the hilbert_space_representation.

source
ExactDiagonalization.scalartypeMethod
scalartype(arg::Type{HilbertSpaceSector{HS, QN}})

Returns the scalar type of the given hilbert space sector type. For HilbertSpaceSector{QN}, it is always Bool.

source
ExactDiagonalization.symmetry_reduceMethod
symmetry_reduce(hsr, lattice, symmetry_irrep_component, complex_type=ComplexF64, tol=√ϵ)

Symmetry-reduce the HilbertSpaceRepresentation using translation group.

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

source
LatticeTools.dimensionMethod
dimension(arg::ReducedHilbertSpaceRepresentation{HSR, BR, C}) -> Int

Dimension of the given reduced hilbert space representation, i.e. number of basis elements.

source