Module Hamiltionians.jl
This module contains definitions of Hamiltonians, in particular specific physical models of interest. These are organised by means of an interface around the abstract type AbstractHamiltonian
, in the spirit of the AbstractArray
interface as discussed in the Julia Documentation.
Rimu.Hamiltonians
— Modulemodule Hamiltonians
This module defines Hamiltonian types, interfaces, and functions for working with Hamiltonians.
Relation to other parts of the Rimu
code
In order to define a specific model Hamiltonian with relevant parameters for the model, instantiate the model like this in the input file:
ham = HubbardReal1D(BoseFS((1,2,0,3)); u=1.0, t=1.0)
In the rest of the Rimu
code, access to properties and matrix elements of the model are then provided by the following methods:
ham[address1, address2]
: indexing of matrix elements (slow - use with caution)ham * dv
,ham(dv::AbstractDVec)
ormul!(dv1, ham, dv2)
: use as linear operatordiagonal_element(ham, add)
: diagonal matrix elementnum_offdiagonals(ham, add)
: number of off-diagonalsget_offdiagonal(ham, add, chosen)
: access off-diagonal matrix elementoffdiagonals(ham, add)
: iterator over off-diagonal matrix elementsrandom_offdiagonal(hops)
: choose random off-diagonaldimension(T, ham)
: dimension of linear spacenear_uniform(ham)
: configuration with particles spread across modesstarting_address(ham)
: address for accessing one of the diagonal elements ofham
Model Hamiltonians
Here is a list of fully implemented model Hamiltonians. So far there are two variants implemented of the one-dimensional Bose-Hubbard model real space as well as a momentum-space Hubbard chain.
Rimu.Hamiltonians.HubbardReal1D
— TypeHubbardReal1D(address; u=1.0, t=1.0)
Implements a one-dimensional Bose Hubbard chain in real space.
\[\hat{H} = -t \sum_{\langle i,j\rangle} a_i^† a_j + \frac{u}{2}\sum_i n_i (n_i-1)\]
Arguments
address
: the starting address, defines number of particles and sites.u
: the interaction parameter.t
: the hopping strength.
See also
Rimu.Hamiltonians.ExtendedHubbardReal1D
— TypeExtendedHubbardReal1D(address; u=1.0, v=1.0, t=1.0)
Implements the extended Hubbard model on a one-dimensional chain in real space.
\[\hat{H} = -t \sum_{\langle i,j\rangle} a_i^† a_j + \frac{u}{2}\sum_i n_i (n_i-1) + v \sum_{\langle i,j\rangle} n_i n_j\]
Arguments
address
: the starting address.u
: on-site interaction parameterv
: the next-neighbor interactiont
: the hopping strength
Rimu.Hamiltonians.HubbardMom1D
— TypeHubbardMom1D(address; u=1.0, t=1.0)
Implements a one-dimensional Bose Hubbard chain in momentum space.
\[\hat{H} = \sum_{k} ϵ_k n_k + \frac{u}{M}\sum_{kpqr} a^†_{r} a^†_{q} a_p a_k δ_{r+q,p+k}\\ ϵ_k = -2t \cos(k)\]
Arguments
address
: the starting address, defines number of particles and sites.u
: the interaction parameter.t
: the hopping strength.
See also
Rimu.Hamiltonians.BoseHubbardReal1D2C
— TypeBoseHubbardReal1D2C(address::BoseFS2C; ua=1.0, ub=1.0, ta=1.0, tb=1.0, v=1.0)
Implements a two-component one-dimensional Bose Hubbard chain in real space.
\[\hat{H} = \hat{H}_a + \hat{H}_b + V\sum_{i} n_{a_i}n_{b_i}\]
Arguments
address
: the starting address, defines number of particles and sites.ua
: the on-site interaction parameter parameter for Hamiltonian a.ub
: the on-site interaction parameter parameter for Hamiltonian b.ta
: the hopping strength for Hamiltonian a.tb
: the hopping strength for Hamiltonian b.v
: the inter-species interaction parameter V.
See also
Rimu.Hamiltonians.BoseHubbardMom1D2C
— TypeBoseHubbardMom1D2C(add::BoseFS2C; ua=1.0, ub=1.0, ta=1.0, tb=1.0, v=1.0)
Implements a two-component one-dimensional Bose Hubbard chain in momentum space.
\[\hat{H} = \hat{H}_a + \hat{H}_b + \frac{V}{M}\sum_{kpqr} b^†_{r} a^†_{q} b_p a_k δ_{r+q,p+k}\]
Arguments
add
: the starting address.ua
: theu
parameter for Hamiltonian a.ub
: theu
parameter for Hamiltonian b.ta
: thet
parameter for Hamiltonian a.tb
: thet
parameter for Hamiltonian b.v
: the inter-species interaction parameter V.
See also
Rimu.Hamiltonians.HubbardRealSpace
— TypeHubbardRealSpace(address; u=ones(C, C), t=ones(C), geometry=PeriodicBoundaries(M,))
Hubbard model in real space. Supports single or multi-component Fock state addresses (with C
components) and various (rectangular) lattice geometries in arbitrary dimensions.
\[ \hat{H} = -\sum_{\langle i,j\rangle,σ} t_σ a^†_{iσ} a_{jσ} + \frac{1}{2}\sum_{i,σ} u_{σ,σ} n_{iσ} (n_{iσ} - 1) + \frac{1}{2}\sum_{i,σ≠τ}u_{σ,τ} n_{iσ} n_{iτ}\]
Address types
BoseFS
: Single-component Bose-Hubbard model.FermiFS
: Single-component Fermi-Hubbard model. This address only provides a single species of (non-interacting) fermions. You probably want to useCompositeFS
.CompositeFS
: For multi-component models.
Geometries
Other parameters
u
: the on-site interaction parameters. Must be a symmetric matrix.u[i, j]
corresponds to the interaction between thei
-th andj
-th component.u[i, i]
corresponds to the interaction of a component with itself. Note thatu[i,i]
must be zero for fermionic components.t
: the hopping strengths. Must be a vector of lengthC
. Thei
-th element of the vector corresponds to the hopping strength of thei
-th component.
Rimu.Hamiltonians.MatrixHamiltonian
— TypeMatrixHamiltonian(
mat::AbstractMatrix{T};
starting_address::Int = starting_address(mat)
) <: AbstractHamiltonian{T}
Wrap an abstract matrix mat
as an AbstractHamiltonian
object for use with regular Vector
s indexed by integers. Works with stochastic methods of lomc!()
. Optionally, a starting_address
can be provided.
Specialised methods are implemented for sparse matrices of type AbstractSparseMatrixCSC
.
Hamiltonians interface
Behind the implementation of a particular model is a more abstract interface for defining hamiltonians. If you want to define a new model you should make use of this interface. The most general form of a model Hamiltonian should subtype to AbstractHamiltonian
and implement the relevant methods.
Rimu.Hamiltonians.AbstractHamiltonian
— TypeAbstractHamiltonian{T}
Supertype that provides an interface for linear operators over a linear space with scalar type T
that are suitable for FCIQMC. Indexing is done with addresses (typically not integers) from an address space that may be large (and will not need to be completely generated).
AbstractHamiltonian
instances operate on vectors of type AbstractDVec
from the module DictVectors
and work well with addresses of type AbstractFockAddress
from the module BitStringAddresses
. The type works well with the external package KrylovKit.jl.
Methods
Provides:
offdiagonals
: iterator over reachable off-diagonal matrix elementsrandom_offdiagonal
: function to generate random off-diagonal matrix elementdimension
: get the dimension of the address space.H[address1, address2]
: indexing withgetindex()
- mostly for testing purposes*(H, v)
: deterministic matrix-vector multiply.H(v)
: equivalent toH * v
.mul!(w, H, v)
: mutating matrix-vector multiply.dot(x, H, v)
: computex⋅(H*v)
minimizing allocations.
Methods that need to be implemented:
num_offdiagonals(::AbstractHamiltonian, address)
get_offdiagonal(::AbstractHamiltonian, address, chosen::Integer)
diagonal_element(::AbstractHamiltonian, address)
starting_address(::AbstractHamiltonian)
Optional methods to implement:
Rimu.Hamiltonians.offdiagonals
— Functionoffdiagonals(h::AbstractHamiltonian, address)
Return an iterator over reachable off-diagonal matrix elements of type <:AbstractOffdiagonals
. Defaults to returning Offdiagonals(h, a)
See also
julia> addr = BoseFS((3,2,1));
julia> H = HubbardReal1D(addr);
julia> h = offdiagonals(H, addr)
6-element Rimu.Hamiltonians.Offdiagonals{BoseFS{6, 3, BitString{8, 1, UInt8}}, Float64, HubbardReal1D{Float64, BoseFS{6, 3, BitString{8, 1, UInt8}}, 1.0, 1.0}}:
(BoseFS{6,3}((2, 3, 1)), -3.0)
(BoseFS{6,3}((2, 2, 2)), -2.449489742783178)
(BoseFS{6,3}((3, 1, 2)), -2.0)
(BoseFS{6,3}((4, 1, 1)), -2.8284271247461903)
(BoseFS{6,3}((4, 2, 0)), -2.0)
(BoseFS{6,3}((3, 3, 0)), -1.7320508075688772)
Rimu.Hamiltonians.diagonal_element
— Functiondiagonal_element(ham, add)
Compute the diagonal matrix element of the linear operator ham
at address add
.
Example
julia> addr = BoseFS((3, 2, 1));
julia> H = HubbardMom1D(addr);
julia> diagonal_element(H, addr)
8.666666666666664
diagonal_element(g::G2Correlator, add::BoseFS2C{NA,NB,M,AA,AB})
The diagonal element in G2Correlator
, where (p-q)=0
, hence it becomes
\[\frac{1}{M}\sum_{k,p=1}^M a^†_{k} b^†_{p} b_p a_k .\]
See also
Rimu.Hamiltonians.starting_address
— Functionstarting_address(h)
Return the starting address for Hamiltonian h
. Part of the AbstractHamiltonian
interface. When called on an AbstractMatrix
return the index of the lowest diagonal element.
Example
julia> addr = BoseFS((3, 2, 1));
julia> H = HubbardMom1D(addr);
julia> addr == starting_address(H)
true
The following functions may be implemented instead of offdiagonals
.
Rimu.Hamiltonians.num_offdiagonals
— Functionnum_offdiagonals(ham, add)
Compute the number of number of reachable configurations from address add
.
Example
julia> addr = BoseFS((3, 2, 1));
julia> H = HubbardMom1D(addr);
julia> num_offdiagonals(H, addr)
10
Rimu.Hamiltonians.get_offdiagonal
— Functionnewadd, me = get_offdiagonal(ham, add, chosen)
Compute value me
and new address newadd
of a single (off-diagonal) matrix element in a Hamiltonian ham
. The off-diagonal element is in the same column as address add
and is indexed by integer index chosen
.
Example
julia> addr = BoseFS((3, 2, 1));
julia> H = HubbardMom1D(addr);
julia> get_offdiagonal(H, addr, 3)
(BoseFS{6,3}((2, 1, 3)), 1.0)
The following functions come with default implementations, but may be customized.
Rimu.Hamiltonians.random_offdiagonal
— Functionrandom_offdiagonal(offdiagonals::AbstractOffdiagonals)
random_offdiagonal(ham::AbstractHamiltonian, add)
Generate a single random excitation, i.e. choose from one of the accessible off-diagonal elements in the column corresponding to address add
of the Hamiltonian matrix represented by ham
. Alternatively, pass as argument an iterator over the accessible matrix elements.
Rimu.Hamiltonians.LOStructure
— TypeHamiltonians.LOStructure(op::AbstractHamiltonian)
Hamiltonians.LOStructure(typeof(op))
LOStructure
speficies properties of the linear operator op
. If a special structure is known this can speed up calculations. Implemented structures are:
Hamiltonians.Hermitian
: The operator is complex and Hermitian or real and symmetric.Hamiltonians.AdjointKnown
: The operator is not Hermitian, but itsadjoint
is implemented.Hamiltonians.AdjointUnknown
:adjoint
for this operator is not implemented.
In order to define this trait for a new linear operator type, define a method for LOStructure(::Type{<:MyNewLOType}) = …
.
Rimu.Hamiltonians.dimension
— Functiondimension(::Type{T}, h)
Return the dimension of Hilbert space as T
. If the result does not fit into T
, return nothing
. If T<:AbstractFloat
, an approximate value computed with the improved Stirling formula may be returned instead.
Examples
julia> dimension(HubbardMom1D(BoseFS((1,2,3))))
28
julia> dimension(HubbardMom1D(near_uniform(BoseFS{200,100})))
julia> dimension(Float64, HubbardMom1D(near_uniform(BoseFS{200,100})))
1.3862737677578234e81
julia> dimension(BigInt, HubbardMom1D(near_uniform(BoseFS{200,100})))
1386083821086188248261127842108801860093488668581216236221011219101585442774669540
Geometry
Rimu.Hamiltonians.LatticeGeometry
— Typeabstract type LatticeGeometry
A LatticeGeometry
controls which sites in an AbstractFockAddress
are considered to be neighbours.
Currently only supported by HubbardRealSpace
.
Available implementations
Interface to implement
Base.size
: return the lattice size.neighbour_site(::LatticeGeometry, ::Int, ::Int)
num_neighbours(::LatticeGeometry)
Rimu.Hamiltonians.PeriodicBoundaries
— TypePeriodicBoundaries(size...) <: LatticeGeometry
Rectangular lattice with periodic boundary conditions of size size
.
The dimension of the lattice is controlled by the number of arguments given to its constructor.
This is the default geometry used by HubbardRealSpace
.
Example
julia> lattice = PeriodicBoundaries(5, 4) # 2D lattice of size 5 × 4
PeriodicBoundaries(5, 4)
julia> num_neighbours(lattice)
4
julia> neighbour_site(lattice, 1, 1)
2
julia> neighbour_site(lattice, 1, 2)
5
julia> neighbour_site(lattice, 1, 3)
6
julia> neighbour_site(lattice, 1, 4)
16
See also
Rimu.Hamiltonians.HardwallBoundaries
— TypeHardwallBoundaries
Rectangular lattice with hard wall boundary conditions of size size
. neighbour_site()
will return 0 for some neighbours of boundary sites.
The dimension of the lattice is controlled by the number of arguments given to its constructor.
Example
julia> lattice = HardwallBoundaries(5) # 1D lattice of size 5
HardwallBoundaries(5)
julia> neighbour_site(lattice, 1, 1)
2
julia> neighbour_site(lattice, 1, 2)
0
See also
Rimu.Hamiltonians.LadderBoundaries
— TypeLadderBoundaries(size...; subgeometry=PeriodicBoundaries) <: LatticeGeometry
Lattice geometry where the first dimension is of size 2 and has hardwall boundary conditions. Using this geometry is more efficient than using HardwallBoundaries
with a size of 2, as it does not generate rejected neighbours.
In other dimensions, it behaves like its subgeometry, which can be any LatticeGeometry
.
Example
julia> lattice = LadderBoundaries(2, 3, 4) # 3D lattice of size 2 × 3 × 4
LadderBoundaries(2, 3, 4)
julia> num_neighbours(lattice)
5
julia> neighbour_site(lattice, 1, 1)
2
julia> neighbour_site(lattice, 1, 2)
3
julia> neighbour_site(lattice, 1, 3)
5
julia> neighbour_site(lattice, 1, 4)
7
julia> neighbour_site(lattice, 1, 5)
19
See also
Rimu.Hamiltonians.num_neighbours
— Functionnum_neighbours(geom::LatticeGeometry)
Return the number of neighbours each lattice site has in this geometry.
Note that for efficiency reasons, all sites are expected to have the same number of neighbours. If some of the neighbours are invalid, this is handled by having neighbour_site
return 0.
Rimu.Hamiltonians.neighbour_site
— Functionneighbour_site(geom::LatticeGeometry, site, i)
Find the i
-th neighbour of site
in the geometry. If the move is illegal, return 0.