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]
Installation
using Pkg; Pkg.add(url="https://github.com/cvsvensson/FermionicHilbertSpaces.jl")
or by adding a registry to your julia environment and then installing the package
using Pkg; Pkg.Registry.add(RegistrySpec(url = "https://github.com/williamesamuelson/PackageRegistry"))
Pkg.add("FermionicHilbertSpaces")
Introduction
The following example demonstrates how to define a fermionic Hilbert space, create fermionic operators, and construct a simple Hamiltonian:
using FermionicHilbertSpaces, LinearAlgebra
N = 2 # number of fermions
spatial_labels = 1:N
internal_labels = (:↑,:↓)
labels = Base.product(spatial_labels, internal_labels)
H = hilbert_space(labels)
16⨯16 SimpleFockHilbertSpace:
modes: [(1, :↑), (2, :↑), (1, :↓), (2, :↓)]
c = fermions(H) # fermionic annihilation operators
OrderedCollections.OrderedDict{Tuple{Int64, Symbol}, SparseArrays.SparseMatrixCSC{Int64, Int64}} with 4 entries:
(1, :↑) => sparse([1, 3, 5, 7, 9, 11, 13, 15], [2, 4, 6, 8, 10, 12, 14, 16], …
(2, :↑) => sparse([1, 2, 5, 6, 9, 10, 13, 14], [3, 4, 7, 8, 11, 12, 15, 16], …
(1, :↓) => sparse([1, 2, 3, 4, 9, 10, 11, 12], [5, 6, 7, 8, 13, 14, 15, 16], …
(2, :↓) => sparse([1, 2, 3, 4, 5, 6, 7, 8], [9, 10, 11, 12, 13, 14, 15, 16], …
Define a simple Hamiltonian from the fermionic operators
H_hopping = c[1,:↑]'c[2,:↑] + c[1,:↓]'c[2,:↓] + hc
H_coulomb = sum(c[n,:↑]'c[n,:↑]c[n,:↓]'c[n,:↓] for n in spatial_labels)
H_hopping + H_coulomb
16×16 SparseArrays.SparseMatrixCSC{Int64, Int64} with 23 stored entries:
⎡⠠⠂⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠰⢂⠑⢄⠀⠀⎥
⎢⠀⠀⠑⢄⠠⢆⠀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠰⢆⎦
Defining a symbolic hamiltonian
@fermions f
matrix_representation(f[1,:↑]'*f[1,:↑], H)
16×16 SparseArrays.SparseMatrixCSC{Int64, Int64} with 8 stored entries:
⎡⠐⢀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠐⢀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠐⢀⠀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠐⢀⎦
Tensor product and partial trace
H1 = hilbert_space(1:2)
H2 = hilbert_space(3:4)
H = tensor_product(H1, H2)
c1,c2,c = fermions(H1), fermions(H2), fermions(H)
c1c3 = tensor_product([c1[1], c2[3]], [H1, H2] => H)
c[1]*c[3] == c1c3
true
partial_trace(tensor_product([c1[1], I/4], [H1, H2] => H), H => H1) == c1[1]
true
Conserved quantum numbers
H = hilbert_space([1,2], ParityConservation())
4⨯4 SymmetricFockHilbertSpace:
modes: [1, 2]
ParityConservation([-1, 1])
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