Internals
LatticeModels.AbstractBonds
LatticeModels.AbstractCurrents
LatticeModels.AbstractLattice
LatticeModels.AbstractSite
LatticeModels.AbstractTranslation
LatticeModels.AdjacencyMatrix
LatticeModels.AdjacencyMatrix
LatticeModels.AdjacencyMatrix
LatticeModels.AdjacencyMatrix
LatticeModels.BallND
LatticeModels.BoundaryConditions
LatticeModels.Box
LatticeModels.Bravais
LatticeModels.BravaisSite
LatticeModels.BravaisTranslation
LatticeModels.BravaisTranslation
LatticeModels.CachedExp
LatticeModels.Currents
LatticeModels.Currents
LatticeModels.DensityCurrents
LatticeModels.DirectedBonds
LatticeModels.Eigensystem
LatticeModels.Evolution
LatticeModels.EvolutionSolver
LatticeModels.FunctionBoundary
LatticeModels.GaugeField
LatticeModels.GenericLattice
LatticeModels.GenericLattice
LatticeModels.GenericLattice
LatticeModels.GenericLattice
LatticeModels.GenericSite
LatticeModels.GreenFunction
LatticeModels.Hamiltonian
LatticeModels.HoneycombLattice
LatticeModels.IncompatibleLattices
LatticeModels.KagomeLattice
LatticeModels.KrylovKitExp
LatticeModels.LandauGauge
LatticeModels.LatticeBasis
LatticeModels.LatticeValue
LatticeModels.LatticeValue
LatticeModels.LatticeValue
LatticeModels.LineIntegralGaugeField
LatticeModels.LocalOperatorCurrents
LatticeModels.LookupTable
LatticeModels.NParticles
LatticeModels.NearestNeighbor
LatticeModels.NearestNeighbor
LatticeModels.NoField
LatticeModels.OperatorBuilder
LatticeModels.OperatorBuilder
LatticeModels.Path
LatticeModels.PointFlux
LatticeModels.PointFlux
LatticeModels.PointFluxes
LatticeModels.PointFluxes
LatticeModels.PointFluxes
LatticeModels.Polygon
LatticeModels.SiteAt
LatticeModels.SiteDistance
LatticeModels.SiteProperty
LatticeModels.SquareLattice
LatticeModels.SubCurrents
LatticeModels.SymmetricGauge
LatticeModels.System
LatticeModels.System
LatticeModels.TimeSequence
LatticeModels.TimeSequence
LatticeModels.TimeSequence
LatticeModels.TimeSequence
LatticeModels.Translation
LatticeModels.TriangularLattice
LatticeModels.TwistedBoundary
LatticeModels.UndefinedLattice
LatticeModels.UnitCell
LatticeModels.Chain
LatticeModels.FermiHubbardSpinSystem
LatticeModels.GrapheneRibbon
LatticeModels.PeriodicBoundary
LatticeModels.adapt_bonds
LatticeModels.addlookuptable
LatticeModels.adjacentsites
LatticeModels.bosehubbard
LatticeModels.brastate
LatticeModels.check_issublattice
LatticeModels.check_samelattice
LatticeModels.check_samesites
LatticeModels.construct_hamiltonian
LatticeModels.construct_operator
LatticeModels.coordoperator
LatticeModels.coordoperators
LatticeModels.coordvalue
LatticeModels.coordvalues
LatticeModels.currentsfrom
LatticeModels.currentsfromto
LatticeModels.densitymatrix
LatticeModels.diagonalelements
LatticeModels.diagonalize
LatticeModels.differentiate
LatticeModels.differentiate!
LatticeModels.dos
LatticeModels.fermihubbard
LatticeModels.fillshapes
LatticeModels.findgroundstate
LatticeModels.getnnbonds
LatticeModels.greenfunction
LatticeModels.greenfunction
LatticeModels.groundstate
LatticeModels.haldane
LatticeModels.hubbard
LatticeModels.integrate
LatticeModels.integrate!
LatticeModels.interaction
LatticeModels.interaction
LatticeModels.kanemele
LatticeModels.ketstate
LatticeModels.lattice
LatticeModels.lattice
LatticeModels.ldos
LatticeModels.line_integral
LatticeModels.localdensity
LatticeModels.localexpect
LatticeModels.periodic_fluxes
LatticeModels.projector
LatticeModels.projector
LatticeModels.qwz
LatticeModels.removedangling!
LatticeModels.sample
LatticeModels.setboundaries
LatticeModels.setnnbonds
LatticeModels.shaperadius
LatticeModels.site_index
LatticeModels.sitedistance
LatticeModels.span_unitcells
LatticeModels.tightbinding_hamiltonian
LatticeModels.timerange
LatticeModels.timestamps
LatticeModels.translate_to_nearest
LatticeModels.vector_potential
LatticeModels.@bravaisdef
Lattice basics
LatticeModels.IncompatibleLattices
— TypeIncompatibleLattices([header, ]lat1, lat2)
An exception thrown when two lattices are incompatible.
LatticeModels.SiteProperty
— TypeSiteProperty
An abstract type for a property of a site.
This interface is used to define various properties of a site. They can be accessed using getsiteproperty
. This interface is used in following places:
lattice[...]
syntax to access sites with specific properties.lattice_value[...]
syntax to access values defined on sites with specific properties.- Functions to generate
LatticeValue
s and operators for specific properties.
Examples
julia> using LatticeModels
julia> l = SquareLattice(3, 3);
julia> l[x = 1, y = 2] # Get site with x = 1 and y = 2
2-dim Bravais lattice site in 2D space at [1.0, 2.0]
julia> l[x = 1] # Get sublattice with x = 1
3-site 2-dim Bravais lattice in 2D space
Unit cell:
Basis site coordinates:
⎡ 0.000⎤
⎣ 0.000⎦
Translation vectors:
⎡ 1.000⎤ ⎡ 0.000⎤
⎣ 0.000⎦ ⎣ 1.000⎦
Lattice type: SquareLattice{2}
Default translations:
:axis1 → Bravais[3, 0]
:axis2 → Bravais[0, 3]
Nearest neighbor hoppings:
1.00000 =>
Bravais[1, 0]
Bravais[0, 1]
1.41421 =>
Bravais[1, -1]
Bravais[1, 1]
2.00000 =>
Bravais[2, 0]
Bravais[0, 2]
Boundary conditions: none
julia> l[x = 1, y = 2, z = 3] # No site with defined z property on a 2D lattice
ERROR: ArgumentError: Invalid axis index 3 of a 2-dim site
[...]
LatticeModels.check_issublattice
— MethodChecks if l1
is sublattice of l2
. Throws an error if not.
LatticeModels.check_samelattice
— MethodChecks if l1
and l2
objects are defined on the same lattice. Throws an error if not.
LatticeModels.check_samesites
— MethodChecks if l1
and l2
objects are defined on the same sites. Throws an error if not.
LatticeModels.lattice
— Methodlattice(any)
Return the lattice of the given object (an operator, LatticeValue
, ...)
LatticeModels.site_index
— Functionsite_index(lat, site[, range])
Return the index of the site
in the lattice lat
. If range
is given, only search in the given range. Return nothing
if the site is not found.
LatticeModels.GenericLattice
— TypeGenericLattice{SiteT}
A generic lattice of SiteT
sites.
Example
julia> using LatticeModels
julia> l = GenericLattice{2}()
0-site GenericLattice{GenericSite{2}} in 2D space
julia> push!(l, GenericSite(0, 0)) # Add a site at (0, 0)
1-site GenericLattice{GenericSite{2}} in 2D space:
Site at [0.0, 0.0]
julia> push!(l, (0, 1)) # Add a site at (0, 1)
2-site GenericLattice{GenericSite{2}} in 2D space:
Site at [0.0, 0.0]
Site at [0.0, 1.0]
julia> push!(l, [1, 0]) # Add a site at (1, 0)
3-site GenericLattice{GenericSite{2}} in 2D space:
Site at [0.0, 0.0]
Site at [0.0, 1.0]
Site at [1.0, 0.0]
julia> l[2]
2-dim GenericSite{2} at [0.0, 1.0]
LatticeModels.GenericLattice
— MethodGenericLattice(lat)
Constructs a GenericLattice
from some other lattice lat
.
LatticeModels.GenericLattice
— MethodGenericLattice{N}()
Constructs an empty N
-dimensional GenericLattice
of GenericSite
s.
LatticeModels.GenericLattice
— MethodGenericLattice{SiteType}()
Constructs an empty GenericLattice
of SiteType
sites.
LatticeModels.GenericSite
— TypeGenericSite{N}
A generic site in an N-dimensional lattice.
LatticeModels.BravaisSite
— TypeBravaisSite{N,NU,B}
A site of a BravaisLattice{N,NU,B}
lattice.
Fields
unitcell
: aUnitCell
object representing the lattice unit cell.latcoords
: aSVector
of sizeN
representing the lattice coordinates of the site.basindex
: anInt
representing the index of the site in the lattice basis.coords
: aSVector
of sizeN
representing the spatial coordinates of the site.
LatticeModels.UnitCell
— TypeUnitCell(translations[, basis; offset, rotate])
Constructs a Bravais lattice unit cell with given translation vectors and locations of basis sites.
Arguments
translations
: anAbstractMatrix
of sizeN×N
representing the translation vectors of the lattice.basis
: anAbstractMatrix
of sizeN×NB
representing the locations of basis sites. If not provided, the lattice basis will consist of one site located in the bottom-left corner of the unit cell.
Keyword arguments
offset
: a keyword argument that specifies how to shift the lattice basis. Possible values::origin
: no shift (default).:center
: shift the lattice so that the center of the basis is at the origin of the unit cell.:centeralign
: shift the lattice so that the center of the basis is at the center of the unit cell.- Also accepts an
AbstractVector
of sizeN
to shift the lattice by a custom vector.
rotate
: a keyword argument that specifies how to rotate the lattice basis. Possible values:nothing
: no rotation (default).- An
AbstractMatrix
of sizeN×N
to rotate the lattice. - A
Real
number to rotate the lattice by this angle in radians. - Also accepts an
AbstractMatrix
of sizeN×N
to rotate the lattice basis.
LatticeModels.span_unitcells
— Methodspan_unitcells([f, ]unitcell, dims...[; boundaries, offset])
Construct a Bravais lattice by spanning unitcell
in dims
dimensions, filtered by f
.
Arguments
f
: a function that defines if the site is included in the lattice. Takes aBravaisSite
, returns aBool
.unitcell
: aUnitCell
object.dims
: a list ofInteger
s orRange
s specifying the size of the lattice in each dimension.
Keyword arguments
default_translations
: a list ofBravaisTranslation
s to add to the lattice as default boundary condition axes.boundaries
: aBoundaryConditions
object specifying the boundary conditions of the lattice.rmdup
: aBool
specifying whether to remove sites that are equivalent after applying the boundary conditions.offset
: the offset of the lattice from the origin. SeeUnitCell
for details.rotate
: a rotation matrix to apply to the lattice. SeeUnitCell
for details.
Keep in mind that the offset and rotation are applied to the unit cell before the lattice is spanned (and f
is applied). To apply them after the lattice is spanned, use the postoffset
and postrotate
keywords.
Examples
julia> using LatticeModels
julia> using LatticeModels
julia> uc = UnitCell([[1, 0] [0, 1]])
1-site Unit cell of a 2-dim Bravais lattice in 2D space:
Basis site coordinates:
⎡ 0.000⎤
⎣ 0.000⎦
Translation vectors:
⎡ 1.000⎤ ⎡ 0.000⎤
⎣ 0.000⎦ ⎣ 1.000⎦
julia> span_unitcells(uc, 3, 3) == SquareLattice(3, 3)
true
Lattice constructors
LatticeModels.HoneycombLattice
— TypeHoneycombLattice
Represents a honeycomb lattice.
Lattice vectors: [1, 0]
and [0.5, √3/2]
, two sites at [0, 0]
and [0.5, √3/6]
in each unit cell.
HoneycombLattice(a, b)
Construct a honeycomb lattice of a×b spanned unit cells.
LatticeModels.KagomeLattice
— TypeKagomeLattice
Represents a kagome lattice.
Lattice vectors: [1, 0]
and [0.5, √3/2]
, three sites at [0, 0]
, [0.5, 0]
and [0.25, √3/4]
in each unit cell.
KagomeLattice(a, b)
Construct a kagome lattice of a×b spanned unit cells.
LatticeModels.SquareLattice
— TypeSquareLattice{N}
Represents a square lattice in N
dimensions.
SquareLattice(sz...)
Construct a square lattice of size sz
.
LatticeModels.TriangularLattice
— TypeTriangularLattice
Represents a triangular lattice. Lattice vectors: [1, 0]
and [0.5, √3/2]
.
TriangularLattice(a, b)
Construct a triangular lattice of a×b spanned unit cells.
LatticeModels.Chain
— MethodChain(sz)
Construct a 1D chain lattice of size sz
.
LatticeModels.GrapheneRibbon
— FunctionGrapheneRibbon(len, wid[, center; kw...])
Construct a graphene ribbon sample with zigzag edges. To get armchair edges, simply rotate the lattice by 90 degrees.
Arguments
len
: the length of the ribbon.wid
: the width of the ribbon.center
: the unit cell coordinates of the bottom-left corner of the ribbon. Default is(0, 0)
.
All other keyword arguments are passed to span_unitcells
(see its documentation for details).
LatticeModels.@bravaisdef
— Macro@bravaisdef MyBravaisLattice UnitCell(...)
@bravaisdef MyBravaisLattice N -> UnitCell(...)
Define a new Bravais lattice type MyBravaisLattice
with a unit cell constructor UnitCell(expr)
. If the notation is N -> UnitCell(expr)
, the unit cell constructor will be dependent on the dimensionality N
. otherwise, the dimensionality will be inferred from the unit cell. N
is the dimensionality of the lattice.
Examples
julia> using LatticeModels
julia> @bravaisdef MyBravaisLattice UnitCell([1 0; 0 1]); # 2D square lattice
julia> MyBravaisLattice(3, 3)
9-site 2-dim Bravais lattice in 2D space
Unit cell:
Basis site coordinates:
⎡ 0.000⎤
⎣ 0.000⎦
Translation vectors:
⎡ 1.000⎤ ⎡ 0.000⎤
⎣ 0.000⎦ ⎣ 1.000⎦
Lattice type: MyBravaisLattice
Default translations:
:axis1 → Bravais[3, 0]
:axis2 → Bravais[0, 3]
Nearest neighbor hoppings:
1.00000 =>
Bravais[1, 0]
Bravais[0, 1]
1.41421 =>
Bravais[1, -1]
Bravais[1, 1]
2.00000 =>
Bravais[2, 0]
Bravais[0, 2]
Boundary conditions: none
LatticeModels.BallND
— TypeBallND{N}([radius, center])
Construct a N
-dimensional ball with a given radius and center. Note the aliases: Circle
and Ball
are BallND{2}
and BallND{3}
respectively.
Arguments
radius
: The radius of the ball.center
: The center of the ball.
LatticeModels.Box
— TypeBox(intervals...)
Construct a box with given horizontal and vertical intervals. Usage: Box(1 .. 3, 2 .. 4)
.
Arguments
intervals
: The intervals for each dimension.
LatticeModels.Path
— TypePath(start, stop)
Construct a path from start
to stop
.
Arguments
start
: The start of the path.stop
: The end of the path.
LatticeModels.Polygon
— TypePolygon{N}([radius, center])
Polygon{N}([center; h])
Construct a regular N
-sided polygon with a given (circumscribed) radius and center. Note the aliases: Triangle
, Square
, and Hexagon
are Polygon{3}
, Polygon{4}
, and Polygon{6}
respectively.
Arguments
radius
: The (circumscribed) radius of the polygon.center
: The center of the polygon.
Keyword Arguments
h
: The distance from the center to the vertices. If given, theradius
is calculated ash / cos(pi / N)
.
LatticeModels.SiteAt
— TypeSiteAt(coords)
Represents a single site at the given coordinates.
LatticeModels.fillshapes
— Methodfillshapes(uc, shapes...[; sites, scale, kw...])
Create a lattice sample with geometry defined by the given shapes. The lattice is filled with sites that are inside the shapes.
Arguments
uc
: TheUnitCell
of the lattice. Might also be a lattice type.shapes
: The shapes to fill the lattice with.
Keyword Arguments
sites
: If given, an attepmt will be made to fill the lattice with the given number of sites. The scaling will be approximate and relying on assumptions that the shapes do not overlap.scale
: The scaling factor for the shapes. Ifsites
is given, the scaling factor will be calculated automatically.
All other keyword arguments are passed to the lattice constructor. See span_unitcells
for more information.
LatticeModels.removedangling!
— Functionremovedangling!(lat[, maxdepth])
Remove dangling sites from the lattice. A site is considered dangling if it has less than 2 neighbors. The function will remove all dangling sites and their neighbors recursively up to maxdepth
levels — the default is Inf
.
LatticeModels.shaperadius
— Methodshaperadius(unitcell, shape, sites)
shaperadius(lat, shape[, sites])
Calculate the radius of a shape such that it contains appriximately sites
sites.
Arguments
unitcell
: TheUnitCell
of the lattice. Might also be a lattice type.lat
: The lattice. It is considered that the lattice was constructed in the same shape.shape
: The shape to calculate the radius for.sites
: The number of sites the shape should contain.
Bonds
LatticeModels.AdjacencyMatrix
— TypeAdjacencyMatrix{LT} where {LT<:Lattice}
Represents the bonds on some lattice.
AdjacencyMatrix(lat[, mat])
Construct an adjacency matrix from the mat
matrix on the lat
lattice.
If mat
is not provided, it is assumed to be a zero matrix.
Example
julia> using LatticeModels
julia> l = SquareLattice(2, 2);
julia> a = AdjacencyMatrix(l)
Adjacency matrix on 4-site SquareLattice in 2D space
Values in a 4×4 SparseArrays.SparseMatrixCSC{Bool, Int64} with 0 stored entries:
⋅ ⋅ ⋅ ⋅
⋅ ⋅ ⋅ ⋅
⋅ ⋅ ⋅ ⋅
⋅ ⋅ ⋅ ⋅
julia> site1, site2, site3, site4 = l;
julia> a[site1, site2] = a[site2, site4] = a[site3, site4] = true;
julia> a
Adjacency matrix on 4-site SquareLattice in 2D space
Values in a 4×4 SparseArrays.SparseMatrixCSC{Bool, Int64} with 6 stored entries:
⋅ 1 ⋅ ⋅
1 ⋅ ⋅ 1
⋅ ⋅ ⋅ 1
⋅ 1 1 ⋅
LatticeModels.AdjacencyMatrix
— MethodAdjacencyMatrix(f, lat)
Constructs an adjacency matrix from the function f
that returns if the sites are connected on the lat
lattice.
LatticeModels.AdjacencyMatrix
— MethodAdjacencyMatrix([lat, ]bonds...)
Constructs an adjacency matrix from the bonds
. If lat
is not provided, it is inferred from the bonds
.
LatticeModels.SiteDistance
— TypeSiteDistance(f, lat)
A bonds type that connects sites based on the distance between them.
Arguments
f
: A function that takes a distance and returns if the distance is allowed.lat
: The lattice where the bonds are defined.
LatticeModels.Translation
— TypeTranslation <: AbstractTranslation
A spatial translation on some lattice.
Fields
lat
: The lattice where the translations are defined.R
: The vector of the translation.
Example
julia> using LatticeModels
julia> gl = GenericLattice([(1, 1), (1, 2), (2, 1), (2, 2)])
4-site GenericLattice{GenericSite{2}} in 2D space:
Site at [1.0, 1.0]
Site at [1.0, 2.0]
Site at [2.0, 1.0]
Site at [2.0, 2.0]
julia> tr = Translation(gl, [1, 0]) # Translation by [1, 0]
Translation by [1.0, 0.0]
on 4-site GenericLattice{GenericSite{2}} in 2D space
julia> site1 = gl[!, x = 1, y = 1] # Site at [1, 1]
2-dim GenericSite{2} at [1.0, 1.0]
julia> site1 + tr # Translated site
2-dim GenericSite{2} at [2.0, 1.0]
julia> site1 - tr # Inverse translation
2-dim GenericSite{2} at [0.0, 1.0]
LatticeModels.UndefinedLattice
— TypeUndefinedLattice
A lattice that is not defined. The bonds can be 'defined' on it in context where the lattice is already defined before, e. g. in construct_operator
.
LatticeModels.adapt_bonds
— Methodadapt_bonds(bonds, lat)
Adapt the bonds to the lattice lat
. The output can be a different type of bonds, more fitting for the concrete type of lattice.
LatticeModels.adjacentsites
— Methodadjacentsites(bonds, site)
Returns the sites that are connected to site
by the bonds
.
LatticeModels.sitedistance
— Methodsitedistance([lat, ]site1, site2)
Returns the distance between two sites on the lat
lattice, taking boundary conditions into account.
Arguments
lat
: The lattice where the sites are defined.site1
andsite2
: The sites to measure the distance between.
LatticeModels.translate_to_nearest
— Methodtranslate_to_nearest(lat, site1, site2)
Translate site2
to its equivalent nearest to site1
in the lattice lat
, taking the boundary conditions into account.
LatticeModels.Bravais
— TypeBravais[ lattice_coords ]
A convenient constructor for a BravaisTranslation
that does not permute sublattices.
LatticeModels.BravaisTranslation
— TypeBravaisTranslation([site_indices, ]translate_uc)
BravaisTranslation(site_indices)
BravaisTranslation([site_indices; ]axis[, dist=1])
A convenient constructor for a BravaisTranslation
object.
Arguments
site_indices
: a::Int => ::Int
pair with indices of sites connected by the bond;
if omitted, the bond connects sites with the same sublattice index.
translate_uc
: The unit cell offset.
Keyword arguments
axis
: The hopping direction axis in terms of unit cell vectors.dist
: The hopping distance in terms of unit cell vectors.
If site_indices
are equal or undefined and translate_uc
is zero, the translation is considered to be a translation of all sites to themselves. An error will be thrown in this case.
LatticeModels.BravaisTranslation
— TypeBravaisTranslation{T, N}
A struct representing bonds in some direction in a lattice.
Note that though the dimension count for the bond is static, it is automatically compatible with higher-dimensional lattices.
LatticeModels.NearestNeighbor
— TypeNearestNeighbor(lat[, N=1])
Returns the nearest neighbor bonds of order N
for the lattice lat
.
Example
julia> using LatticeModels
julia> lat = HoneycombLattice(5, 5);
julia> NearestNeighbor(lat)
BravaisSiteMapping with 3 translations:
1 => 2, [0, -1]
1 => 2, [-1, 0]
1 => 2, [0, 0]
on 50-site HoneycombLattice in 2D space
julia> lat = SquareLattice(3, 3, 3, 3);
julia> NearestNeighbor(lat, 4)
BravaisSiteMapping with 12 translations:
Bravais[1, -1, -1, -1]
Bravais[1, 1, -1, -1]
Bravais[1, -1, 1, -1]
Bravais[1, 1, 1, -1]
Bravais[2, 0, 0, 0]
Bravais[0, 2, 0, 0]
Bravais[0, 0, 2, 0]
Bravais[1, -1, -1, 1]
Bravais[1, 1, -1, 1]
⋮
on 81-site SquareLattice in 4D space
LatticeModels.NearestNeighbor
— TypeNearestNeighbor{N}
A bonds type that connects sites that are nearest neighbors of order N
on some lattice.
LatticeModels.getnnbonds
— Methodgetnnbonds(lat)
Returns the nearest neighbor bonds of the lattice lat
.
LatticeModels.setnnbonds
— Methodsetnnbonds(lat, args...; overwrite=false)
Adds the nearest neighbor bonds args
to the lattice lat
. If overwrite
is true
, the default nearest neighbor bonds are replaced by args
. Otherwise, the new bonds are merged with the default.
Each args
can be a bonds type or a distance-bonds pair.
Example
julia> using LatticeModels
julia> lat = SquareLattice(3, 3);
julia> lat2 = setnnbonds(lat, SiteDistance(0 .. 1), SiteDistance(1 .. 2));
julia> lat2.nnbonds
Nearest neighbor hoppings:
#1 =>
SiteDistance(0 .. 1)
#2 =>
SiteDistance(1 .. 2)
Boundary conditions
LatticeModels.BoundaryConditions
— TypeBoundaryConditions
A collection of boundary conditions for a lattice.
Fields
bcs
: A tuple of boundary conditions.depth
: The upper limit of the depth of the boundary conditions (used for routing).
LatticeModels.FunctionBoundary
— TypeFunctionBoundary <: Boundary
A boundary condition with a function that returns the phase factor for a given site. The boundary condition is encoded in form $ψ(x + R) = f(x)ψ(x)$, where $f(x)$ is the function and $R$ is the translation vector.
FunctionBoundary(f, translation)
Construct a FunctionBoundary
with a given function and translation.
Arguments
f
: The function that returns the phase factor for a given site.translation
: The translation vector of the boundary representad asAbstractTranslation
. If an array is passed, it is converted toTranslation
automatically.
LatticeModels.TwistedBoundary
— TypeTwistedBoundary <: Boundary
A boundary condition with a phase twist. A PeriodicBoundary
is a special case of TwistedBoundary
with zero twist.
TwistedBoundary(translation, Θ)
Construct a TwistedBoundary
with a given translation and twist angle.
Arguments
translation
: The translation vector of the boundary representad asAbstractTranslation
. If an array is passed, it is converted toTranslation
automatically.Θ
: The twist angle in radians.
LatticeModels.PeriodicBoundary
— MethodPeriodicBoundary(translation)
Construct a PeriodicBoundary
with a given translation.
LatticeModels.setboundaries
— Methodsetboundaries(lat, boundaries...[; checkboundaries=true, rmdup=false])
Set the boundary conditions for the lattice lat
.
Arguments
lat
: The lattice.boundaries
: The boundary conditions. It can be a singleBoundary
or aTuple
ofBoundary
objects.
Keyword arguments
checkboundaries
: Iftrue
, check if the boundary conditions overlap within the lattice sites.rmdup
: Iftrue
, remove duplicate sites from the lattice.
Example
julia> using LatticeModels
julia> l = SquareLattice(4, 4);
julia> l2 = setboundaries(l, [4, 0] => true, [0, 4] => pi);
julia> l2.boundaries
Boundary conditions (depth = 1):
Bravais[4, 0] → periodic
Bravais[0, 4] → twist θ = 3.14
LatticeValue
LatticeModels.LatticeValue
— TypeLatticeValue{T, LT}
Represents a value of type T
on a LT
lattice.
Fields
- lattice: the
AbstractLattice
object the value is defined on - values: the values on different sites
LatticeModels.LatticeValue
— MethodLatticeValue(lat, values)
Constructs a LatticeValue object.
Arguments
lat
: the lattice the value is defined on.values
: anAbstractVector
of values on the lattice.
Example
julia> using LatticeModels
julia> l = SquareLattice(2, 2);
julia> LatticeValue(l, [1, 2, 3, 4]) # Custom values
LatticeValue{Int64} on a 4-site SquareLattice in 2D space
Values stored in a Vector{Int64}:
[1, 2, 3, 4]
julia> x = LatticeValue(l, :x) # x-coordinate
LatticeValue{Float64} on a 4-site SquareLattice in 2D space
Values stored in a Vector{Float64}:
[1.0, 1.0, 2.0, 2.0]
julia> l2 = SquareLattice(3, 3);
julia> y = LatticeValue(l2, :y) # y-coordinate
LatticeValue{Float64} on a 9-site SquareLattice in 2D space
Values stored in a Vector{Float64}:
[1.0, 2.0, 3.0, 1.0, 2.0, 3.0, 1.0, 2.0, 3.0]
julia> x + y
ERROR: Matching lattices expected.
Got following:
#1: 4-site SquareLattice in 2D space
#2: 9-site SquareLattice in 2D space
[...]
LatticeModels.LatticeValue
— MethodLatticeValue(lat, prop)
Generates a LatticeValue
representing the site property prop
on lattice lat
.
Arguments
lat
: the lattice the value is defined on.prop
: theSiteProperty
to be represented. Can be aSiteProperty
or aSymbol
defining it.
LatticeModels.brastate
— Methodbrastate(lv)
Converts a LatticeValue
to a Bra
wavefunction vector.
LatticeModels.coordvalue
— Methodcoordvalue(lat, coord)
Generates a LatticeValue
representing the spatial coordinate coord
on lattice lat
.
Arguments
lat
: the lattice the value is defined on.coord
: the coordinate to be represented. Can be an integer or a symbol.
LatticeModels.coordvalues
— Methodcoordvalues(lat)
Generates a tuple of LatticeValue
s representing spatial coordinates on lattice lat
.
LatticeModels.ketstate
— Methodketstate(lv)
Converts a LatticeValue
to a Ket
wavefunction vector.
Operators and observables
LatticeModels.LatticeBasis
— TypeLatticeBasis <: QuantumOpticsBase.Basis
Basis for a lattice system.
LatticeModels.coordoperator
— Methodcoordoperator(sys, crd)
coordoperator(basis, crd)
coordoperator(lat[, internal], crd)
Generate a coordinate operator for the given lattice.
Arguments
sys
: aSystem
for which the coordinate operators are to be generated.basis
: a one-particleBasis
for which the coordinate operators are to be generated.lat
: a lattice for which the coordinate operators are to be generated.internal
: The basis for the internal degrees of freedom.crd
: The coordinate to generate the operator for. Must be an integer representing the coordinate index (e. g.1
for x,2
for y, etc.) or a symbol (e. g.:x
,:y
, etc.).
LatticeModels.coordoperators
— Methodcoordoperators(sys)
coordoperators(basis)
coordoperators(lat[, internal])
Generate a Tuple
of coordinate operators for the given lattice.
Arguments
sys
: aSystem
for which the coordinate operators are to be generated.basis
: a one-particleBasis
for which the coordinate operators are to be generated.lat
: a lattice for which the coordinate operators are to be generated.internal
: The basis for the internal degrees of freedom.
LatticeModels.interaction
— Methodinteraction(f, [T, ]sys)
Create an two-site interaction operator for a given NParticles
system. The function f
takes two arguments, which are the two sites, and returns the interaction energy.
LatticeModels.interaction
— Methodinteraction(f, [T, ]sys, K)
Create an 2K
-site interaction operator for a given many-body system. The function f
takes two K
-tuples of integer numbers, which are site indices for creation and annihilation operators, and returns the interaction energy.
If the system sys
has internal degrees of freedom, the function f
should take four K
-tuples: first two are site & internal indices for creation operators, and the last two are the same for annihilation operators.
LatticeModels.AdjacencyMatrix
— MethodAdjacencyMatrix(op::Operator)
Generates an AdjacencyMatrix
for the provided operator.
LatticeModels.localdensity
— Methodlocaldensity(state)
Compute the local density of given state
. The result is a LatticeValue
with the same lattice as the input state.
Arguments
state
: AKet
orBra
representing the wavefunction or anOperator
representing the density matrix.
LatticeModels.localexpect
— Methodlocalexpect(op, state)
Compute the expectation values of the operator op
acting on the local Hilbert space of state
. The result is a LatticeValue
with the same lattice as the input state.
Arguments
op
: ADataOperator
representing the operator.state
: AKet
orBra
representing the wavefunction or anOperator
representing the density matrix.
Hamiltonians
LatticeModels.Hamiltonian
— TypeHamiltonian <: QuantumOpticsBase.DataOperator
A wrapper for a Hamiltonian operator. Contains the operator matrix and the system it acts on.
Hamiltonian(sys, op)
Create a Hamiltonian operator for a given system and a given operator.
Arguments
sys
: the system the Hamiltonian acts on.op
: the operator matrix.
Example
julia> using LatticeModels
julia> l = SquareLattice(4, 4);
julia> H = tightbinding_hamiltonian(l)
Hamiltonian(dim=16x16)
System: One particle on 16-site SquareLattice in 2D space
16×16 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 48 stored entries:
⎡⠪⡢⠑⢄⠀⠀⠀⠀⎤
⎢⠑⢄⠪⡢⠑⢄⠀⠀⎥
⎢⠀⠀⠑⢄⠪⡢⠑⢄⎥
⎣⠀⠀⠀⠀⠑⢄⠪⡢⎦
julia> l2 = SquareLattice(5, 5);
julia> H2 = tightbinding_hamiltonian(l2)
Hamiltonian(dim=25x25)
System: One particle on 25-site SquareLattice in 2D space
25×25 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 80 stored entries:
⎡⠪⡢⡈⠢⡀⠀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠢⡈⠠⡢⡈⠢⡀⠀⠀⠀⠀⠀⠀⎥
⎢⠀⠈⠢⡈⠊⡠⡈⠢⡀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠈⠢⡈⠪⠂⡈⠢⡀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠈⠢⡈⠪⡢⠈⠢⡀⎥
⎢⠀⠀⠀⠀⠀⠀⠀⠈⠢⡀⠪⡢⡀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠀⠈⠀⎦
julia> H + H2
ERROR: Incompatible Hamiltonians:
#1: One particle on 16-site SquareLattice in 2D space
#2: One particle on 25-site SquareLattice in 2D space
[...]
LatticeModels.NParticles
— MethodNParticles(lat[, internal], N[; T=0, statistics=FermiDirac, occupations_type])
NParticles(sys, N[; T=0, statistics=FermiDirac, occupations_type])
Create a manybody system with a given lattice and a given number of particles.
Arguments
lat
: the lattice of the system.internal
: The basis for the internal degrees of freedom.sys
: a one-particle system.N
: the number of particles in the system.
Keyword Arguments
T
: the temperature of the system. Default is0
.statistics
: the statistics of the particles. Default isFermiDirac
.occupations_type
: The occupations type for the many-body operator. By default, the occupation numbers are stored in vectors, but you can use, for example, set it toFermionBitstring
s for better performance on fermion systems.
Example
julia> using LatticeModels
julia> lat = SquareLattice(3, 3);
julia> NParticles(lat, 4, statistics=BoseEinstein)
NParticles(4 bosons) on 9-site SquareLattice in 2D space
LatticeModels.System
— MethodSystem(lat[, internal; T, μ, N, statistics])
Create a system with a given lattice and optionally internal degrees of freedom.
Arguments
lat
: the lattice of the system.internal
: The basis for the internal degrees of freedom.
Keyword Arguments
T
: the temperature of the system. Default is0
.μ
: the chemical potential of the system. Usemu
synonym if Unicode input is not available.N
: the number of particles in the system.statistics
: the statistics of the particles. Default isFermiDirac
.
Example
julia> using LatticeModels
julia> lat = SquareLattice(3, 3);
julia> System(lat)
One particle on 9-site SquareLattice in 2D space
julia> System(lat, N=4, statistics=BoseEinstein)
4 non-interacting bosons on 9-site SquareLattice in 2D space
julia> System(lat, mu=0, statistics=BoseEinstein)
Non-interactng bosons with fixed μ=0.0 on 9-site SquareLattice in 2D space
LatticeModels.System
— MethodSystem(mb[; T])
Create a system with a given many-body basis mb
.
This function is used to create a many-body system from an arbitrary many-body basis with a lattice.
Example
julia> using LatticeModels
julia> lat = SquareLattice(3, 3);
julia> bas = SpinBasis(1//2) ⊗ LatticeBasis(lat);
julia> mbas = ManyBodyBasis(bas, fermionstates(bas, 2));
julia> System(mbas, T=2)
Many-body system on (9-site SquareLattice in 2D space) ⊗ Spin(1/2) (153 states, T=2.0)
LatticeModels.sample
— MethodReturns the Sample
of the object.
Define this function for your type to implement Sample
API.
This function can be considered stable internal API. Feel free to use it in your packages.
LatticeModels.OperatorBuilder
— TypeOperatorBuilder
A helper struct for building custom operators. This struct is used to build operators for a given system or lattice.
LatticeModels.OperatorBuilder
— MethodOperatorBuilder([T, ]sys, [; field, auto_hermitian, auto_pbc_field])
OperatorBuilder([T, ]lat, [internal; field, auto_hermitian])
Construct an OperatorBuilder
for a given system or lattice.
Arguments
T
: The type of the matrix elements. Defaults toComplexF64
.sys
: ASystem
object representing the system.lat
: The lattice on which the operator is defined.internal
: The basis for the internal degrees of freedom.
Keyword arguments
field
: The gauge field to use for the hopping operators. Defaults toNoField()
, which corresponds to zero magnetic field.auto_hermitian
: Whether to automatically add the hermitian conjugate of the operator. Defaults tofalse
.auto_pbc_field
: Whether to automatically adapt the field to the periodic boundary conditions of the lattice. Defaults totrue
.
Example
julia> using LatticeModels
julia> l = SquareLattice(5, 5);
julia> builder = OperatorBuilder(l, field=LandauGauge(0.1), auto_hermitian=true)
OperatorBuilder(field=LandauGauge(0.1), auto_hermitian=true)
System: One particle on 25-site SquareLattice in 2D space
julia> hx = Bravais[1, 0]; hy = Bravais[0, 1];
julia> for site in l; builder[site, site + hx] = builder[site, site + hy] = 1; end
julia> H = Hamiltonian(builder)
Hamiltonian(dim=25x25)
System: One particle on 25-site SquareLattice in 2D space
25×25 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 80 stored entries:
⎡⠪⡢⡈⠢⡀⠀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠢⡈⠠⡢⡈⠢⡀⠀⠀⠀⠀⠀⠀⎥
⎢⠀⠈⠢⡈⠊⡠⡈⠢⡀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠈⠢⡈⠪⠂⡈⠢⡀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠈⠢⡈⠪⡢⠈⠢⡀⎥
⎢⠀⠀⠀⠀⠀⠀⠀⠈⠢⡀⠪⡢⡀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠀⠈⠀⎦
julia> H == tightbinding_hamiltonian(l, field=LandauGauge(0.1))
true
LatticeModels.construct_hamiltonian
— Methodconstruct_hamiltonian([T, ]sys, terms...[; field, auto_pbc_field])
construct_hamiltonian([T, ]lat[, internal, terms...; field, auto_pbc_field])
Construct a Hamiltonian for the given system. Does the same as construct_operator
, but wraps the result in a Hamiltonian
type.
LatticeModels.construct_operator
— Methodconstruct_operator([T, ]sys, terms...[; field, auto_pbc_field])
construct_operator([T, ]lat[, internal, terms...; field])
Construct an operator for the given system.
Each of the terms
describes a term of the Hamiltonian. The term can be given in several ways:
- A
DataOperator
on the lattice, internal or composite basis (will be matched automatically). - A
Pair
of an "internal" and an "on-lattice" part (e.g.int_p => lat_p
):- The "internal" part can be a
DataOperator
, a matrix or a number. - The "on-lattice" part can be a
LatticeValue
(represents a diagonal term), a site (represents a local on-site potential), a bond (represents a hopping term) or asite1 => site2
pair (represents a single hopping). - Identity "internal" or "on-lattice" parts can be omitted.
- The "internal" part can be a
See documentation for more details.
Arguments
T
: The element type of the Hamiltonian. Default isComplexF64
.sys
: TheSystem
for which the Hamiltonian is constructed.lat
: The lattice for which the Hamiltonian is constructed.internal
: The basis for the internal degrees of freedom.
Keyword Arguments
field
: The gauge field to use for the bond operators. Default isNoField()
.auto_pbc_field
: Whether to automatically adapt the field to the periodic boundary conditions of the lattice. Defaults totrue
.
LatticeModels.tightbinding_hamiltonian
— Methodtightbinding_hamiltonian([T, ]sys[, args...; t1=1, t2=0, t3=0, field, auto_pbc_field])
tightbinding_hamiltonian([T, ]lat[, internal, args...; t1=1, t2=0, t3=0, field, auto_pbc_field])
Construct a tight-binding Hamiltonian for the given system.
Arguments
T
: The element type of the Hamiltonian. Default isComplexF64
.sys
: TheSystem
for which the Hamiltonian is constructed.lat
: The lattice for which the Hamiltonian is constructed.internal
: The basis for the internal degrees of freedom.
All other arguments are interpreted as terms of the Hamiltonian and passed to construct_hamiltonian
.
Keyword Arguments
t1
,t2
,t3
: The hopping amplitudes for the nearest, next-nearest, and next-next-nearest neighbors, respectively.field
: The gauge field to use for the bond operators. Default isNoField()
.auto_pbc_field
: Whether to automatically adapt the field to the periodic boundary conditions of the lattice. Defaults totrue
.
LatticeModels.GaugeField
— TypeGaugeField <: AbstractField
A gauge field defined by a vector potential function.
GaugeField(func; n)
Create a gauge field with a given vector potential function func
.
Arguments
func
: a function that takes a point in space and returns the vector potential at this point as aSVector
orTuple
.n
: the number of steps to use in the trapezoidal rule integration.
Example
field = GaugeField(n = 10) do p
(-0.5 * p[2], 0.5 * p[1], 0)
end
LatticeModels.LineIntegralGaugeField
— TypeLineIntegralGaugeField <: AbstractField
A gauge field defined by a line integral function.
LineIntegralGaugeField(func)
Create a gauge field with a given line integral function func
. The function should take two points in space and return the line integral of the vector potential between these points.
Example
field = LineIntegralGaugeField() do p1, p2
0.5 * (p1[1] * p2[2] - p1[2] * p2[1]) # A = [-y/2, x/2, 0]
end
LatticeModels.NoField
— TypeNoField <: AbstractField
A stub object representing zero magnetic field. Use it as a default magnetic field argument in functions — this will not cause any performance overhead.
LatticeModels.line_integral
— Methodline_integral(field, p1, p2[, n_steps=1])
Calculates the $\int_{p1}^{p2} \overrightarrow{A} \cdot \overrightarrow{dl}$ integral using the trapezoidal rule. Increase n_steps
to improve accuracy (note that for linear field gauges like Landau or symmetrical the formula is already pefrectly accurate). If needed, redefine this function for specific field types — this is likely to boost accuracy and performance.
LatticeModels.vector_potential
— Methodvector_potential(field, point)
Returns vector potential $\overrightarrow{A}$ for field
in location point
.
This function should be defined for new field types, but it is not necessary unless you want to use built-in trapezoidal rule integrating.
LatticeModels.LandauGauge
— TypeLandauGauge <: AbstractField
An object representing Landau gauge uniform magnetic field along z-axis.
Fields
B
: The magnetic field value
LatticeModels.PointFlux
— TypePointFlux(flux, [point; gauge])
Construct a PointFlux
object with given flux and point.
The optional gauge
argument can be used to specify the gauge of the field. Possible values are :axial
($A(r) = B \times \frac{r}{|r|}$) and :singular
(the the phase changes if the particle passes below the point). The default is :axial
.
LatticeModels.PointFlux
— TypePointFlux{GaugeT} <: AbstractField
An object representing a small magnetic flux through given point. The field is directed along z-axis.
Fields
flux
: The magnetic flux value.point
: ATuple
of x and y coordinates of the point.
LatticeModels.PointFluxes
— TypePointFluxes{GaugeT} <: AbstractField
An object representing a collection of small magnetic fluxes through given points. The field is directed along z-axis.
Fields
fluxes
: The magnetic flux values.points
: A vector ofNTuple{2,Float64}
points.
LatticeModels.PointFluxes
— MethodPointFluxes([fluxes, points; offset=(0, 0), gauge=:axial])
Construct a PointFluxes
object with given fluxes and points.
The optional gauge
argument can be used to specify the gauge of the field.
Arguments
fluxes
: A vector of flux values. Also can be a single value, in which case it will be broadcasted to all points.points
: A vector of points or aAbstractLattice
object. In the latter case the sites will be interpreted as points.
If both arguments are omitted, an empty field is created. You can add fluxes to it later using push!
or append!
.
Keyword arguments
offset
: An offset to add to the points, default is(0, 0)
. Valid only ifpoints
is aAbstractLattice
, otherwise an error is thrown.gauge
: The gauge of the field. Possible values are:axial
and:singular
. The default is:axial
.
LatticeModels.PointFluxes
— MethodPointFluxes(fields[; gauge])
Construct a PointFluxes
object from a collection of PointFlux
objects.
An error is thrown if the gauges of the fields are inconsistent. You can specify the gauge of the field explicitly using the gauge
keyword argument.
See also: PointFlux
.
LatticeModels.SymmetricGauge
— TypeSymmetricGauge <: AbstractField
An object representing symmetrically gauged uniform magnetic field along z-axis.
Fields
B
: The magnetic field value
Base.append!
— Methodappend!(fields, fields2)
Add a PointFluxes
object to another PointFluxes
object. An error is thrown if the gauges of the fields are inconsistent.
Base.push!
— Methodpush!(fields, field)
Add a PointFlux
object to a PointFluxes
object. An error is thrown if the gauges of the fields are inconsistent.
LatticeModels.periodic_fluxes
— Methodperiodic_fluxes(l, fl)
Construct a PointFluxes
object by periodic replication of a point flux over a Bravais lattice.
Arguments
l
: A Bravais lattice.fl
: APointFlux
object.
Built-in models
LatticeModels.FermiHubbardSpinSystem
— MethodFermiHubbardSpinSystem(l, N_up, N_down[; occupations_type])
Generates a Fermi-Hubbard model basis system on given lattice l
. The number of particles is N_up
and N_down
for spin up and down respectively. The occupations_type
keyword argument defines the type of the occupation basis.
LatticeModels.bosehubbard
— Methodbosehubbard([type, ]lat, N[; U, T, t1, t2, t3, field])
$\hat{H} = \sum_{i,j}^\text{sites} t_{ij} c^\dagger_i c_j + \sum_i^\text{sites} \frac{U}{2} \hat{n}_i (\hat{n}_i - 1)$
Generates a Bose-Hubbard model hamiltonian on given lattice lat
.
Arguments
type
: The element type of the resulting operator. Default isComplexF64
.N
: The number of particles.
Keyword arguments
t1
,t2
,t3
denote the coefficient on first, second and third hoppings respectively. By defaultt1
is equal to one, the rest are zero.U
: The interaction strength. Default is zero.T
: The temperature of the system. Default is zero.field
: The magnetic field. Default isNoField()
.
LatticeModels.fermihubbard
— Methodfermihubbard([type, ]lat, N[; U, T, t1, t2, t3, field])
$\hat{H} = \sum_{i,j}^\text{sites} t_{ij} c^\dagger_i c_j + \sum_i^\text{sites} U \hat{n}_i^{\uparrow} \hat{n}_i^{\downarrow}$
Generates a Fermi-Hubbard model hamiltonian on given lattice lat
.
Arguments
type
: The element type of the resulting operator. Default isComplexF64
.N
: The number of particles.
Keyword arguments
t1
,t2
,t3
denote the coefficient on first, second and third hoppings respectively. By defaultt1
is equal to one, the rest are zero.U
: The interaction strength. Default is zero.T
: The temperature of the system. Default is zero.field
: The magnetic field. Default isNoField()
.
LatticeModels.haldane
— Functionhaldane(lat, t1, t2[, m=0; T, μ, field, statistics])
$\hat{H} = \sum_i^\text{sublattice A} m c^\dagger_i c_i + \sum_j^\text{sublattice B} m c^\dagger_j c_j + \sum_{i, j}^\text{adjacent} \left( t_1 c^\dagger_i c_j + h. c. \right) + \sum_{i, j}^\text{2-connected, counter-clockwise} \left( i \cdot t_2 c^\dagger_i c_j + h. c. \right)$
Generates a Haldane topological insulator hamiltonian operator on given lattice lat
.
Keyword arguments
T
: The temperature of the system. Default is zero.μ
: The chemical potential. Usemu
as a synonym if Unicode input is not available.field
: The magnetic field. Default isNoField()
.statistics
defines the particle statistics, eitherFermiDirac
orBoseEinstein
.
LatticeModels.hubbard
— Methodhubbard([type, ]sys[; U, T, t1, t2, t3, field])
Generates a Hubbard model hamiltonian on given manybody system sys
.
See fermihubbard
, bosehubbard
for more specific models.
LatticeModels.kanemele
— Methodkanemele(lat, t1, t2[; T, μ, field, statistics])
$\hat{H} = \sum_{i, j}^\text{adjacent} \left( t_1 c^\dagger_i c_j + h. c. \right) + \sum_{i, j}^\text{2-connected, counter-clockwise} \left( i \cdot t_2 c^\dagger_i σ_z c_j + h. c. \right)$
Generates a Kane-Mele hamiltonian operator on given lattice lat
.
Keyword arguments
T
: The temperature of the system. Default is zero.μ
: The chemical potential. Use keywordmu
as a synonym if Unicode input is not available.field
: The magnetic field. Default isNoField()
.statistics
defines the particle statistics, eitherFermiDirac
orBoseEinstein
.
LatticeModels.qwz
— Methodqwz(m[; T, μ, field, statistics])
qwz(lat[, m; T, μ, field, statistics])
$\hat{H} = \sum_i^\text{sites} m_i c^\dagger_i \sigma_z c_i + \sum_i^\text{sites} \left( c^\dagger_{i + \hat{x}} \frac{\sigma_z - i \sigma_x}{2} c_i + c^\dagger_{i + \hat{y}} \frac{\sigma_z - i \sigma_y}{2} c_i + h. c. \right)$
Generates a QWZ model hamiltonian operator on given square lattice lat
.
Arguments
m
(either aLatticeValue
or a number) defines the $m_i$ factors
Keyword arguments
T
: The temperature of the system. Default is zero.μ
: The chemical potential. Usemu
as a synonym if Unicode input is not available.field
: The magnetic field. Default isNoField()
.statistics
defines the particle statistics, eitherFermiDirac
orBoseEinstein
.
Diagonalization
LatticeModels.Eigensystem
— TypeEigensystem{LT, MT} where {LT<:AbstractLattice, MT<:AbstractMatrix}
Eigenvalues and eigenvectors for some operator.
LatticeModels.densitymatrix
— Methoddensitymatrix(eig::Eigensystem[; T=0, μ, N, statistics, info=true])
Creates an Operator
representing a equilibrium density matrix, given the eigensystem eig
of the Hamiltonian.
The resulting distribution will be Fermi-Dirac or Bose-Einstein if the statistics
is specified, otherwise the Gibbs distribution will be used.
Keyword arguments
T
is the temperature of the system. Default is zero.μ
is the chemical potential. Usemu
as a synonym if Unicode input is not available.N
is the number of particles. If specified, the chemical potential is found automatically.statistics
defines the particle statistics, eitherFermiDirac
orBoseEinstein
.info
is a boolean flag to enable/disable logging. Default istrue
.
Note that if eig
is a diagonalized Hamiltonian
, the μ
, N
and statistics
parameters are inserted automatically.
LatticeModels.diagonalize
— Methoddiagonalize(op::DataOperator[, routine; params...])
Finds eigenvalues and eigenvectors for a Operator
and stores them in an Eigensystem
.
Two routines are available:
:lapack
uses theeigen
function from the standardLinearAlgebra
package.:krylovkit
uses the Lanczos algorithm from theKrylovKit
package. Accepts following parameters:v0
is the starting vector. Default isrand(ComplexF64, size(op.data, 1))
.n
is the target number of eigenvectors. Default is 10.
KrylovKit.eigsolve
function. See its documentation for details.:auto
automatically selects the routine based on the size of the operator.
The default routine is :lapack
for dense operators. If the operator matrix is less than 5000×5000, it is automatically converted to a dense operator. In other cases :krylovkit
is used.
Example
julia> using LatticeModels
julia> l = SquareLattice(4, 4);
julia> H = tightbinding_hamiltonian(l)
Hamiltonian(dim=16x16)
System: One particle on 16-site SquareLattice in 2D space
16×16 SparseArrays.SparseMatrixCSC{ComplexF64, Int64} with 48 stored entries:
⎡⠪⡢⠑⢄⠀⠀⠀⠀⎤
⎢⠑⢄⠪⡢⠑⢄⠀⠀⎥
⎢⠀⠀⠑⢄⠪⡢⠑⢄⎥
⎣⠀⠀⠀⠀⠑⢄⠪⡢⎦
julia> eig = diagonalize(H)
Diagonalized Hamiltonian (16 eigenvectors)
Eigenvalues in range -3.23607 .. 3.23607
System: One particle on 16-site SquareLattice in 2D space
LatticeModels.findgroundstate
— Methodfindgroundstate(eig::HamiltonianEigensystem)
findgroundstate(ham::Hamiltonian)
Finds the ground state of a Hamiltonian. Returns the energy and the state.
Example
eig = diagonalize(ham)
E, ψ = findgroundstate(eig)
LatticeModels.groundstate
— Methodgroundstate(eig::HamiltonianEigensystem)
groundstate(ham::Hamiltonian)
Finds the ground state of a Hamiltonian. Returns the state.
Example
eig = diagonalize(ham)
ψ = groundstate(eig)
LatticeModels.projector
— Methodprojector(f, eig::Eigensystem)
Returns an Operator
representing a function applied to the diagonalized operator defined by the formula below:
$\hat{\mathcal{P}} = \sum_i f(A_i) |\psi_i⟩⟨\psi_i|$
LatticeModels.projector
— Methodprojector(eig::Eigensystem)
returns an Operator
that projects onto the eigenvectors of the spectrum, defined by the formula below.
$\hat{\mathcal{P}} = \sum_i |\psi_i⟩⟨\psi_i|$
Green's function
LatticeModels.GreenFunction
— TypeGreenFunction
A Green's function for a given lattice and Hamiltonian.
LatticeModels.diagonalelements
— Methoddiagonalelements(gf::GreenFunctionEval)
Return the diagonal elements of the Green's function as a LatticeValue
.
LatticeModels.dos
— Methoddos(eig[, E; broaden])
dos(gf[, E; broaden])
Calculates the DOS (density of states) for a given eigensystem at energy E
. If E
is not specified, a function that calculates the DOS at a given energy is returned.
Arguments
eig
is anEigensystem
orHamiltonianEigensystem
.gf
is aGreenFunction
.E
is the energy at which the DOS is calculated.
Keyword arguments
broaden
is the broadening factor for the energy levels, default is0.1
.
LatticeModels.greenfunction
— Methodgreenfunction(psi0, hamp, hamm[; E₀, tol, showprogress, kw...])
Calculates the Green's function for a many-body system with a given initial state psi0
.
Arguments
psi0
is the initial state.hamp
is the Hamiltonian for the subspace with one more particle than inpsi0
.hamm
is the Hamiltonian for the subspace with one less particle than inpsi0
.
Keyword arguments
E₀
is the energy shift for the Green's function. Default is0
. UseE0
as a synonym if Unicode input is not available.tol
is the tolerance for the new eigenvectors. Default is1e-5
.showprogress
is a flag to show the progress bar. Default istrue
.
All other keyword arguments are passed to the diagonalize
function. See its documentation for details.
LatticeModels.greenfunction
— Methodgreenfunction(ham_eig::HamiltonianEigensystem)
Creates a Green's function for a given one-body Hamiltonian eigensystem.
LatticeModels.ldos
— Methodldos(gf::GreenFunction, E[, site; broaden])
ldos(gf::GreenFunction, site[; broaden])
Calculates the LDOS (local density of states) for a given Green's function at energy E
. broaden
is the broadening factor for the energy levels, default is 0.1
.
If site
is not specified, the LDOS for all sites is calculated and returned as a LatticeValue
. Otherwise, the LDOS for the given site is returned as a Real
value.
If E
is not specified, a function that calculates the LDOS at site site
for given energy is returned.
Currents
LatticeModels.Currents
— TypeCurrents <: AbstractCurrents
A AbstractCurrents
instance that stores values for all currents explicitly.
LatticeModels.Currents
— MethodCurrents(currents[, adjacency_matrix])
Creates a Currents
instance for currents
.
Arguments
currents
: TheAbstractCurrents
object to be turned intoCurrents
. That might be time-consuming, because this requires evaluation of the current between all pairs.adjacency_matrix
: If provided, the current will be evaluated only between adjacent sites.
Examples
julia> using LatticeModels
julia> lat = SquareLattice(4, 4); site1, site2 = lat[1:2];
julia> H0 = tightbinding_hamiltonian(lat); psi = groundstate(H0);
julia> H1 = tightbinding_hamiltonian(lat, field=LandauGauge(0.1));
julia> currents = DensityCurrents(H1, psi)
Density currents for system:
One particle on 16-site SquareLattice in 2D space
julia> c2 = Currents(currents)
Currents{SparseArrays.SparseMatrixCSC{Float64, Int64}} on 16-site SquareLattice in 2D space
julia> c2[site1, site2] ≈ currents[site1, site2]
true
LatticeModels.DensityCurrents
— TypeDensityCurrents <: AbstractCurrents
Density currents for given state and given hamiltonian.
DensityCurrents(hamiltonian, state)
Constructs a DensityCurrents
object for given hamiltonian
and state
.
Arguments
hamiltonian
: AHamiltonian
object representing the Hamiltonian of the system.state
: AKet
orBra
representing the wavefunction or anOperator
representing the density matrix.
Example
Examples
julia> using LatticeModels
julia> lat = SquareLattice(4, 4);
julia> H0 = tightbinding_hamiltonian(lat); psi = groundstate(H0);
julia> H1 = tightbinding_hamiltonian(lat, field=LandauGauge(0.1));
julia> currents = DensityCurrents(H1, psi)
Density currents for system:
One particle on 16-site SquareLattice in 2D space
LatticeModels.LocalOperatorCurrents
— TypeLocalOperatorCurrents <: AbstractCurrents
Local operator (e. g. spin) currents for given state and given hamiltonian.
LocalOperatorCurrents(hamiltonian, state, op)
Constructs a DensityCurrents
object for given hamiltonian
and state
.
Arguments
hamiltonian
: AHamiltonian
object representing the Hamiltonian of the system.state
: AKet
orBra
representing the wavefunction or anOperator
representing the density matrix.op
: A local (on-site) operator; either anOperator
or a matrix of such.
Example
julia> using LatticeModels
julia> lat = SquareLattice(4, 4); site1, site2 = lat[1:2];
julia> H0 = qwz(lat); psi = groundstate(H0);
julia> H1 = qwz(lat, field=LandauGauge(0.1));
julia> op = [1 0; 0 -1]; # Spin operator
julia> currents = LocalOperatorCurrents(H1, psi, op)
Currents of Operator(Spin(1/2))
1 0
0 -1
For system:
One particle on (16-site SquareLattice in 2D space) ⊗ Spin(1/2)
julia> op2 = one(SpinBasis(3//2)); # Invalid operator
julia> LocalOperatorCurrents(H1, psi, op2)
ERROR: ArgumentError: Operator must be defined on the internal basis of the Hamiltonian.
[...]
LatticeModels.SubCurrents
— TypeSubCurrents{CT<:AbstractCurrents} <: AbstractCurrents
A lazy wrapper for a Currents
object representing the same currents but on a smaller lattice.
LatticeModels.currentsfrom
— Methodcurrentsfrom(currents, src)
Create a LatticeValue
object with the currents from src
region to all other sites.
Arguments
currents
: TheAbstractCurrents
object to process.src
: The source region. Can be a site/collection of sites or aLatticeValue{Bool}
mask.
LatticeModels.currentsfromto
— Functioncurrentsfromto(currents, src[, dst])
Finds the total current from src
to dst
regions. If dst
is not provided, the current from src
to all other sites is returned.
Arguments
currents
: TheAbstractCurrents
object to process.src
: The source region.dst
: The destination region.
Both src
and dst
can be a site/collection of sites or a LatticeValue{Bool}
mask.
LatticeModels.lattice
— MethodGets the lattice where the given AbstractCurrents
object is defined.
Evolution
LatticeModels.CachedExp
— TypeCachedExp([ham; threshold=1e-10, nztol=1e-14])
A EvolutionSolver
that finds the matrix exponential of the Hamiltonian and caches it. The matrix exponential is computed using a scaling and squaring method, so this solver works well with sparse or GPU arrays.
Arguments
ham
: The Hamiltonian of the system. It can be anOperator
or its matrix.threshold
: The threshold for the error in the matrix exponential.nztol
: The tolerance for dropping small elements in the matrix exponential if it is sparse.
LatticeModels.Evolution
— TypeEvolution([solver, ]hamiltonian, states...; timedomain, namedstates...)
Create an Evolution
object that can be used to evolve states in time according to the Schrödinger equation.
Arguments
solver
: AEvolutionSolver
object that will be used to evolve the states. If omitted, aCachedExp
solver will be created.hamiltonian
: The Hamiltonian of the system. It can be a matrix, a time-dependent operator or a function that returns the Hamiltonian at a given time.states
andnamedstates
: The states to be evolved. They can beKet
wavefunctions orDataOperator
density matrices.timedomain
: The time domain to be used for the evolution. If omitted, the non-iterableEvolution
object will be returned, and you will be able to call it with the time domain later.
See EvolutionSolver
for more information about solvers.
Please note that the Evolution
object is a stateful iterator. This means that it keeps track of the current time and the states as they evolve. You can perform evolution multiple times, but the timeline will be kept and the states will be updated in place.
Also do not edit the states in place, as this will affect the evolution. If you need to modify the states or save them, make a copy of them first.
LatticeModels.KrylovKitExp
— TypeKrylovKitExp([ham; kw...])
A EvolutionSolver
that uses the exponentiate
function from KrylovKit.jl to evolve the wavefunction vectors. This solver is useful for large, sparse, time-dependent Hamiltonians.
Arguments
ham
: The Hamiltonian of the system. It can be anOperator
or its matrix.kw...
: Keyword arguments to be passed toexponentiate
.
LatticeModels.TimeSequence
— TypeTimeSequence{ET}
A time-ordered sequence of values.
LatticeModels.TimeSequence
— MethodTimeSequence(times, values)
Constructs a TimeSequence
with the given times
and values
.
LatticeModels.TimeSequence
— MethodTimeSequence(f, ev, times)
TimeSequence(f, ev_iter)
Constructs a TimeSequence by iterating the evolution iterator ev
and applying the function f
to each moment.
Arguments
f
: A function that takes a moment and returns a value. The function is applied to each moment in the evolution.ev
: AnEvolution
object.times
: A range of times to evaluate the function at.ev_iter
: An evoltuion iterator that yields moments. Think of it asev(times)
.
LatticeModels.TimeSequence
— MethodTimeSequence{ET}()
Constructs an empty TimeSequence
with eltype ET
.
LatticeModels.differentiate!
— Methoddifferentiate!(tseq::TimeSequence)
Differentiate the values stored in tseq
by time using the symmetric difference formula. The new values are written into tseq
.
Example
julia> using LatticeModels
julia> tseq = TimeSequence(0:0.1:10, 0:0.1:10) # f(t) = t
TimeSequence{Float64} with 101 entry
Timestamps in range 0.0 .. 10.0:
0.0 => 0.0
0.1 => 0.1
0.2 => 0.2
0.3 => 0.3
0.4 => 0.4
0.5 => 0.5
0.6 => 0.6
0.7 => 0.7
0.8 => 0.8
⋮
julia> differentiate!(tseq) # f'(t) = 1
TimeSequence{Float64} with 100 entries
Timestamps in range 0.05 .. 9.95:
0.05 => 1.0
0.15 => 1.0
0.25 => 1.0
0.35 => 1.0
0.45 => 1.0
0.55 => 1.0
0.65 => 1.0
0.75 => 1.0
0.85 => 1.0
⋮
LatticeModels.differentiate
— Methoddifferentiate(tseq::TimeSequence)
Differentiate the values stored in tseq
and create a copy; see differentiate!
.
LatticeModels.integrate!
— Methodintegrate!(tseq::TimeSequence)
Integrate the values stored in tseq
over time using the trapezoidal rule. The new values are written into tseq
. The first value is set to zero.
Example
julia> using LatticeModels
julia> tseq = TimeSequence(0:0.1:10, 0:0.1:10) # f(t) = t
TimeSequence{Float64} with 101 entry
Timestamps in range 0.0 .. 10.0:
0.0 => 0.0
0.1 => 0.1
0.2 => 0.2
0.3 => 0.3
0.4 => 0.4
0.5 => 0.5
0.6 => 0.6
0.7 => 0.7
0.8 => 0.8
⋮
julia> integrate!(tseq) # F(t) = t^2 / 2
TimeSequence{Float64} with 101 entry
Timestamps in range 0.0 .. 10.0:
0.0 => 0.0
0.1 => 0.005
0.2 => 0.02
0.3 => 0.045
0.4 => 0.08
0.5 => 0.125
0.6 => 0.18
0.7 => 0.245
0.8 => 0.32
⋮
LatticeModels.integrate
— Methodintegrate(tseq::TimeSequence)
Integrate the values stored in tseq
and create a copy; see integrate!
.
LatticeModels.timerange
— Methodtimerange(tseq::TimeSequence)
Returns the range of the timestamps of the TimeSequence
.
LatticeModels.timestamps
— Methodtimestamps(tseq::TimeSequence)
Returns the timestamps of the TimeSequence
.
Internals
LatticeModels.AbstractLattice
— TypeAbstractLattice{SiteT}
An abstract type for a lattice of SiteT
sites.
Methods for subtypes to implement
length(l::AbstractLattice)
: Return the number of sites in the lattice.site_index(l::AbstractLattice, site::SiteT)
: Return the index of the site in the lattice.getindex(l::AbstractLattice, i::Int)
: Return the site with the given index.getindex(l::AbstractLattice, is::AbstractVector{Int})
: Return anAbstractLattice
with the sites at the given indices.
Optional methods for mutable lattices
emptymutable(l::AbstractLattice, ::Type{SiteT})
: Return an empty mutable instance of lattice.copymutable(l::AbstractLattice)
: Return a mutable copy of the lattice.push!(l::AbstractLattice, site::SiteT)
: Add a site to the lattice.deleteat!(l::AbstractLattice, is::AbstractVector{Int})
: Remove the sites with the given indices from the lattice.
LatticeModels.AbstractSite
— TypeAbstractSite{N}
An abstract type for a site of a N
-dimensional lattice.
Fields
coords
: ASVector
of sizeN
representing the spatial coordinates of the site. All subtypes are expected to have this field.
LatticeModels.AbstractBonds
— TypeAbstractBonds{LT}
An abstract type for bonds on some lattice.
Methods for subtypes to implement
lattice(bonds::AbstractBonds)
: Returns the lattice where the bonds are defined.isadjacent(bonds::AbstractBonds, site1::AbstractSite, site2::AbstractSite)
: Returns if the sites are connected by the bonds.
Optional methods for subtypes to implement
adapt_bonds(bonds::AbstractBonds, l::AbstractLattice)
LatticeModels.DirectedBonds
— TypeDirectedBonds{LT} <: AbstractBonds{LT}
An abstract type for bonds on some lattice that have a direction.
Methods for subtypes to implement
lattice(bonds::DirectionalBonds)
: Returns the lattice where the bonds are defined.destinations(bonds::DirectionalBonds, site::AbstractSite)
: Returns the sites where the
site
is connected to, accounting for the direction of the bonds.
LatticeModels.AbstractTranslation
— TypeAbstractTranslation{LT}
An abstract type for translations on some lattice.
Methods for subtypes to implement
lattice(bonds::AbstractTranslation)
: Returns the lattice where the translations are defined.destination(bonds::AbstractTranslation, site::AbstractSite)
: Returns the site where thesite
is translated to.
Optional methods for subtypes to implement
adapt_bonds(bonds::AbstractTranslation, l::AbstractLattice)
: Adapt the translation to the latticel
. The output can be a different type of translation, more fitting for the concrete type of lattice.inv(bonds::AbstractTranslation)
: Returns the inverse of the translation, if any.
LatticeModels.AbstractCurrents
— TypeAbstractCurrents
Supertype for all type representing currents-like values on a lattice. Subtypes must implement Base.getindex(::Int, ::Int)
and lattice
functions.
LatticeModels.LookupTable
— TypeLookupTable
A helper data structure to quickly find the index of a site in a lattice.
Relies on sitekey(site)
and secondarykey(function)
functions to determine the index of a site in the lattice.
Works well under following assumptions:
- The
sitekey
is some integer property of the sites. - The sites in the lattice are ordered by
sitekey
. - The numbering is mostly contiguous, i.e. there are no (or few) gaps in the numbering.
- The
secondarykey
is also integer, mostly contiguous, ordered and unique for all sites with the samesitekey
.
Set them to nothing
to disable usage (this is the default behavior).
LatticeModels.EvolutionSolver
— TypeEvolutionSolver
Abstract type for solvers that can be used to evolve states in time according to the Schrödinger equation.
See also concrete implementations: CachedExp
, KrylovKitExp
.
Methods to implement
update_solver!(solver, hamiltonian, dt, force=false)
: Update the solver to evolve the states according to the given Hamiltonian and time step. Ifforce
istrue
, the solver should always update, even if the Hamiltonian and time step are seemingly the same as the previous ones.step!(solver, state, cache)
: Evolve the given state in time using the solver. Thecache
argument is used to store intermediate results and can benothing
if the solver does not need it.evolution_cache(solver, state)
: Return a cache object that can be used to store intermediate results for the given state. Returnsnothing
if the solver does not need a cache for the given state (this is the default implementation).
LatticeModels.addlookuptable
— Functionaddlookuptable(lat)
Adds a lookup table to the lattice lat
and returns the lattice with the lookup table.
Make sure you add the lookup table to the lattice after you stop making changes to it. Otherwise the results may be unpredictable.
This operation is not in-place.