API Reference
API Reference
Complete function and type documentation
Complete reference for all SHTnsKit.jl functions and types.
Configuration Management
Configuration Creation
SHTnsKit.create_config — Function
create_config(lmax::Int; mmax=lmax, mres=1, nlat=lmax+2, nlon=max(2*lmax+1,4),
norm::Symbol=:orthonormal, cs_phase::Bool=true,
real_norm::Bool=false, robert_form::Bool=false,
grid_type::Symbol=:gauss) -> SHTConfigCompatibility wrapper for configuration creation used in some docs/snippets. Supports Gauss–Legendre (grid_type = :gauss) and regular equiangular (grid_type = :regular or :regular_poles) grids, forwarding to the appropriate creator. nlat/nlon defaults are adjusted to satisfy accuracy constraints for the chosen grid.
create_config(::Type{T}, lmax::Int, nlat::Int, mres::Int=1; kwargs...) -> SHTConfigCompatibility method that ignores the element type T and calls create_config. Provided to match older example signatures like create_config(Float64, 20, 22, 1).
SHTnsKit.create_gauss_config — Function
create_gauss_config(lmax::Int, nlat::Int; mmax::Int=lmax, nlon::Int=max(2*lmax+1, 4)) -> SHTConfigCreate a Gauss–Legendre based SHT configuration. Constraints:
nlat ≥ lmax+1for exactness up tolmaxin θ integration.nlon ≥ 2*mmax+1to resolve azimuthal orders up tommax.
SHTnsKit.create_regular_config — Function
create_regular_config(lmax::Int, nlat::Int; mmax::Int=lmax, mres::Int=1,
nlon::Int=max(2*lmax+1, 4), norm::Symbol=:orthonormal,
cs_phase::Bool=true, real_norm::Bool=false,
robert_form::Bool=false, include_poles::Bool=false,
precompute_plm::Bool=true, use_dh_weights::Bool=false) -> SHTConfigCreate an equiangular (regular) grid configuration. Regular grids use simple θ = (i+0.5)π/nlat nodes by default; set include_poles=true to place nodes directly on the poles. By default associated Legendre tables are precomputed, which mirrors SHTns' regular-grid behaviour and improves performance.
Driscoll-Healy Quadrature
Set use_dh_weights=true to use Driscoll-Healy quadrature for exact spherical harmonic transforms. This requires:
include_poles=truenlat = 2*(lmax+1)for exactness up to degree lmaxnlatmust be even
The DH grid uses θ = πj/n for j=0,...,n-1 (includes north pole, excludes south pole) and provides exact quadrature via specially computed weights.
Reference: Driscoll & Healy (1994), "Computing Fourier transforms and convolutions on the 2-sphere", Adv. Appl. Math., 15, 202-250.
SHTnsKit.destroy_config — Function
destroy_config(cfg::SHTConfig)No-op placeholder for API symmetry with libraries that require explicit teardown.
Index Utilities
SHTnsKit.nlm_calc — Function
nlm_calc(lmax::Integer, mmax::Integer, mres::Integer) -> IntCalculate the total number of packed spectral coefficients for a real scalar field. Real fields only need m ≥ 0 modes due to Hermitian symmetry: Yl^{-m} = (-1)^m (Yl^m)*
The packing scheme groups modes by m-order, then by l-degree within each m block. This layout is optimized for vectorized operations and cache efficiency.
SHTnsKit.LM_index — Function
LM_index(lmax::Int, mres::Int, l::Int, m::Int) -> IntConvert 2D spherical harmonic indices (l,m) to a 1D packed array index. This implements the SHTns packing convention for efficient memory layout.
The packing scheme groups coefficients by m-order blocks, with each block containing all valid l-degrees for that m. This enables vectorized operations on modes with the same azimuthal symmetry.
Matches SHTns C macro: LM(shtns, l, m)
Scalar Field Transforms
Forward Transform (Synthesis)
SHTnsKit.synthesis — Function
synthesis(cfg::SHTConfig, alm::AbstractMatrix; real_output::Bool=true) -> MatrixInverse transform back to a grid (nlat, nlon). If real_output=true, Hermitian symmetry is enforced before IFFT. Optional fft_scratch lets you reuse a preallocated (nlat,nlon) complex buffer for lower allocations.
Parallelizes over m-modes using static scheduling for consistent performance.
SHTnsKit.synthesis! — Function
synthesis!(plan::SHTPlan, f_out::AbstractMatrix, alm::AbstractMatrix; real_output=true)In-place inverse scalar SHT writing spatial field into f_out. Streams m→k directly without building a (θ×m) intermediate.
synthesis!(cfg::SHTConfig, f_out::AbstractMatrix, alm::AbstractMatrix;
real_output::Bool=true, use_fused_loops::Bool=true)In-place inverse transform. Writes the spatial field into f_out to reduce allocations. f_out must be (nlat, nlon) and real if real_output=true. You may pass a complex fft_scratch buffer of size (nlat,nlon) to reuse FFT workspace.
Backward Transform (Analysis)
SHTnsKit.analysis — Function
analysis(cfg::SHTConfig, f::AbstractMatrix) -> Matrix{ComplexF64}Forward transform on Gauss–Legendre × equiangular grid. Returns coefficients alm[l+1, m+1] with orthonormal normalization.
Parallelizes over m-modes using static scheduling for consistent performance.
SHTnsKit.analysis! — Function
analysis!(plan::SHTPlan, alm_out::AbstractMatrix, f::AbstractMatrix)In-place forward scalar SHT writing coefficients into alm_out.
analysis!(cfg::SHTConfig, alm_out::AbstractMatrix, f::AbstractMatrix; use_fused_loops=true)In-place forward transform. Writes coefficients into alm_out to avoid allocating the output matrix each call. alm_out must be size (lmax+1, mmax+1).
Complex Field Transforms
SHTnsKit.synthesis_packed_cplx — Function
synthesis_packed_cplx(cfg::SHTConfig, alm_packed::AbstractVector{<:Complex}) -> Matrix{ComplexF64}Synthesize complex spatial field from packed complex coefficients (LM_cplx order). Returns an nlat × nlon complex array.
SHTnsKit.analysis_packed_cplx — Function
analysis_packed_cplx(cfg::SHTConfig, z::AbstractMatrix{<:Complex}) -> Vector{ComplexF64}Analyze complex spatial field into packed complex coefficients (LM_cplx order). Input z must be nlat × nlon complex.
Vector Field Transforms
Vector fields on the sphere are decomposed into spheroidal and toroidal components:
- Spheroidal: Poloidal component (has radial component)
- Toroidal: Azimuthal component (purely horizontal)
Vector Synthesis
SHTnsKit.synthesis_sphtor — Function
synthesis_sphtor(cfg, Slm, Tlm; real_output=true) -> (Vt, Vp)Transform spheroidal/toroidal coefficients to horizontal vector field components. Returns colatitude (Vt) and azimuthal (Vp) components on the spatial grid.
Vector Analysis
SHTnsKit.analysis_sphtor — Function
analysis_sphtor(cfg, Vt, Vp) -> (Slm, Tlm)Transform horizontal vector field components to spheroidal/toroidal coefficients. Input: colatitude (Vt) and azimuthal (Vp) components on spatial grid. The returned coefficients satisfy δlm = −l(l+1) Slm (divergence) and ζlm = −l(l+1) Tlm (vorticity) when expressed in the internal normalization.
QST Transforms (3D Vector Fields)
SHTnsKit.synthesis_qst — Function
synthesis_qst(cfg, Qlm, Slm, Tlm; real_output=true) -> (Vr, Vt, Vp)Transform QST spectral coefficients to 3D spatial vector field components. Returns radial (Vr), colatitude (Vt), and azimuthal (Vp) components.
SHTnsKit.analysis_qst — Function
analysis_qst(cfg, Vr, Vt, Vp) -> (Qlm, Slm, Tlm)Transform 3D spatial vector field to QST spectral coefficients. Input: radial (Vr), colatitude (Vt), and azimuthal (Vp) components.
Rotations
SHTnsKit.SH_Zrotate — Function
SH_Zrotate(cfg::SHTConfig, Qlm::AbstractVector{<:Complex}, alpha::Real, Rlm::AbstractVector{<:Complex})Rotate a real-field SH expansion around the Z-axis by angle alpha. Input and output are packed Qlm vectors (LM order, m ≥ 0). In-place supported if Rlm === Qlm.
SHTnsKit.SH_Yrotate — Function
SH_Yrotate(cfg::SHTConfig, Qlm::AbstractVector{<:Complex}, alpha::Real, Rlm::AbstractVector{<:Complex})Rotate a real-field SH expansion around the Y-axis by angle alpha. Uses Wigner-d mixing per l; dispatches to the general rotation engine.
SHTnsKit.SH_Yrotate90 — Function
SH_Yrotate90(cfg::SHTConfig, Qlm::AbstractVector{<:Complex}, Rlm::AbstractVector{<:Complex})SHTnsKit.SH_Xrotate90 — Function
SH_Xrotate90(cfg::SHTConfig, Qlm::AbstractVector{<:Complex}, Rlm::AbstractVector{<:Complex})Rotate around X-axis by 90 degrees using ZYZ equivalence: Rz(π/2)·Ry(π/2)·Rz(-π/2).
Energy/Power Spectrum Analysis
SHTnsKit.energy_scalar_l_spectrum — Function
energy_scalar_l_spectrum(cfg, alm; real_field=true) -> Vector{Float64}Compute energy spectrum as a function of spherical harmonic degree l. Returns E(l) = Σₘ |a_lm|² for each l = 0..lmax.
SHTnsKit.energy_vector_l_spectrum — Function
energy_vector_l_spectrum(cfg, Slm, Tlm; real_field=true) -> Vector{Float64}Compute kinetic energy spectrum by degree l for vector fields. Returns KE(l) = Σₘ l(l+1)[|Slm|² + |Tlm|²] for each l = 1..lmax.
SHTnsKit.energy_scalar — Function
energy_scalar(cfg, alm; real_field=true) -> Float64Compute the total energy of a scalar field from its spherical harmonic coefficients.
For a scalar field f(θ,φ) = Σ alm Yl^m(θ,φ), the energy is defined as: E = (1/2) ∫ |f|² dΩ = (1/2) Σ |a_lm|²
This represents the L² norm of the field, which is conserved under orthonormal spherical harmonic transforms (Parseval's identity).
SHTnsKit.energy_vector — Function
energy_vector(cfg, Slm, Tlm; real_field=true) -> Float64Compute the total kinetic energy of a vector field from spheroidal/toroidal coefficients.
For a vector field V = ∇×(T Yl^m êᵣ) + ∇ₕ(S Yl^m), the kinetic energy is: KE = (1/2) ∫ |V|² dΩ = (1/2) Σ [l(l+1)|Slm|² + l(l+1)|Tlm|²]
SHTnsKit.enstrophy — Function
enstrophy(cfg, Tlm; real_field=true) -> Float64Compute the enstrophy (kinetic energy of vorticity) from toroidal coefficients.
For a toroidal field with coefficients Tlm, the enstrophy is: Z = (1/2) ∫ |∇×V|² dΩ = (1/2) Σ l²(l+1)²|Tlm|²
This is a measure of the small-scale intensity of rotational motion.
Threading Control
SHTnsKit.shtns_use_threads — Function
shtnsusethreads(num_threads) -> Int
Buffer Helpers
SHTnsKit.scratch_spatial — Function
scratch_spatial(cfg::SHTConfig, T::Type=Float64) -> Matrix{T}Allocate a spatial grid buffer of size (nlat, nlon).
SHTnsKit.scratch_fft — Function
scratch_fft(cfg::SHTConfig, T::Type=ComplexF64) -> Matrix{T}Allocate a complex buffer sized for φ-FFTs (nlat, nlon). Reuse this across analysis!/synthesis! calls via the fft_scratch keyword to avoid per-call allocations.
Distributed Transforms (MPI)
When using MPI with PencilArrays, the following functions are available via the parallel extension:
dist_analysis(cfg, fθφ)- Distributed spatial to spectral transformdist_synthesis(cfg, Alm; prototype_θφ, real_output)- Distributed spectral to spatial transformdist_analysis_sphtor(cfg, Vθ, Vφ)- Distributed vector analysisdist_synthesis_sphtor(cfg, Slm, Tlm; prototype_θφ)- Distributed vector synthesis
See the Distributed Guide for detailed usage.
Gradient and Differential Operators
SHTnsKit.synthesis_grad — Function
synthesis_grad(cfg::SHTConfig, Slm::AbstractMatrix; real_output::Bool=true)Gradient synthesis alias for compatibility with SHTns.
SHTnsKit.divergence_from_spheroidal — Function
divergence_from_spheroidal(cfg, Slm) -> MatrixReturn divergence spectral coefficients δlm corresponding to spheroidal modes via δlm = −l(l+1) Slm (`δ{00} = 0`).
SHTnsKit.vorticity_from_toroidal — Function
vorticity_from_toroidal(cfg, Tlm) -> MatrixReturn vorticity spectral coefficients ζlm = −l(l+1) Tlm.
Usage Examples
Basic Transform
using SHTnsKit
lmax = 16
nlat = lmax + 2
nlon = 2*lmax + 1
cfg = create_gauss_config(lmax, nlat; nlon=nlon)
# Create test pattern
spatial = zeros(cfg.nlat, cfg.nlon)
for i in 1:cfg.nlat
x = cfg.x[i]
spatial[i, :] .= (3*x^2 - 1)/2
end
# Transform roundtrip
Alm = analysis(cfg, spatial)
recovered = synthesis(cfg, Alm)
destroy_config(cfg)Vector Field Transform
using SHTnsKit
lmax = 32
nlat = lmax + 2
nlon = 2*lmax + 1
cfg = create_gauss_config(lmax, nlat; nlon=nlon)
# Create velocity field
Vθ = zeros(cfg.nlat, cfg.nlon)
Vφ = zeros(cfg.nlat, cfg.nlon)
for i in 1:cfg.nlat, j in 1:cfg.nlon
θ = cfg.θ[i]
φ = cfg.φ[j]
Vθ[i,j] = cos(θ) * sin(φ)
Vφ[i,j] = cos(φ)
end
# Decompose into spheroidal/toroidal
Slm, Tlm = analysis_sphtor(cfg, Vθ, Vφ)
# Reconstruct
Vθ_rec, Vφ_rec = synthesis_sphtor(cfg, Slm, Tlm)
destroy_config(cfg)Error Handling
Common Errors
BoundsError: Invalid lmax/mmax valuesAssertionError/DimensionMismatch: Array size mismatches
Best Practices
# Always check array sizes
@assert size(Alm) == (cfg.lmax+1, cfg.mmax+1) "Wrong spectral array size"
@assert size(spatial) == (cfg.nlat, cfg.nlon) "Wrong spatial array size"
# Always destroy configurations
cfg = create_gauss_config(32, 34; nlon=65)
# ... work with cfg ...
destroy_config(cfg)