FermionicHilbertSpaces.jl

FermionicHilbertSpaces.jl is a Julia package for defining fermionic Hilbert spaces and operators. The central features are fermionic tensor products and partial traces, which differ from the standard tensor product since fermions anticommute. [1]

Introduction

The following example demonstrates how to define a fermionic Hilbert space and construct a simple Hamiltonian.

using FermionicHilbertSpaces, LinearAlgebra
@fermions c
labels = [(i, σ) for i in 1:2 for σ in (:↑,:↓)]
H = hilbert_space(c, labels)
16-dimensional FermionicSpace
Modes: (c[(1, :↑), (1, :↓), (2, :↑), (2, :↓)])

We now have a Hilbert space representing 2N fermions. To construct operators on this space we first call @fermions c which makes c represent symbolic fermions. Indexing into c returns an operator representing an annihilation operator e.g. c[1,:↑]. The creation operator is given by the adjoint c[1,:↑]'. These can be multiplied and added together. Here is a simple Hamiltonian with hopping and Coulomb interaction.

hopping = c[1,:↑]'c[2,:↑] + c[1,:↓]'c[2,:↓] + hc
coulomb = sum(c[n,:↑]'c[n,:↑]c[n,:↓]'c[n,:↓] for n in 1:2)
ham = hopping + coulomb
Sum with 6 terms: 
c†[2, :↑]*c[1, :↑]+c†[2, :↓]*c[1, :↓]-c†[2, :↑]*c†[2, :↓]*c[2, :↑]*c[2, :↓] + ...

To get the matrix representation of this operator on the Hilbert space, do

mat = matrix_representation(ham, H)
16×16 SparseArrays.SparseMatrixCSC{Int64, Int64} with 23 stored entries:
⎡⠀⢀⠂⡀⢄⠀⠀⠀⎤
⎢⠈⠠⠀⢀⠀⠀⢄⠀⎥
⎢⠀⠑⠀⠀⠀⢀⠂⡀⎥
⎣⠀⠀⠀⠑⠈⠠⠑⢄⎦

Tensor product and partial trace

This package also includes functionality for combining Hilbert spaces and operators on them, and taking partial traces, in a way that is consistent with fermionic anticommutation relations.

H1 = hilbert_space(c, 1:2)
H2 = hilbert_space(c, 3:4)
H = tensor_product(H1, H2)
c1 = matrix_representation(c[1], H1)
c3 = matrix_representation(c[3], H2)
c1c3 = matrix_representation(c[1] * c[3], H)
# Use embed to embed operators into a larger space
embed(c1, H1 => H) * embed(c3, H2 => H) ≈ c1c3 #true
# Or call tensor_product to combine operators from different spaces
tensor_product((c1, c3), (H1, H2) => H) ≈ c1c3
true

Let's partial trace to sites 1 and 3. Let's get a Hilbert space for those sites by using the function subregion, and then we can partial trace to that space with partial_trace.

Hsub = subregion([c[1], c[3]], H)
dim(Hsub) / dim(H) * partial_trace(c1c3, H => Hsub) ≈ matrix_representation(c[1] * c[3], Hsub)
true

Conserved quantum numbers

This package also includes some functionality for working with conserved quantum numbers. If we have for example number conservation, we might want to get a block structure of the hamiltonian. Here's how one can do that:

H = hilbert_space(c, labels, NumberConservation())
matrix_representation(ham, H)
16×16 SparseArrays.SparseMatrixCSC{Int64, Int64} with 23 stored entries:
⎡⢀⠐⠄⠀⠀⠀⠀⠀⎤
⎢⠀⠁⢐⠐⠂⡀⠀⠀⎥
⎢⠀⠀⠈⠠⠄⢅⢀⠀⎥
⎣⠀⠀⠀⠀⠀⠐⠕⢅⎦

This has a block structure corresponding to the different sectors. To only look at some sectors, for example the sectors with 0, 2 and 4 particles, use

H = hilbert_space(c, labels, NumberConservation([0, 2, 4]))
matrix_representation(ham, H)
8×8 SparseArrays.SparseMatrixCSC{Int64, Int64} with 11 stored entries:
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  ⋅
 ⋅   1  ⋅  -1  1  ⋅   ⋅  ⋅
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  ⋅
 ⋅  -1  ⋅   ⋅  ⋅  ⋅  -1  ⋅
 ⋅   1  ⋅   ⋅  ⋅  ⋅   1  ⋅
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  ⋅
 ⋅   ⋅  ⋅  -1  1  ⋅   1  ⋅
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  2

Those sectors have even fermionic parity, which can alternatively be specified with ParityConservation.

H = hilbert_space(c, labels, ParityConservation(1))
matrix_representation(ham, H)
8×8 SparseArrays.SparseMatrixCSC{Int64, Int64} with 11 stored entries:
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  ⋅
 ⋅   1  ⋅  -1  1  ⋅   ⋅  ⋅
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  ⋅
 ⋅  -1  ⋅   ⋅  ⋅  ⋅  -1  ⋅
 ⋅   1  ⋅   ⋅  ⋅  ⋅   1  ⋅
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  ⋅
 ⋅   ⋅  ⋅  -1  1  ⋅   1  ⋅
 ⋅   ⋅  ⋅   ⋅  ⋅  ⋅   ⋅  2

References

[1] Szalay, Szilárd, et al. "Fermionic systems for quantum information people." Journal of Physics A: Mathematical and Theoretical 54.39 (2021): 393001, arXiv:2006.03087