API
LatticeUtilities.Bond
LatticeUtilities.Bond
LatticeUtilities.Lattice
LatticeUtilities.Lattice
LatticeUtilities.UnitCell
LatticeUtilities.UnitCell
LatticeUtilities.bond_to_vec
LatticeUtilities.bond_to_vec!
LatticeUtilities.build_neighbor_table
LatticeUtilities.build_neighbor_table
LatticeUtilities.calc_k_point
LatticeUtilities.calc_k_point!
LatticeUtilities.calc_k_points
LatticeUtilities.calc_k_points!
LatticeUtilities.displacement_to_vec
LatticeUtilities.displacement_to_vec!
LatticeUtilities.loc_to_pos
LatticeUtilities.loc_to_pos
LatticeUtilities.loc_to_pos!
LatticeUtilities.loc_to_pos!
LatticeUtilities.loc_to_site
LatticeUtilities.loc_to_site
LatticeUtilities.loc_to_unitcell
LatticeUtilities.map_neighbor_table
LatticeUtilities.norbits
LatticeUtilities.nsites
LatticeUtilities.pbc!
LatticeUtilities.simplify
LatticeUtilities.simplify!
LatticeUtilities.site_to_loc
LatticeUtilities.site_to_loc!
LatticeUtilities.site_to_orbital
LatticeUtilities.site_to_site
LatticeUtilities.site_to_unitcell
LatticeUtilities.sites_to_bond
LatticeUtilities.sites_to_displacement
LatticeUtilities.sites_to_displacement!
LatticeUtilities.translational_avg!
LatticeUtilities.unitcell_to_loc
LatticeUtilities.unitcell_to_loc!
LatticeUtilities.valid_loc
LatticeUtilities.valid_site
LatticeUtilities.Bond
— TypeBond{D}
Defines a bond in a D
dimensional lattice.
Fields
orbitals::NTuple{2,Int}
: Orbital pair associated with bond. Bond goes fromorbitals[1]
toorbitals[2]
.displacement::SVector{D,Int}
: Displacement in unit cells in the direction of each lattice vector associated with bond.
LatticeUtilities.Bond
— MethodBond(orbitals, displacement)
Bond(; orbitals, displacement)
Construct a Bond
LatticeUtilities.Lattice
— TypeLattice{D}
A type defining a finite lattice in D
dimensions.
Fields
N::Int
: Number of unit cells in finite lattice.L::SVector{D,Int}
: Size of finite lattice in the direction of each lattice vector in unit cells.periodic::SVector{D,Bool}
: Specifies whether each lattice vector direction hosts periodic or open boundary conditions.lvec::MVector{D,Int}
: (private) Temporary storage vector containing intermediate location and displacement vectors.
LatticeUtilities.Lattice
— MethodLattice(L, periodic)
Lattice(; L, periodic)
Constructs a Lattice
.
LatticeUtilities.UnitCell
— TypeUnitCell{D,T<:AbstractFloat,N}
A type defining a unit cell in D
spatial dimensions.
Fields
n::Int
: Number of orbitals in the unit cell.lattice_vecs::SMatrix{D,D,T,N}
: Matrix where columns give the lattice vectors.reciprocal_vecs::SMatrix{D,D,T,N}
: Matrix where columns give the reciprocal lattice vectors.basis_vecs::Vector{SVector{D,T}}
: Vector of basis vectors specifying the location of each orbital in unit cell.
LatticeUtilities.UnitCell
— MethodUnitCell(lattice_vecs::AbstractMatrix{T}, basis_vecs::AbstractMatrix{T}) where {T<:AbstractFloat}
UnitCell(lattice_vecs::AbstractVector{Vector{T}}, basis_vecs::AbstractVector{Vector{T}}) where {T<:AbstractFloat}
UnitCell(; lattice_vecs, basis_vecs)
Returns an instance of the type UnitCell
.
Base.:==
— MethodBase.:(==)(b1::Bond{D1}, b2::Bond{D2}) where {D1, D2}
Tests if two bonds b1
and b2
are equivalent.
Base.ndims
— MethodBase.ndims(unit_cell::UnitCell{D}) where {D}
Base.ndims(lattice::Lattice{D}) where {D}
Base.ndims(bond::Bond{D}) where {D}
Return the spatial dimensions of a UnitCell
, Lattice
or Bond
.
Base.show
— MethodBase.show(io::IO, bond::Bond{D}) where {D}
Base.show(io::IO, ::MIME"text/plain", bond::Bond{D}) where {D}
Show lattice.
Base.show
— MethodBase.show(io::IO, lattice::Lattice{D}) where {D}
Base.show(io::IO, ::MIME"text/plain", lattice::Lattice{D}) where {D}
Show lattice.
Base.show
— MethodBase.show(io::IO, uc::UnitCell{D,T}) where {D,T}
Base.show(io::IO, ::MIME"text/plain", uc::UnitCell{D,T}) where {D,T}
Show unit cell.
LatticeUtilities.bond_to_vec!
— Methodbond_to_vec!(Δr::AbstractVector{T}, bond::Bond{D}, unit_cell::UnitCell{D,T}) where {D,T}
Calculate the displacement vector associated with a bond
.
LatticeUtilities.bond_to_vec
— Methodbond_to_vec(bond::Bond{D}, unit_cell::UnitCell{D,T}) where {D,T}
Return the displacement vector associated with a bond
as SVector{D,T}
.
LatticeUtilities.build_neighbor_table
— Methodbuild_neighbor_table(bonds, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
Construct the neighbor table corresponding bonds
, a tuple/vector of Bond
defintions.
LatticeUtilities.build_neighbor_table
— Methodbuild_neighbor_table(bond::Bond{D}, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
Construct the neighbor table corresponding to bond
.
LatticeUtilities.calc_k_point!
— Methodcalc_k_point!(k_point::AbstractVector{T}, k_loc,
unit_cell::UnitCell{D,T}, lattice::Lattice{D}) where {D,T}
Calculate the k-point k_point
corresponding to the k-point location k_loc
where k_loc
is a vector of tuple of integers of length D
.
LatticeUtilities.calc_k_point
— Methodcalc_k_point(k_loc, unit_cell::UnitCell{D,T}, lattice::Lattice{D}) where {D,T}
Return the k-point corresponding to the k-point location k_loc
.
LatticeUtilities.calc_k_points!
— Methodcalc_k_points!(k_points::AbstractArray{T}, unit_cell::UnitCell{D,T}, lattice::Lattice{D}) where {D,T}
Calculate the k-point grid k_points
assicated with a finite lattice.
LatticeUtilities.calc_k_points
— Methodcalc_k_points(unit_cell::UnitCell{D,T}, lattice::Lattice{D}) where {D,T}
Return the k-points associated with a finite lattice as a vector of static vectors, Vector{SVector{D,T}}
. If the system has open boundary conditions in a given direction, it will treat the linear extent of the system in that direction as equalling L=1
for the purposes of calculating the k-points.
LatticeUtilities.displacement_to_vec!
— Methoddisplacement_to_vec!(Δr::AbstractVector{T}, Δl, o_init::Int, o_final::Int, unit_cell::UnitCell{D,T}) where {D,T}
Computes the position space displacement vector Δr
corresponding to a displacement definition given by initial and final orbitals o₁
and o₂
in the unit cell respectively, along with a displacement in unit cells Δl
.
Arguments
Δr::AbstractVector{T}
: displacement vector in position space.Δl
: displacement in unit cells.o_init::Int
: initial orbital in unit cell.o_final::Int
: final orbital in unit cell.unit_cell::UnitCell{D,T}
: unit cell.
LatticeUtilities.displacement_to_vec
— Methoddisplacement_to_vec(Δl, o_init::Int, o_final::Int, unit_cell::UnitCell{D,T})::SVector{D,T} where {D,T}
Returns the position space displacement vector Δr
corresponding to a displacement definition given by initial and final orbitals o_init
and o_final
in the unit cell respectively, along with a displacement in unit cells Δl
.
LatticeUtilities.loc_to_pos!
— Methodloc_to_pos!(r::AbstractVector{T}, l, o::Int, unit_cell::UnitCell{D,T}) where {D,T}
Calculate the position r
of an orbital o
at location l
.
LatticeUtilities.loc_to_pos!
— Methodloc_to_pos!(r::AbstractVector{T}, l, unit_cell::UnitCell{D,T}) where {D,T}
Calculate the position r
of a unit cell at location l
.
LatticeUtilities.loc_to_pos
— Methodloc_to_pos(l, s::Int, unit_cell::UnitCell{D,T}) where {D,T}
Return the position r
of a orbital o
at location l
as a vector or type SVector{D,T}
.
LatticeUtilities.loc_to_pos
— Methodloc_to_pos(l, unit_cell::UnitCell{T}) where {T}
Return the position r
of a unit cell at location l
.
LatticeUtilities.loc_to_site
— Methodloc_to_site(l, o::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
Given a unit cell location l
and orbital species o
, return the corresponding site s
in the lattice. If the location is not valid owing to open boundary conditions then return s = 0
.
LatticeUtilities.loc_to_site
— Methodloc_to_site(u::Int, o::Int, unit_cell::UnitCell{D}) where {D}
Given a unit cell index u
and orbital o
, return the correspond site s
.
LatticeUtilities.loc_to_unitcell
— Methodloc_to_unitcell(l, lattice::Lattice{D}) where {D}
Return the unit cell found at location l
in the lattice.
LatticeUtilities.map_neighbor_table
— Methodmap_neighbor_table(neighbor_table::Matrix{Int})
For a given neighbor table, return a dictionary that reports the bonds and neighbors associated with each site in the lattice. If neighbor_table
is modified, then a new map must be constructed.
LatticeUtilities.norbits
— Methodnorbits(unit_cell::UnitCell)
Return the number of orbitals in a unit cell.
LatticeUtilities.nsites
— Methodnsites(unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
nsites(; unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
Return the total number of sites/orbitals in a finite lattice.
LatticeUtilities.pbc!
— Methodpbc!(l::AbstractVector{Int}, lattice::Lattice{D}) where {D}
Apply periodic boundary to unit cell location l
.
LatticeUtilities.simplify!
— Methodsimplify!(Δl::AbstractVector{Int}, lattice::Lattice{D}) where {D}
Simplify displacement Δl
so that it is as short as possible accounting for periodic boundary conditions where necessary.
LatticeUtilities.simplify
— Methodsimplify(bond::Bond{D}, lattice::Lattice{D}) where {D}
Simplify a bond so that the displacement is the shortest possible accounting for periodic boundary conditions where necessary, returning the new bond.
LatticeUtilities.site_to_loc!
— Methodsite_to_loc!(l::AbstractVector{Int}, s::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
site_to_loc!(; l::AbstractVector{Int}, s::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
For a given site s
in the lattice, calculate the location l
of the unit cell it is in and return the orbital species o
of the site.
LatticeUtilities.site_to_loc
— Methodsite_to_loc(s::Int, unit_cell::UnitCell{D,T}, lattice::Lattice{D}) where {D,T}
For a given site s
in the lattice, return the location l
of the unit cell it is in and the orbital species o
.
LatticeUtilities.site_to_orbital
— Methodsite_to_orbital(s::Int, unit_cell::UnitCell)
Return the orbtial species of site s
.
LatticeUtilities.site_to_site
— Methodsite_to_site(s::Int, Δl, o::Int, unit_cell::UnitCell, lattice::Lattice)
Given an initial site s
, and a displacement in unit cells Δl
and a terminating orbital species o
, return the resulting site s′
in the lattice. If the displacement is not allowed as a result of open boundary conditions, then s′=0
is returned.
LatticeUtilities.site_to_unitcell
— Methodsite_to_unitcell(s::Int, unit_cell::UnitCell)
Return the unit cell u
containing lattice site s
.
LatticeUtilities.sites_to_bond
— Methodsites_to_bond(s::Int, s′::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
Return the Bond
associated with getting displaced from site s
to s′
.
LatticeUtilities.sites_to_displacement!
— Methodfunction sites_to_displacement!(Δl::AbstractVector{Int}, s::Int, s′::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
When getting displaced from site s
to s′
, calculate the corresponding displacement in unit cells Δl
.
LatticeUtilities.sites_to_displacement
— Methodfunction sites_to_displacement(s::Int, s′::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
When getting displaced from site s
to s′
, return the corresponding displacement in unit cells Δl::MVector{D,Int}
.
LatticeUtilities.translational_avg!
— Methodfunction translational_avg!(fg::AbstractArray{Complex{T}}, f::AbstractArray{Complex{T}}, g::AbstractArray{Complex{T}};
restore::Bool=true) where {T<:AbstractFloat}
Let f[i]
and g[j]
be two distinct multi-dimensional arrays, where i
and j
represent an index into them and also correspond to the position of a unit cell in a periodic finite lattice This method then computes in-place the product (f⋅g)[i-j]
that is averaged over translation symmetry. If restore = true
then f
and g
are left unchanged, otherwise they will be left modified.
LatticeUtilities.unitcell_to_loc!
— Methodunitcell_to_loc!(l::AbstractVector{Int}, u::Int, lattice::Lattice{D}) where {D}
Calculate the location l
of a unit cell u
.
LatticeUtilities.unitcell_to_loc
— Methodunitcell_to_loc(u::Int, lattice::Lattice{D})::SVector{D,Int} where {D}
Return the location of unit cell u
as an instance of type SVector{D,Int}
.
LatticeUtilities.valid_loc
— Methodvalid_loc(l, lattice::Lattice{D})::Bool where {D}
Determine if l
describes a valid location in the lattice.
LatticeUtilities.valid_site
— Methodvalid_site(s::Int, unit_cell::UnitCell{D}, lattice::Lattice{D}) where {D}
valid_site(; s, unit_cell, lattice) = valid_site(s, unit_cell, lattice)
Return whether s
is a valid site index.