Solver API

Wave Types

Abstract Types

2D Wave Types

PhoXonic.TEWaveType
TEWave <: PhotonicWave

Transverse Electric polarization (H_z mode) for 2D photonic crystals.

Electric field is in the xy-plane, magnetic field H is along z. The eigenvalue equation is: ∇×(ε⁻¹∇×Hz) = (ω/c)² μ Hz

Field Components

  • Solved: H_z (scalar)
  • Derived: Ex, Ey from H_z

Example

solver = Solver(TEWave(), geo, (64, 64))
solver = Solver(TEWave(), geo, (64, 64), KrylovKitMethod())
solver = Solver(TEWave(), geo, (64, 64), LOBPCGMethod())

See also: TMWave, FullVectorEM

source
PhoXonic.TMWaveType
TMWave <: PhotonicWave

Transverse Magnetic polarization (E_z mode) for 2D photonic crystals.

Magnetic field is in the xy-plane, electric field E is along z. The eigenvalue equation is: ∇×(μ⁻¹∇×Ez) = (ω/c)² ε Ez

Field Components

  • Solved: E_z (scalar)
  • Derived: Hx, Hy from E_z

Example

solver = Solver(TMWave(), geo, (64, 64))
solver = Solver(TMWave(), geo, (64, 64), KrylovKitMethod())
solver = Solver(TMWave(), geo, (64, 64), LOBPCGMethod())

See also: TEWave, FullVectorEM

source
PhoXonic.SHWaveType
SHWave <: PhononicWave

Shear Horizontal (out-of-plane, anti-plane) elastic wave for 2D phononic crystals.

Displacement uz is perpendicular to the xy-plane. The eigenvalue equation is: ∇·(C₄₄∇uz) = -ρω² u_z

Field Components

  • Solved: u_z (scalar)
  • Decoupled from in-plane (P-SV) modes

Notes

  • Eigenvalues ω² can be O(10¹⁰) for typical materials
  • KrylovKitMethod uses automatic scaling for numerical stability
  • LOBPCGMethod works without explicit scaling

Example

steel = IsotropicElastic(ρ=7800.0, λ=115e9, μ=82e9)
epoxy = IsotropicElastic(ρ=1180.0, λ=4.43e9, μ=1.59e9)
geo = Geometry(lat, epoxy, [(Circle([0,0], 0.3), steel)])

solver = Solver(SHWave(), geo, (64, 64))
solver = Solver(SHWave(), geo, (64, 64), KrylovKitMethod())
solver = Solver(SHWave(), geo, (64, 64), LOBPCGMethod())

See also: PSVWave, FullElastic

source
PhoXonic.PSVWaveType
PSVWave <: PhononicWave

P-SV (in-plane) elastic wave for 2D phononic crystals.

Displacement (ux, uy) is in the xy-plane, coupling P (longitudinal) and SV (shear vertical) polarizations.

Field Components

  • Solved: ux, uy (2 components)
  • Produces 2 bands per k-point (quasi-P and quasi-SV)

Notes

  • Eigenvalues ω² can be O(10¹⁰) for typical materials
  • KrylovKitMethod uses automatic scaling for numerical stability
  • LOBPCGMethod works without explicit scaling

Example

steel = IsotropicElastic(ρ=7800.0, λ=115e9, μ=82e9)
epoxy = IsotropicElastic(ρ=1180.0, λ=4.43e9, μ=1.59e9)
geo = Geometry(lat, epoxy, [(Circle([0,0], 0.3), steel)])

solver = Solver(PSVWave(), geo, (64, 64))
solver = Solver(PSVWave(), geo, (64, 64), KrylovKitMethod())
solver = Solver(PSVWave(), geo, (64, 64), LOBPCGMethod())

See also: SHWave, FullElastic

source

1D Wave Types

PhoXonic.Photonic1DType
Photonic1D <: PhotonicWave

Scalar electromagnetic wave for 1D photonic structures (Bragg reflectors, etc.).

The eigenvalue equation is: -d/dx(ε⁻¹ d/dx E) = (ω/c)² E

Example

lat = lattice_1d(1.0)
geo = Geometry(lat, mat1, [(Segment(0.0, 0.5), mat2)])
solver = Solver(Photonic1D(), geo, 128; cutoff=20)
solver = Solver(Photonic1D(), geo, 128, KrylovKitMethod(); cutoff=20)
solver = Solver(Photonic1D(), geo, 128, LOBPCGMethod(); cutoff=20)

See also: TEWave, TMWave, Longitudinal1D

source
PhoXonic.Longitudinal1DType
Longitudinal1D <: PhononicWave

Longitudinal elastic wave for 1D phononic structures (superlattices, etc.).

The eigenvalue equation is: d/dx(C₁₁ du/dx) = -ρω² u

Notes

  • Eigenvalues ω² can be O(10¹⁰) for typical materials
  • KrylovKitMethod uses automatic scaling for numerical stability
  • LOBPCGMethod works without explicit scaling

Example

lat = lattice_1d(1.0)
steel = IsotropicElastic(ρ=7800.0, λ=115e9, μ=82e9)
epoxy = IsotropicElastic(ρ=1180.0, λ=4.43e9, μ=1.59e9)
geo = Geometry(lat, epoxy, [(Segment(0.0, 0.5), steel)])

solver = Solver(Longitudinal1D(), geo, 128; cutoff=20)
solver = Solver(Longitudinal1D(), geo, 128, KrylovKitMethod(); cutoff=20)
solver = Solver(Longitudinal1D(), geo, 128, LOBPCGMethod(); cutoff=20)

See also: SHWave, PSVWave, Photonic1D

source

3D Wave Types

PhoXonic.TransverseEMType
TransverseEM <: PhotonicWave

Transverse electromagnetic wave for 3D photonic crystals (recommended).

Uses H-field formulation with transverse projection, ensuring ∇·H = 0. This is the MPB-compatible formulation for accurate band structure calculations.

Field Components

  • Basis: 2 polarization vectors (e₁, e₂) perpendicular to k+G
  • Matrix size: 2N×2N (vs 3N×3N for FullVectorEM)
  • No spurious longitudinal modes

Mathematical Formulation

For each plane wave G, the H field is expanded as: H_G = h₁(G) e₁(G) + h₂(G) e₂(G) where e₁, e₂ are orthonormal vectors perpendicular to k+G.

Example

lat = fcc_lattice(1.0)
geo = Geometry(lat, air, [(Sphere([0,0,0], 0.25), dielectric)])

# Recommended for band structure calculations
solver = Solver(TransverseEM(), geo, (16, 16, 16); cutoff=5)
bands = compute_bands(solver, kpath; bands=1:10)

See also: FullVectorEM, TEWave, TMWave

source
PhoXonic.FullVectorEMType
FullVectorEM <: PhotonicWave

Full vector electromagnetic wave for 3D photonic crystals.

Uses H-field formulation: ∇×(ε⁻¹∇×H) = (ω/c)² μ H

Field Components

  • Solved: Hx, Hy, H_z (3 components)
  • Physical modes: 2 transverse modes per k-point

Notes

  • Produces spurious longitudinal modes at ω ≈ 0 (unphysical)
  • Use shift parameter with iterative solvers to skip spurious modes

Example

lat = cubic_lattice(1.0)
geo = Geometry(lat, air, [(Sphere([0,0,0], 0.3), rod)])
solver = Solver(FullVectorEM(), geo, (16, 16, 16); cutoff=3)

# Iterative solvers require shift to skip spurious modes
solver = Solver(FullVectorEM(), geo, (16, 16, 16), KrylovKitMethod(shift=0.01); cutoff=3)
solver = Solver(FullVectorEM(), geo, (16, 16, 16), LOBPCGMethod(shift=0.01); cutoff=3)

See also: TEWave, TMWave, FullElastic

source
PhoXonic.FullElasticType
FullElastic <: PhononicWave

Full elastic wave for 3D phononic crystals.

Solves the elastodynamic equation: ∇·σ = -ρω²u where σ = C:ε

Field Components

  • Solved: ux, uy, u_z (3 components)
  • Produces 3 bands per k-point (1 quasi-P + 2 quasi-S)

Notes

  • Eigenvalues ω² can be O(10¹⁰) for typical materials
  • Use shift parameter with iterative solvers for 3D

Example

lat = cubic_lattice(1.0)
steel = IsotropicElastic(ρ=7800.0, λ=115e9, μ=82e9)
epoxy = IsotropicElastic(ρ=1180.0, λ=4.43e9, μ=1.59e9)
geo = Geometry(lat, epoxy, [(Sphere([0,0,0], 0.3), steel)])

solver = Solver(FullElastic(), geo, (16, 16, 16); cutoff=3)

# Iterative solvers may require shift for 3D
solver = Solver(FullElastic(), geo, (16, 16, 16), KrylovKitMethod(shift=0.01); cutoff=3)
solver = Solver(FullElastic(), geo, (16, 16, 16), LOBPCGMethod(shift=0.01); cutoff=3)

See also: SHWave, PSVWave, FullVectorEM

source

Solver

Solver Types

PhoXonic.AbstractSolverType
AbstractSolver

Abstract base type for all solvers in PhoXonic.jl.

Subtypes:

  • Solver: Plane wave expansion (PWE) solver

This abstract type allows for future extensions (e.g., GME solver).

source
PhoXonic.SolverType
Solver{D<:Dimension, W<:WaveType, M<:SolverMethod} <: AbstractSolver

Solver for computing band structures of photonic/phononic crystals using the plane wave expansion (PWE) method.

Fields

  • wave: Wave type (TE, TM, SH, etc.)
  • geometry: Crystal geometry
  • basis: Plane wave basis
  • resolution: Grid resolution for discretization
  • material_arrays: Pre-computed material arrays on the grid
  • method: Solver method (DenseMethod(), BasicRSCG(), etc.)
source

High-level API

PhoXonic.solveFunction
solve(solver::Solver, k; bands=1:10)

Solve the eigenvalue problem at wave vector k.

Arguments

  • solver: The solver
  • k: Wave vector (in units of reciprocal lattice vectors or absolute)
  • bands: Which bands to return (default: 1:10)

Returns

  • frequencies: Eigenfrequencies (sorted, positive)
  • modes: Corresponding eigenvectors
source
PhoXonic.solve_at_kFunction
solve_at_k(solver, k, method; bands=1:10, X0=nothing, P=nothing)

Solve eigenvalue problem at a single k-point with explicit control over method, initial vectors, and preconditioner. Returns only frequencies.

For eigenvectors, use solve_at_k_with_vectors instead.

Arguments

  • solver::Solver: Solver instance
  • k: Wave vector (2D: Vector{Float64}, 1D: Float64)
  • method::SolverMethod: Solver method (DenseMethod, LOBPCGMethod, etc.)

Keyword Arguments

  • bands: Range of bands to compute (default: 1:10)
  • X0: Initial guess for eigenvectors (default: nothing = auto)
    • Size: dim × nev where dim = matrix_dimension(solver) and nev = maximum(bands)
    • Type: Matrix{ComplexF64}
    • If nothing, random orthonormal vectors are used
  • P: Preconditioner (default: nothing = use method's preconditioner setting)
    • Must implement ldiv!(y, P, x)

Returns

  • Vector{Float64}: Frequencies for the requested bands

Example

solver = Solver(PSVWave(), geo, (64, 64); cutoff=20)
dim = matrix_dimension(solver)  # 2514 for cutoff=20

# Frequencies only
freqs = solve_at_k(solver, k, LOBPCGMethod(); bands=1:20)

# With eigenvectors - use solve_at_k_with_vectors
freqs, vecs = solve_at_k_with_vectors(solver, k, DenseMethod(); bands=1:20)

# Custom preconditioner
using LinearAlgebra
LHS, _ = build_matrices(solver, k)
P = Diagonal(1.0 ./ diag(LHS))
freqs = solve_at_k(solver, k, LOBPCGMethod(); bands=1:20, P=P)

See also: solve_at_k_with_vectors, solve, matrix_dimension

source
solve_at_k(solver, k, method; kwargs...)

Solve eigenvalue problem at a single k-point.

See concrete method signatures for detailed documentation and keyword arguments.

source

Mid-level API

For custom algorithms and advanced analysis:

PhoXonic.solve_at_k_with_vectorsFunction
solve_at_k_with_vectors(solver, k, method; bands=1:10, X0=nothing, P=nothing)

Solve eigenvalue problem at a single k-point and return both frequencies and eigenvectors.

This function always returns eigenvectors, unlike solve_at_k which returns only frequencies by default. Use this when you need the eigenvectors for further analysis (e.g., computing overlaps, mode profiles, or topological invariants).

Arguments

  • solver::AbstractSolver: The solver object
  • k: Wave vector (Real for 1D, AbstractVector for 2D/3D)
  • method::SolverMethod: Solver method (DenseMethod(), KrylovKitMethod(), LOBPCGMethod())

Keyword Arguments

  • bands: Range of bands to compute (default: 1:10)
  • X0: Initial guess for eigenvectors (for LOBPCG, optional)
  • P: Preconditioner (for LOBPCG, optional)

Returns

  • (frequencies, eigenvectors): Tuple of frequencies (Vector{Float64}) and eigenvectors (Matrix{ComplexF64})

Example

ω, vecs = solve_at_k_with_vectors(solver, [0.1, 0.2], DenseMethod(); bands=1:4)
W = get_weight_matrix(solver)
overlap = vecs' * W * vecs  # Should be ≈ I (orthonormal)

See also: solve_at_k, get_weight_matrix, build_matrices

source
solve_at_k_with_vectors(solver, k, method; kwargs...)

Solve eigenvalue problem and return both frequencies and eigenvectors.

See concrete method signatures for detailed documentation and keyword arguments.

source
PhoXonic.build_matricesFunction
build_matrices(solver::Solver, k)

Build the eigenvalue problem matrices LHS * ψ = ω² * RHS * ψ for wave vector k.

source
build_matrices(solver::Solver{Dim3, FullVectorEM}, k)

Build 3N×3N matrices for 3D photonic crystal (H-field formulation).

Formulation: curl × ε⁻¹ × curl H = (ω²/c²) μ H

Storage order: [Hx; Hy; H_z] (component-order for FFT efficiency)

source
build_matrices(solver::Solver{Dim3, TransverseEM}, k)

Build 2N×2N matrices for 3D photonic crystal with transverse projection.

This is the MPB-compatible formulation that enforces ∇·H = 0 by projecting the 3N×3N FullVectorEM matrices onto the 2N-dimensional transverse subspace.

Algorithm

  1. Build 3N×3N matrices (LHS3N, RHS3N) as in FullVectorEM
  2. Construct polarization basis P (3N × 2N) where each column is a transverse polarization vector perpendicular to k+G
  3. Project: LHS = P' * LHS3N * P, RHS = P' * RHS3N * P

Advantages over FullVectorEM

  • Matrix size: 2N×2N (vs 3N×3N) → ~2.25× faster solve
  • No spurious longitudinal modes (ω ≈ 0)
  • Correct physical results matching MPB
source
build_matrices(solver::Solver{Dim3, FullElastic}, k)

Build 3N×3N matrices for 3D phononic crystal.

Formulation: -∇·(C:∇u) = ω² ρ u

Storage order: [ux; uy; u_z] (component-order)

source
build_matrices(solver, k)

Build the LHS and RHS matrices for the generalized eigenvalue problem.

See concrete method signatures for detailed documentation.

source
PhoXonic.get_weight_matrixFunction
get_weight_matrix(solver::Solver{Dim3, TransverseEM})

Not directly available for TransverseEM due to k-dependency.

The weight matrix for TransverseEM is the projected μ matrix: W = P' * μ_3N * P, where P is the polarization basis that depends on k. Since this varies with each k-point, a single k-independent weight matrix cannot be provided.

Alternative

Use the RHS matrix returned by build_matrices(solver, k) instead:

LHS, RHS = build_matrices(solver, k)
# RHS is the 2N×2N projected weight matrix for this k-point

See also

source
get_weight_matrix(solver)

Return the weight matrix W for computing inner products of eigenvectors.

The weight matrix satisfies: vecs' * W * vecs ≈ I for orthonormal eigenvectors. This is the RHS matrix of the generalized eigenvalue problem, which defines the inner product in the function space.

See concrete method signatures for detailed documentation.

source

Utilities

PhoXonic.matrix_dimensionFunction
matrix_dimension(solver::Solver) -> Int

Return the dimension of the eigenvalue problem matrix.

This is the size of the matrices A and B in the generalized eigenvalue problem A x = λ B x. Use this to determine the size of initial vectors for iterative solvers like LOBPCG.

Returns

  • num_pw for scalar waves (TE, TM, SH, 1D waves)
  • 2 * num_pw for 2D vector waves (P-SV)
  • 3 * num_pw for 3D vector waves (FullVectorEM, FullElastic)

Example

solver = Solver(PSVWave(), geo, (64, 64); cutoff=20)
dim = matrix_dimension(solver)  # 2 * num_pw

# Create initial vectors for LOBPCG
X0 = randn(ComplexF64, dim, 20)
X0, _ = qr(X0)
X0 = Matrix(X0)
source
matrix_dimension(solver)

Return the matrix dimension for the eigenvalue problem.

See concrete method signatures for detailed documentation.

source
PhoXonic.group_velocityFunction
group_velocity(solver::Solver, k; bands=1:10, δk=1e-5)

Compute group velocity v_g = ∂ω/∂k at wave vector k.

Arguments

  • solver: The solver
  • k: Wave vector
  • bands: Which bands to compute (default: 1:10)
  • δk: Finite difference step size (default: 1e-5)

Returns

  • v_g: Group velocity vectors for each band
    • 2D: Vector of SVector{2} (vgx, vgy)
    • 1D: Vector of Float64
source

Solver Method Types

Abstract Types

PhoXonic.RSCGMethodType
RSCGMethod <: IterativeMethod

Abstract type for Reduced Shifted Conjugate Gradient methods. Used for Green's function calculations.

source

Concrete Methods

PhoXonic.DenseMethodType
DenseMethod <: SolverMethod

Dense matrix eigenvalue solver using LAPACK's eigen function.

This is the default method, suitable for small to medium-sized systems. It builds the full N×N matrices and computes all eigenvalues/eigenvectors.

Fields

  • shift::Float64: Minimum eigenvalue cutoff (default: 0.0)
    • Eigenvalues with ω² < shift are filtered out after solving
    • Use shift > 0 for 3D FullVectorEM to skip spurious longitudinal modes (λ ≈ 0)
    • Recommended: shift = 0.01 for 3D photonic crystals
    • Note: Unlike iterative methods, this is post-hoc filtering, not shift-and-invert

Complexity

  • Memory: O(N²)
  • Time: O(N³)

where N = numplanewaves × ncomponents(wave)

Recommended For

  • 1D systems (any resolution)
  • 2D systems with N < ~1000
  • Development and debugging (exact results)

Example

# Default (implicit)
solver = Solver(TEWave(), geo, (64, 64))

# Explicit
solver = Solver(TEWave(), geo, (64, 64), DenseMethod())

# 3D with shift to skip spurious modes
solver = Solver(FullVectorEM(), geo, (12, 12, 12), DenseMethod(shift=0.01))

See also: KrylovKitMethod, LOBPCGMethod

source
PhoXonic.BasicRSCGType
BasicRSCG <: RSCGMethod

Basic RSCG method for Green's function computation. Suitable for DOS/LDOS calculations in large systems.

source
PhoXonic.KrylovKitMethodType
KrylovKitMethod <: IterativeMethod

Iterative eigenvalue solver using KrylovKit.jl. Suitable for large systems where dense matrix storage is impractical.

Uses shift-and-invert spectral transformation when shift > 0:

  • Original problem: A x = λ B x
  • Transformed: (A - σB)⁻¹ B x = μ x where μ = 1/(λ - σ)

This is particularly useful for 3D calculations where the H-field formulation produces N spurious longitudinal modes with λ ≈ 0. Setting shift > 0 (e.g., 0.01) effectively filters out these unphysical modes and focuses on the transverse (physical) modes.

Fields

  • tol::Float64: Convergence tolerance (default: 1e-8)
  • maxiter::Int: Maximum iterations (default: 300)
  • krylovdim::Int: Krylov subspace dimension (default: 30)
  • verbosity::Int: Output verbosity level (0=silent, 1=warn, 2=info)
  • shift::Float64: Spectral shift σ for shift-and-invert (default: 0.0)
    • shift = 0: Standard generalized eigenvalue problem (no transformation)
    • shift > 0: Shift-and-invert, finds eigenvalues closest to σ
    • Recommended: shift = 0.01 for 3D photonic crystals to skip longitudinal modes
  • matrix_free::Bool: Use matrix-free operators for shift-and-invert (default: false)
    • false: Build dense matrices (faster for N < 2000)
    • true: Use O(N) memory with iterative inner solver (for large 3D problems)

Example

# 2D calculation (no shift needed)
method = KrylovKitMethod()

# 3D calculation (use shift to skip longitudinal modes)
method = KrylovKitMethod(shift=0.01)

# Target specific frequency range
method = KrylovKitMethod(shift=1.5)  # Find modes near ω² = 1.5

See also: DenseMethod, Solver

source
PhoXonic.LOBPCGMethodType
LOBPCGMethod <: IterativeMethod

Locally Optimal Block Preconditioned Conjugate Gradient method. Based on Knyazev (2001), SIAM J. Sci. Comput. Vol.23, No.2, pp.517-541.

Solves the symmetric generalized eigenvalue problem A x = λ B x where both A and B are Hermitian and B is positive definite. This method is particularly effective for:

  • Large-scale band structure calculations with compute_bands
  • Phononic crystal calculations (steel/epoxy, silicon/air, etc.)

With warm start enabled (default), LOBPCG reuses eigenvectors from previous k-points as initial guesses, achieving significant speedup while maintaining accuracy. The first k-point uses Dense (if first_dense=true) to provide accurate initial eigenvectors.

Fields

  • tol::Float64: Convergence tolerance (default: 1e-3)
  • maxiter::Int: Maximum iterations (default: 100)
  • shift::Float64: Spectral shift for shift-and-invert transformation (default: 0.0)
  • warm_start::Bool: Use previous eigenvectors as initial guess (default: true)
  • scale::Bool: Scale matrix A by max|A| (default: false, scaling can hurt phononic problems)
  • first_dense::Bool: Solve first k-point with Dense for accurate warm start (default: true)
  • preconditioner: Preconditioner type (default: :none)

Notes

  • Requires symmetric/Hermitian matrices A and B
  • B must be positive definite
  • Works directly with dense matrices (no matrix-free support yet)
  • Block method: computes multiple eigenvalues simultaneously
  • For phononic problems, leave scale=false and preconditioner=:none (defaults)
  • When shift > 0, uses shift-and-invert to skip eigenvalues near zero (useful for 3D H-field formulation where spurious modes exist at λ ≈ 0)
  • Note: With shift > 0, the method falls back to dense eigen solver internally because the shifted matrix (A - σB) is not positive definite as required by LOBPCG.

Example

# Default with warm start (recommended for band structure)
solver = Solver(PSVWave(), geo, (64, 64), LOBPCGMethod(); cutoff=20)

# For 3D FullVectorEM (skip spurious modes at λ ≈ 0)
method = LOBPCGMethod(shift=0.01)

See also: KrylovKitMethod, DenseMethod, matrix_dimension

source

Band Structure

PhoXonic.BandStructureType
BandStructure

Results of a band structure calculation.

Fields

  • kpoints: K-points used (as vector of SVectors)
  • distances: Cumulative distances along the path
  • frequencies: Matrix of frequencies (nk × nbands)
  • labels: High-symmetry point labels (index, label)
source
PhoXonic.compute_bandsFunction
compute_bands(solver::Solver, kpath; bands=1:10, verbose=false, track_bands=false)

Compute band structure along a k-point path.

When the solver uses LOBPCGMethod with warm_start=true, automatically:

  1. Solves first k-point with Dense (if first_dense=true)
  2. Uses previous eigenvectors as initial guess for subsequent k-points
  3. Applies matrix scaling (if scale=true)

This can achieve up to 38x speedup for large problems.

Arguments

  • solver: The solver
  • kpath: K-point path (SimpleKPath, KPathInterpolant, or Vector)
  • bands: Which bands to compute (default: 1:10)
  • verbose: Print progress (default: false)
  • track_bands: Track bands across k-points using eigenvector overlap (default: false). When true, bands are reordered at each k-point to maintain continuity, preventing spurious crossings from appearing as discontinuities.

Returns

A BandStructure object containing frequencies for each k-point and band.

source
compute_bands(solver::Solver{Dim3}, kpath::SimpleKPath{3}; bands=1:10, verbose=false, track_bands=false)

Compute 3D band structure along a k-point path.

Arguments

  • solver: 3D solver (FullVectorEM or FullElastic)
  • kpath: K-point path from simple_kpath_cubic, simple_kpath_fcc, etc.
  • bands: Which bands to compute (default: 1:10)
  • verbose: Print progress (default: false)
  • track_bands: Track bands across k-points using eigenvector overlap (default: false). When true, bands are reordered at each k-point to maintain continuity.

Returns

A BandStructure{3} object containing frequencies for each k-point and band.

Example

lat = fcc_lattice(1.0)
geo = Geometry(lat, Dielectric(1.0), [(Sphere([0,0,0], 0.25), Dielectric(12.0))])
solver = Solver(FullVectorEM(), geo, (12,12,12), KrylovKitMethod(shift=0.01); cutoff=3)
kpath = simple_kpath_fcc(a=1.0, npoints=20)
bands = compute_bands(solver, kpath; bands=1:6, track_bands=true)

Note

At Γ point (k=0), the lowest transverse modes also have ω→0, which can cause anomalous values. Consider using a small k offset or skipping the Γ point.

source
compute_bands(solver, kpath; kwargs...)

Compute band structure along a k-path.

See concrete method signatures for detailed documentation and keyword arguments.

source
PhoXonic.find_bandgapFunction
find_bandgap(bs::BandStructure, band1::Int, band2::Int)

Find the band gap between band1 and band2.

Returns

Named tuple with:

  • gap: Band gap width (0 if bands overlap)
  • gap_ratio: Gap-to-midgap ratio
  • min_upper: Minimum of upper band
  • max_lower: Maximum of lower band
source
PhoXonic.find_all_gapsFunction
find_all_gaps(bs::BandStructure; threshold=0.0)

Find all band gaps in the band structure.

Returns

Vector of named tuples with gap information.

source

K-path

Simple K-paths

PhoXonic.SimpleKPathType
SimpleKPath{D}

A simple k-path structure for basic band structure calculations. Use this when Brillouin.jl's full functionality is not needed.

source
PhoXonic.simple_kpath_squareFunction
simple_kpath_square(; a=1.0, npoints=50)

Create a simple k-path for square lattice without Brillouin.jl. Γ → X → M → Γ

source
PhoXonic.simple_kpath_cubicFunction
simple_kpath_cubic(; a=1.0, npoints=50)

Create a simple k-path for simple cubic (SC) lattice.

Path: Γ → X → M → Γ → R → X

High-symmetry points (in units of 2π/a):

  • Γ = (0, 0, 0)
  • X = (1/2, 0, 0)
  • M = (1/2, 1/2, 0)
  • R = (1/2, 1/2, 1/2)

Arguments

  • a: Lattice constant (default: 1.0)
  • npoints: Number of points per segment (default: 50)

Returns

A SimpleKPath{3} object that can be passed to compute_bands.

Example

lat = cubic_lattice(1.0)
geo = Geometry(lat, Dielectric(1.0), [(Sphere([0,0,0], 0.3), Dielectric(12.0))])
solver = Solver(FullVectorEM(), geo, (12,12,12), KrylovKitMethod(shift=0.01); cutoff=3)
kpath = simple_kpath_cubic(a=1.0, npoints=30)
bands = compute_bands(solver, kpath; bands=1:6)

See also: kpath_cubic for Brillouin.jl-based k-path, Brillouin.jl documentation

source
PhoXonic.simple_kpath_fccFunction
simple_kpath_fcc(; a=1.0, npoints=50)

Create a simple k-path for face-centered cubic (FCC) lattice.

Path: Γ → X → W → L → Γ → K

High-symmetry points (values shown are multiplied by 2π/a to give Cartesian coordinates):

  • Γ = (0, 0, 0)
  • X = (0, 1, 0) × 2π/a
  • W = (1/2, 1, 0) × 2π/a
  • L = (1/2, 1/2, 1/2) × 2π/a
  • K = (3/4, 3/4, 0) × 2π/a

Arguments

  • a: Conventional lattice constant (default: 1.0)
  • npoints: Number of points per segment (default: 50)

Returns

A SimpleKPath{3} object that can be passed to compute_bands.

Example

lat = fcc_lattice(1.0)
geo = Geometry(lat, Dielectric(1.0), [(Sphere([0,0,0], 0.25), Dielectric(12.0))])
solver = Solver(FullVectorEM(), geo, (12,12,12), KrylovKitMethod(shift=0.01); cutoff=3)
kpath = simple_kpath_fcc(a=1.0, npoints=30)
bands = compute_bands(solver, kpath; bands=1:6)

See also: kpath_fcc for Brillouin.jl-based k-path, Brillouin.jl documentation

source
PhoXonic.simple_kpath_bccFunction
simple_kpath_bcc(; a=1.0, npoints=50)

Create a simple k-path for body-centered cubic (BCC) lattice.

Path: Γ → H → N → Γ → P → H

High-symmetry points (values shown are multiplied by 2π/a to give Cartesian coordinates):

  • Γ = (0, 0, 0)
  • H = (0, 0, 1) × 2π/a
  • N = (0, 1/2, 1/2) × 2π/a
  • P = (1/4, 1/4, 1/4) × 2π/a

Arguments

  • a: Conventional lattice constant (default: 1.0)
  • npoints: Number of points per segment (default: 50)

Returns

A SimpleKPath{3} object that can be passed to compute_bands.

Example

kpath = simple_kpath_bcc(a=1.0, npoints=30)
bands = compute_bands(solver, kpath; bands=1:6)

See also: kpath_bcc for Brillouin.jl-based k-path, Brillouin.jl documentation

source

Brillouin.jl Integration

PhoXonic.kpath_from_brillouinFunction
kpath_from_brillouin(sgnum::Int, Rs; N=100)

Create a k-path using Brillouin.jl's irrfbz_path.

Arguments

  • sgnum: Space group number (e.g., 1 for P1, 227 for diamond)
  • Rs: Direct lattice vectors as a matrix or vector of vectors
  • N: Number of interpolation points (default: 100)

Returns

A KPathInterpolant that can be iterated over.

Example

# Square lattice (space group 1, primitive cell)
Rs = [[1.0, 0.0], [0.0, 1.0]]
kpi = kpath_from_brillouin(1, Rs; N=100)

# Iterate over k-points
for k in kpi
    # k is an SVector
end
source
PhoXonic.kpath_squareFunction
kpath_square(; a=1.0, N=100)

Create a standard k-path for a 2D square lattice: Γ → X → M → Γ

Uses Brillouin.jl with space group 1 (P1).

source
PhoXonic.kpath_hexagonalFunction
kpath_hexagonal(; a=1.0, N=100)

Create a standard k-path for a 2D hexagonal lattice: Γ → M → K → Γ

Uses Brillouin.jl with appropriate space group.

source
PhoXonic.kpath_cubicFunction
kpath_cubic(; a=1.0, N=100)

Create a k-path for simple cubic (SC) lattice using Brillouin.jl.

Path: Γ → X → M → Γ → R → X | R → M

Uses space group 221 (Pm-3m).

Arguments

  • a: Lattice constant (default: 1.0)
  • N: Number of interpolation points (default: 100)

Returns

A KPathInterpolant from Brillouin.jl that can be passed to compute_bands.

Example

lat = cubic_lattice(1.0)
geo = Geometry(lat, Dielectric(1.0), [(Sphere([0,0,0], 0.3), Dielectric(12.0))])
solver = Solver(FullVectorEM(), geo, (12,12,12), KrylovKitMethod(shift=0.01); cutoff=3)
kpath = kpath_cubic(a=1.0, N=100)
bands = compute_bands(solver, kpath; bands=1:6)

Reference: Setyawan & Curtarolo, Comp. Mat. Sci. 49, 299 (2010). DOI:10.1016/j.commatsci.2010.05.010

See also: Brillouin.jl irrfbz_path

source
PhoXonic.kpath_fccFunction
kpath_fcc(; a=1.0, N=100)

Create a k-path for face-centered cubic (FCC) lattice using Brillouin.jl.

Path: Γ → X → U | K → Γ → L → W → X

Uses space group 225 (Fm-3m).

Arguments

  • a: Conventional lattice constant (default: 1.0)
  • N: Number of interpolation points (default: 100)

Returns

A KPathInterpolant from Brillouin.jl that can be passed to compute_bands.

Example

lat = fcc_lattice(1.0)
geo = Geometry(lat, Dielectric(1.0), [(Sphere([0,0,0], 0.25), Dielectric(12.0))])
solver = Solver(FullVectorEM(), geo, (12,12,12), KrylovKitMethod(shift=0.01); cutoff=3)
kpath = kpath_fcc(a=1.0, N=100)
bands = compute_bands(solver, kpath; bands=1:6)

Reference: Setyawan & Curtarolo, Comp. Mat. Sci. 49, 299 (2010). DOI:10.1016/j.commatsci.2010.05.010

See also: Brillouin.jl irrfbz_path

source
PhoXonic.kpath_bccFunction
kpath_bcc(; a=1.0, N=100)

Create a k-path for body-centered cubic (BCC) lattice using Brillouin.jl.

Path: Γ → H → N → Γ → P → H | P → N

Uses space group 229 (Im-3m).

Arguments

  • a: Conventional lattice constant (default: 1.0)
  • N: Number of interpolation points (default: 100)

Returns

A KPathInterpolant from Brillouin.jl that can be passed to compute_bands.

Example

kpath = kpath_bcc(a=1.0, N=100)
bands = compute_bands(solver, kpath; bands=1:6)

Reference: Setyawan & Curtarolo, Comp. Mat. Sci. 49, 299 (2010). DOI:10.1016/j.commatsci.2010.05.010

See also: Brillouin.jl irrfbz_path

source