API
Configuration
SHTns.SHTnsCfg — Typemutable struct SHTnsCfg{TR<:Union{Real,Complex}, N<:SHTnsNorm, T<:SHTnsType}cfg::Ptr{shtns_info}norm::Nshtype::Trobert_form::Boolnlm::Intlmax::Intmmax::Intmres::Intnlat_2::Intnlat::Intnphi::Intnspat::Intli::Vector{Int}mi::Vector{Int}ct::Vector{Float64}st::Vector{Float64}nlat_padded::Intnlm_cplx::Inthowmany::Int
Configuration of spherical harmonic transform.
Transform types (SHTnsType)
SHTns.SHTnsType — Typeabstract type SHTnsTypeSHTnsType is an abstract type for the spherical harmonic transform types. All subtypes contain the following keyword arguments:
contiguous_lat::Bool=truecontiguous_phi::Bool=falsepadding::Bool=falsegpu::Bool=falsesouthpolefirst::Bool=falsefloat32::Bool=false
SHTns.Gauss — Typestruct SHTns.Gauss <: SHTnsTypeSHTns.GaussFly — Typestruct SHTns.GaussFly <: SHTnsTypeSHTns.QuickInit — Typestruct SHTns.QuickInit <: SHTnsTypeSHTns.RegDCT — Typestruct SHTns.RegDCT <: SHTnsTypeSHTns.RegFast — Typestruct SHTns.RegFast <: SHTnsTypeSHTns.RegPoles — Typestruct SHTns.RegPoles <: SHTnsTypeSpherical harmonic normalizations (SHTnsNorm)
SHTns.ForRotations — Typestruct SHTns.ForRotations <: SHTnsNormSHTns.FourPi — Typestruct SHTns.FourPi <: SHTnsNormSHTns.Orthonormal — Typestruct SHTns.Orthonormal <: SHTnsNormSHTns.Schmidt — Typestruct SHTns.Schmidt <: SHTnsNormTransforms
The transforms from spatial space to spectral space, also called analysis, are available through SHTns.analys and SHTns.analys!. The transforms from spectral space to spatial space, also called synthesis, are available through SHTns.synth and SHTns.synth!.
SHTns.analys — Functionanalys(cfg::SHTnsCfg, v)
analys(cfg::SHTnsCfg, utheta, uphi)
analys(cfg::SHTnsCfg, ur, utheta, uphi)Transforms the spatial data into spherical harmonics coefficients qlm; slm and tlm; qlm, slm and tlm for scalar; 2D; 3D fields, respectively.
SHTns.analys! — Functionanalys!(cfg::SHTnsCfg, v, qlm)
analys!(cfg::SHTnsCfg, utheta, uphi, slm, tlm)
analys!(cfg::SHTnsCfg, ur, utheta, uphi, qlm, slm, tlm)In-place transforms of the spatial data into spherical harmonics coefficients for scalar, 2D or 3D fields.
This function modifies the input arrays v, ur, utheta and uphi.
SHTns.synth — Functionsynth(cfg::SHTnsCfg, qlm)
synth(cfg::SHTnsCfg, slm, tlm)
synth(cfg::SHTnsCfg, qlm, slm, tlm)Transforms spherical harmonics coefficients qlm into spatial data v; slm and tlm into spatial data utheta and uphi; qlm, slm and tlm into spatial data ur, utheta and uphi.
SHTns.synth! — Functionsynth!(cfg::SHTnsCfg, qlm, v)
synth!(cfg::SHTnsCfg, slm, tlm, utheta, uphi)
synth!(cfg::SHTnsCfg, qlm, slm, tlm, ur, utheta, uphi)In-place transforms of the spherical harmonics coefficients into spatial data for scalar, 2D or 3D fields.
Miscellaneous
SHTns.LM — FunctionLM(cfg::SHTnsCfg{Real,T,N}, l, m)Returns index corresponding to l and m in real spherical harmonic expansion.
LM(cfg::SHTnsCfg{Complex,T,N}, l, m)Returns index corresponding to l and m in complex spherical harmonic expansion.
SHTns.gauss_weights — Functiongauss_weights(cfg::SHTnsCfg)Returns Gauss quadrature weights of length cfg.nlat (shtns_gauss_wts only returns half of the symmetric weights).
SHTns.grid — Functiongrid(cfg::SHTnsCfg)Returns latitudes lat::Vector{Float64} and longitudes lon::Vector{Float64}.