Functions

Docstrings

FermionicHilbertSpaces.AdditiveConstraintType
AdditiveConstraint(allowed_values, subspaces=missing, functions)

Constraint enforcing that the sum of user-specified per-subspace contributions lies in allowed_values.

The constraint is evaluated on the factors used during state generation. For a composite Hilbert space this means atomic_factors(space).

source
FermionicHilbertSpaces.LazyOperatorType
LazyOperator{O,S,T,C,P}

A matrix-free (lazy) representation of a symbolic operator on a Hilbert space. Acts directly on vectors and matrices via without constructing a sparse matrix.

Constructed via matrix_representation(op, space, :lazy). LazyOperator conforms to the SciMLOperators.AbstractSciMLOperator interface.

source
FermionicHilbertSpaces.NumberConservationType
NumberConservation(total=missing, subspaces=missing, weights=missing)

Constraint enforcing conservation of a (possibly weighted) particle number. total can be a single value or collection of allowed values.

source
FermionicHilbertSpaces.PassthroughType

A source factor that goes unchanged into exactly one target slot. The trivial inner mapper is kept only for phase-factor computation; split/combine bypass it entirely.

source
FermionicHilbertSpaces.PieceIndexType
PieceIndex{block, slot}

A routing index carried in the type domain. getpiece(blocks, ref) extracts blocks[block][slot] with compile-time-known indices, so heterogeneous tuples can be indexed without type instability.

source
FermionicHilbertSpaces.ProductSpaceMapperType
ProductSpaceMapper(source::ProductSpace, targets)

Mapper between source and the partition (or sub-collection) targets.

Fields:

  • leaves : per source factor, an Uncovered, Passthrough or Fractured leaf describing how that factor fractures.
  • factor_routes : per source factor, one PieceIndex{target, group_slot} per piece (in the factor's piece order). Drives combine.
  • target_routes : per target, one PieceIndex{factor, piece_slot} per group (in the target's group order). Drives split.
  • target_spaces : the targets themselves.

factor_routes and target_routes are the two directional views of the same piece table; they are exact inverses of each other by construction.

source
FermionicHilbertSpaces._combined_sector_constraintMethod
_combined_sector_constraint(spaces)

Given a collection of Hilbert spaces, construct a combined constraint from the stored constraints of any SectorHilbertSpace inputs (localized to their respective subspaces). Returns missing if no usable sector information exists.

source
FermionicHilbertSpaces.apply_local_operatorMethod
apply_local_operator(op::SpinSym, state::SpinState) -> (newstate, amplitude)

Apply a single spin operator to a spin state. Returns (newstate, amplitude) where amplitude is 0 if the operation is not allowed, and a nonzero value otherwise.

source
FermionicHilbertSpaces.bdg_hilbert_spaceMethod
bdg_hilbert_space(labels)

This hilbert space uses Nambu states to describe non-interacting systems with superconductive pairing. Matrix representations of quadratic operators take the form
[H Δ; -Δ* -H*]
where H is hermitian and Δ is antisymmetric.

source
FermionicHilbertSpaces.constrain_spaceFunction
constrain_space(space, constraint; kwargs...)
constrain_space(space, states)

Build a constrained Hilbert space from space, either by applying a constraint or by explicitly providing a list of allowed basis states.

source
FermionicHilbertSpaces.embedMethod
embed(m, Hsub => H; complement=complementary_subsystem(H, Hsub), kwargs...)

Compute the embedding of a matrix m in the basis Hsub into the basis H. Fermionic phase factors are included if the two spaces are fermionic Hilbert spaces.

source
FermionicHilbertSpaces.embedMethod
embed(Hsub => H; kwargs...)

Compute the embedding map from Hsub into H. Fermionic phase factors are included if the two spaces are fermionic Hilbert spaces.

source
FermionicHilbertSpaces.fermion_sparse_matrixMethod
fermion_sparse_matrix(fermion_number, H::AbstractHilbertSpace)

Constructs a sparse matrix of size representing a fermionic annihilation operator at bit position fermion_number on the Hilbert space H.

source
FermionicHilbertSpaces.generalized_kronFunction
generalized_kron(ms, Hs, H::AbstractHilbertSpace=tensor_product(Hs))

Compute the tensor product of matrices or vectors in ms with respect to the spaces Hs, respectively. Return a matrix in the space H, which defaults to the tensor_product product of Hs.

source
FermionicHilbertSpaces.generate_statesMethod
generate_states(spaces, constraints; partial_processor=nothing, process_result=(state, space) -> state)

Generate all tensor product states from spaces satisfying constraint. Uses backtracking with pruning via valid_branch.

partial_processor(partial_state, depth, spaces) is called whenever a branch is accepted. process_result(full_state, spaces) can transform each completed state before storing it.

source
FermionicHilbertSpaces.localize_constraintMethod
localize_constraint(c, H::SectorHilbertSpace)

Return a version of constraint c that is bound to H as its subspace, so that when sector_function is later called on a larger combined space it knows to extract the contribution from H only. Already-localized (non-Missing subspace) or unsupported constraints are returned unchanged.

source
FermionicHilbertSpaces.matrix_representationFunction
matrix_representation(op, space::AbstractHilbertSpace; kwargs...)

Return the matrix representation of symbolic operator op in Hilbert space space.

Keyword arguments include projection (default false) which if true will ignore the error thrown when the operator maps outside the space. Use chunking to control NCAdd construction strategy: NoChunking(), TermChunking(scheduler), or StateChunking(scheduler).

Examples

@fermions f
H = hilbert_space(f, 1:2)
op = f[1]' * f[2] + hc
M = matrix_representation(op, H)
size(M) == (dim(H), dim(H))
source
FermionicHilbertSpaces.partial_trace!Function
partial_trace!(mout, m, H::AbstractHilbertSpace, Hsub::AbstractHilbertSpace, complement, extend_state=StateExtender((Hsub, complement), H); skipmissing=true, phase_factors=true)

Compute the partial trace of m from H to Hsub.

source
FermionicHilbertSpaces.partial_traceMethod
partial_trace(m, H => Hsub; complement=complementary_subsystem(H, Hsub))

Compute the partial trace of m from H to Hsub. Fermionic phase factors are included if both H and Hsub are Fermionic, unless specified otherwise in kwargs.

source
FermionicHilbertSpaces.partial_traceMethod
partial_trace(H => Hsub; kwargs...)

Compute the partial trace map from H to Hsub, represented by a sparse matrix of dimension dim(Hsub)^2 x dim(H)^2 that can be multiplied with a vectorized density matrix.

source
FermionicHilbertSpaces.permutation_projectorMethod
permutation_projector(H, Hs, perms, weights=nothing, T = Float64; normalize=true)

Build a group-averaged projector/operator from permutation operators: P = sum(weights[g] * R_g for g).

Hs may be either a full partition of H (covering all atomic factors) or a proper subsystem list (covering only a subset of atomic factors). In the subsystem case the projector is first constructed on Hsub = tensor_product(Hs...) and then embedded back into H via embed(Psub, Hsub => H).

source
FermionicHilbertSpaces.removefermionMethod
removefermion(digitposition, f::FockNumber)

Return (newfocknbr, fermionstatistics) where newfocknbr is the state obtained by removing a fermion at digitposition from f and fermionstatistics is the phase from the Jordan-Wigner string, or 0 if the operation is not allowed.

source
FermionicHilbertSpaces.representationMethod
representation(op_or_state, space::AbstractHilbertSpace; kwargs...)

Return a concrete representation of op_or_state in Hilbert space space.

  • For symbolic operators, returns a sparse matrix (or lazy operator when lazy=true).
  • For AbstractBasisState, returns a one-hot sparse column vector.
  • For SymbolicState (ket/bra), returns the corresponding column/row vector.
source
FermionicHilbertSpaces.sectorFunction
sector(qn, H)

Return the sector of H corresponding to quantum number qn. For qn === nothing, this returns H for non-sector spaces.

source
FermionicHilbertSpaces.single_particle_hilbert_spaceMethod
single_particle_hilbert_space(labels)

A hilbert space suitable for non-interacting systems with fermion number conservation. Matrix representations of symbolic operators give the single particle hamiltonian, without any contribution from the identity matrix.

source
FermionicHilbertSpaces.space_indicesMethod
indices(Hsub, H)
indices(qn, H)

Return the basis-state indices in H belonging to a given sector, specified either by a sector Hilbert space Hsub or by a quantum number qn.

source
FermionicHilbertSpaces.split_stateMethod
split_state(state, mapper)

Split a state using mapper.

Contract: returns a tuple with one entry per target subsystem. Each entry is a weighted collection ((substate, weight), ...) for that target.

source
FermionicHilbertSpaces.state_mapperMethod
state_mapper(H, Hs)

Create a mapper object for decomposing states in H into the target subsystems Hs. The returned object must subtype AbstractStateMapper.

source
FermionicHilbertSpaces.subregionMethod
subregion(Hs, H::AbstractHilbertSpace)

Return the subsystem of H spanned by the factors Hs.

This is primarily used for product spaces where the subsystem is specified by a list/tuple of factor spaces.

Examples

H1 = hilbert_space(1:1)
H2 = hilbert_space(2:2)
H3 = hilbert_space(3:3)
H = tensor_product((H1, H2, H3))
Hsub = subregion((H1, H3), H)
source
FermionicHilbertSpaces.symmetric_sectorMethod
symmetric_sector(H, Hs, sector=:symmetric, T=Float64; normalize=true, cutoff=0.9, atol=1e-10)

Construct the projector onto symmetric sectors of H under permutations of the partition Hs.

If Combinatorics.jl is loaded, sector can be :symmetric or :antisymmetric to automatically generate the full set of permutations and weights for the symmetric or antisymmetric representation. Otherwise, sector must be a tuple of (perms, weights) to specify the permutations and their corresponding weights directly.

Hs may be either a full partition of H or a proper subsystem list; in the latter case the sector is constructed on tensor_product(Hs...) and embedded back into H.

source
FermionicHilbertSpaces.tensor_productMethod
tensor_product(ms, Hs, H::AbstractHilbertSpace; kwargs...)

Compute the ordered product of the fermionic embeddings of the matrices ms in the spaces Hs into the space H. kwargs can be passed a bool phase_factors.

source
FermionicHilbertSpaces.togglefermionsMethod
togglefermions(digitpositions, daggers, f::FockNumber)

Return (newfocknbr, fermionstatistics) where newfocknbr is the state obtained by toggling fermions at digitpositions with daggers in the Fock state f, and fermionstatistics is the phase from the Jordan-Wigner string. If the operation puts two fermions one the same site, the resulting state is undefined.

source
FermionicHilbertSpaces.valid_branchMethod
valid_branch(constraint, partial_state, remaining_spaces) -> Bool

Return `true` if the branch should be explored, `false` to prune. By default this calls `constraint.f(partial_state, remaining_spaces)`.
source
FermionicHilbertSpaces.vector_representationMethod
vector_representation(op::NCMul, space) -> vector or adjoint-vector

For a product of ket SymbolicStates, apply each ket's operator action in sequence and return a sparse column vector. For bras, return the adjoint row vector.

source
FermionicHilbertSpaces.@bosonsMacro
@bosons b

Create a bosonic basis b. Index into it to get bosonic mode operators: b[i] returns the annihilation operator for mode i, b[i]' for the creation operator.

source
FermionicHilbertSpaces.@fermionsMacro
@fermions a b ...

Create one or more fermion species with the given names. Indexing into fermions species gives a concrete fermion. Fermions in one @fermions block anticommute with each other, and commute with fermions in other @fermions blocks.

Examples:

  • @fermions a b creates two species of fermions that anticommute:
    • a[1]' * a[1] + a[1] * a[1]' == 1
    • a[1]' * b[1] + b[1] * a[1]' == 0
  • @fermions a; @fermions b creates two species of fermions that commute with each other:
    • a[1]' * a[1] + a[1] * a[1]' == 1
    • a[1] * b[1] - b[1] * a[1] == 0

See also @majoranas.

source
FermionicHilbertSpaces.@majoranasMacro
@majoranas a b ...

Create one or more Majorana species with the given names. Indexing into Majorana species gives a concrete Majorana. Majoranas in one @majoranas block anticommute with each other, and commute with Majoranas in other @majoranas blocks.

Examples:

  • @majoranas a b creates two species of Majoranas that anticommute:
    • a[1] * a[1] + a[1] * a[1] == 1
    • a[1] * b[1] + b[1] * a[1] == 0
  • @majoranas a; @majoranas b creates two species of Majoranas that commute with each other:
    • a[1] * a[1] + a[1] * a[1] == 1
    • a[1] * b[1] - b[1] * a[1] == 0

See also @fermions.

source