Momentum space API
See Momentum space guide for the narrative version and Reciprocal lattice and Brillouin zone for the mathematics.
Abstract types
LatticeCore.AbstractMomentumLattice — Type
AbstractMomentumLattice{D, T} <: AbstractLattice{D, T}Abstract supertype for k-space lattices. A momentum lattice is itself an AbstractLattice — its "sites" are k-points, its "positions" are k-vectors in the reciprocal basis — so 02's lattice interface (num_sites, position, traits, test suite) can be re-used for k-space work. See dev/note/04_architecture/05_momentum_space.
Required interface for concrete subtypes
num_k_points(ml)::Intk_point(ml, i::Int)::SVector{D, T}reciprocal_basis(ml)::SMatrix{D, D, T}
The graph-neighbour side of AbstractLattice (neighbors, bonds) is deliberately vacuous for momentum lattices: k-space is not a graph in the same sense. Concrete momentum lattices return an empty neighbour list; MC code should never walk them as a graph.
Concrete momentum lattices
LatticeCore.PeriodicMomentumLattice — Type
PeriodicMomentumLattice{D, T}Concrete momentum lattice for a Bravais periodic lattice. Stores the reciprocal basis, the mesh dimensions, the Γ-offset, and an explicit list of k-points (eager construction).
Construct via monkhorst_pack or gamma_centered rather than calling the struct directly.
Interface
LatticeCore.num_k_points — Function
num_k_points(ml::AbstractMomentumLattice) → IntNumber of k-points stored in the lattice.
LatticeCore.k_point — Function
k_point(ml::AbstractMomentumLattice, i::Int) → SVector{D, T}The i-th k-vector in real (Cartesian) units — i.e. already multiplied through the reciprocal basis.
LatticeCore.reciprocal_basis — Function
reciprocal_basis(ml::AbstractMomentumLattice) → SMatrix{D, D, T}Reciprocal-space basis matrix. Columns are the primitive reciprocal lattice vectors.
Mesh constructors
LatticeCore.monkhorst_pack — Function
monkhorst_pack(basis::SMatrix{D, D, T}, mesh::NTuple{D, Int})
→ PeriodicMomentumLattice{D, T}Construct a Monkhorst–Pack half-shifted mesh over the Brillouin zone spanned by basis. For each axis the fractional k-coordinates are
frac_i ∈ { (n + 0.5) / N - 0.5 | n = 0, …, N - 1 }so the mesh is centred at Γ and avoids the BZ edges.
LatticeCore.gamma_centered — Function
gamma_centered(basis::SMatrix{D, D, T}, mesh::NTuple{D, Int})
→ PeriodicMomentumLattice{D, T}Construct a Γ-centred regular mesh. Fractional k-coordinates are n / N with n = 0, …, N - 1, so the mesh includes Γ and walks to but not through the BZ edge.
Entry points
LatticeCore.reciprocal_lattice — Function
reciprocal_lattice(lat::AbstractLattice) → PeriodicMomentumLatticeConstruct the reciprocal lattice for a Bravais-like lat. Concrete lattices are expected to specialise this method; the fallback throws.
LatticeCore.fourier_module — Function
fourier_module(lat::AbstractLattice) → AbstractMomentumLatticeQuasicrystal-side entry point: construct the discrete Fourier module (Bragg peak set) for a cut-and-project lattice. Concrete quasicrystal lattices implement this in their own package; the fallback throws.
LatticeCore.momentum_lattice — Function
momentum_lattice(lat::AbstractLattice) → AbstractMomentumLatticeUnified trait-dispatched entry point. Returns the reciprocal lattice for Bravais-like structures, the Fourier module for quasicrystals, and throws for lattices without any k-space representation.
Structure factor
LatticeCore.structure_factor — Function
structure_factor(lat, state, k::SVector) → Float64Scalar structure factor at a single k-point:
S(k) = (1/N) |⟨ Σ_i s_i exp(-i k · r_i) ⟩|²computed on a single snapshot (no thermal averaging). The default implementation is naive O(N) per k-point; FFT-based specialisations on regular meshes can be added later without changing the API.
structure_factor(lat, state, ml::AbstractMomentumLattice) → Vector{Float64}Evaluate the structure factor at every k-point of ml. Naive O(N · M) where M = num_k_points(ml).
Quasicrystal Fourier-module skeletons
LatticeCore.AcceptanceWindow — Type
AcceptanceWindowAbstract supertype for acceptance windows in the internal (perp) space of a cut-and-project quasicrystal. Concrete windows (PolygonalWindow, IntervalWindow, etc.) live in QuasiCrystal.jl.
LatticeCore.HyperReciprocalLattice — Type
HyperReciprocalLattice{DPhys, DHyper, T}Infinite-abstract representation of a quasicrystal's reciprocal structure: the higher-dimensional reciprocal basis, the parallel and perpendicular projections, and the acceptance window used to decide peak intensities.
Concrete construction (building the projections, sampling peaks) lives in QuasiCrystal.jl.
LatticeCore.BraggPeakSet — Type
BraggPeakSet{DPhys, T}Finite materialisation of a HyperReciprocalLattice up to a cutoff radius: a list of physical-space k-positions with intensities and a back-pointer to the higher-dimensional indices that generated them.
Being a subtype of AbstractMomentumLattice, a BraggPeakSet can be passed straight to structure_factor or to a StructureFactorObserver so the same code paths work for periodic and quasiperiodic lattices.