SHTns 3.6.1
Data Structures | Macros | Typedefs | Enumerations | Functions
shtns.h File Reference

shtns.h is the definition file for SHTns : include this file in your source code to use SHTns. More...

#include <complex.h>

Go to the source code of this file.

Data Structures

struct  shtns_info
 

Macros

#define SHTNS_INTERFACE   0x30602
 SHTns interface version (loosely follow versions) allowing simple version checks: (major << 16) | (minor << 8) | patchlevel should be increased at least each time this file changes.
 
#define SHT_NO_CS_PHASE   (256*4)
 don't include Condon-Shortley phase (add to last argument of shtns_create)
 
#define SHT_REAL_NORM   (256*8)
 use a "real" normalization. (add to last argument of shtns_create)
 
#define SHT_NATIVE_LAYOUT   0
 Tells shtns_init to use Native layout.
 
#define SHT_THETA_CONTIGUOUS   256
 use Contiguous latitudes
 
#define SHT_PHI_CONTIGUOUS   (256*2)
 use Contiguous longitudes
 
#define SHT_SOUTH_POLE_FIRST   (256*32)
 latitudinal data are stored starting from south pole.
 
#define SHT_SCALAR_ONLY   (256*16)
 don't compute vector matrices. (add to flags in shtns_set_grid)
 
#define SHT_LOAD_SAVE_CFG   (256*64)
 try to load and save the config. (add to flags in shtns_set_grid)
 
#define SHT_ALLOW_GPU   (256*128)
 allows to use a GPU. This needs special care because the same plan cannot be used simultaneously by different threads anymore.
 
#define SHT_ALLOW_PADDING   (256*256)
 allows SHTns to add extra space between lines of the spatial array to avoid cache bank conflicts.
 
#define SHT_ROBERT_FORM   (256*512)
 use Robert form for vector transforms. See also shtns_robert_form
 
#define SHT_FP32   (256*1024)
 float support. Currently only on GPU, where it then replaces double (fp64) support which is likely to segfault.
 
#define NSPAT_ALLOC(shtns)   (shtns->nspat)
 total number of 'doubles' required for a spatial field (includes FFTW reserved space).
 
#define PHI_DEG(shtns, ip)   (360./(shtns->nphi*shtns->mres))*(ip)
 phi angle value in degrees for given index ip.
 
#define PHI_RAD(shtns, ip)   (2.*M_PI/(shtns->nphi*shtns->mres))*(ip)
 phi angle value in radians for given index ip.
 
#define SH_to_grad_spat(shtns, S, Gt, Gp)   SHsph_to_spat(shtns, S, Gt, Gp)
 Compute the spatial representation of the gradient of a scalar SH field.
 
#define SH_to_grad_spat_l(shtns, S, Gt, Gp, ltr)   SHsph_to_spat_l(shtns, S, Gt, Gp, ltr)
 Compute the spatial representation of the gradient of a scalar SH field.
 
#define SH_to_grad_spat_ml(shtns, im, S, Gt, Gp, ltr)   SHsph_to_spat_ml(shtns, im, S, Gt, Gp, ltr)
 Compute the spatial representation of the gradient of a scalar SH field.
 
Access to spherical harmonic components

The following macros give access to single spherical harmonic coefficient or perform loops spanning all of them.

#define LiM(shtns, l, im)   ( (((im)*(2*shtns->lmax + 2 - ((im)+1)*shtns->mres))>>1) + (l) )
 LiM(shtns, l,im) : macro returning array index for given l and im, corresponding to config shtns.
 
#define LM(shtns, l, m)   ( (((((unsigned short)(m))/shtns->mres)*(2*shtns->lmax + 2 - ((m)+shtns->mres)))>>1) + (l) )
 LM(shtns, l,m) : macro returning array index for given l and m, corresponding to config shtns.
 
#define LM_LOOP(shtns, action)   { int lm=0; do { action } while(++lm < shtns->nlm); }
 LM_LOOP( shtns, action ) : macro that performs "action" for every (l,m), with lm set, but neither l, m nor im.
 
#define LM_L_LOOP(shtns, action)   { int lm=0; do { int l=shtns->li[lm]; action } while(++lm < shtns->nlm); }
 LM_L_LOOP : loop over all (l,im) and perform "action" : l and lm are defined (but NOT m and im).
 
#define LM_cplx(shtns, l, m)   ( ((l) <= shtns->mmax) ? (l)*((l)+1)+(m) : shtns->mmax*(2*(l) - shtns->mmax) + (l)+(m) )
 LM_cplx(shtns, l,m) : macro returning array index for SH representation of complex-valued fields for given l and m. (mres=1 only!)
 

Typedefs

typedef _Complex double cplx
 double precision complex number data type
 
typedef _Complex float cplx_f
 single precision (float) complex number data type
 
typedef struct shtns_infoshtns_cfg
 pointer to data structure describing an SHT, returned by shtns_init() or shtns_create().
 
typedef struct shtns_rot_ * shtns_rot
 pointer to data structure describing a rotation, returned by shtns_rotation_create().
 

Enumerations

enum  shtns_norm { sht_orthonormal , sht_fourpi , sht_schmidt , sht_for_rotations }
 Different Spherical Harmonic normalizations. See also section Normalization for details. More...
 
enum  shtns_type {
  sht_gauss , sht_auto , sht_reg_fast , sht_reg_dct ,
  sht_quick_init , sht_reg_poles , sht_gauss_fly
}
 different SHT types and algorithms More...
 

Functions

long nlm_calc (long lmax, long mmax, long mres)
 compute number of spherical harmonics modes (l,m) to represent a real-valued spatial field for given size parameters. Does not require any previous setup.
 
long nlm_cplx_calc (long lmax, long mmax, long mres)
 compute number of spherical harmonics modes (l,m) to represent a complex-valued spatial field, for given size parameters. Does not require any previous setup.
 
void shtns_verbose (int)
 controls output during initialization: 0=no output (default), 1=some output, 2=degug (if compiled in)
 
void shtns_print_version (void)
 print version information to stdout.
 
const char * shtns_get_build_info (void)
 get version and build information as a string (the same as the one printed by shtns_print_version())
 
void shtns_print_cfg (shtns_cfg)
 print information about given config to stdout.
 
initialization
shtns_cfg shtns_init (enum shtns_type flags, int lmax, int mmax, int mres, int nlat, int nphi)
 Simple initialization of the spherical harmonic transforms of given size. Calls shtns_create and shtns_set_grid_auto.
 
shtns_cfg shtns_create (int lmax, int mmax, int mres, enum shtns_norm norm)
 Defines the sizes of the spectral description. Use for advanced initialization.
 
int shtns_set_grid (shtns_cfg, enum shtns_type flags, double eps, int nlat, int nphi)
 Precompute everything for a given spatial grid. Use for advanced initialization, after shtns_create.
 
int shtns_set_grid_auto (shtns_cfg, enum shtns_type flags, double eps, int nl_order, int *nlat, int *nphi)
 Precompute everything and choose the optimal nlat and nphi for a given non-linear order.
 
shtns_cfg shtns_create_with_grid (shtns_cfg, int mmax, int nofft)
 Copy a given config but allow a different (smaller) mmax and the possibility to enable/disable fft.
 
int shtns_use_threads (int num_threads)
 Enables multi-thread transform using OpenMP with num_threads (if available). Returns number of threads that will be used.
 
int shtns_set_batch (shtns_cfg shtns, int howmany, long spec_dist)
 beta: Modify plan to perform several transforms together (batch). This is useful to get good performance for small transforms on GPU, but also works on CPU.
 
void shtns_reset (void)
 destroy all configs, free memory, and go back to initial state.
 
void shtns_destroy (shtns_cfg)
 free memory of given config, which cannot be used afterwards.
 
void shtns_unset_grid (shtns_cfg)
 unset the grid.
 
void shtns_robert_form (shtns_cfg, int robert)
 set the use of Robert form. If robert != 0, the vector synthesis returns a field multiplied by sin(theta), while the analysis divides by sin(theta) before the transform.
 
void * shtns_malloc (size_t bytes)
 allocate appropriate CPU memory (pinned for gpu, aligned for avx, ...). Use shtns_free to free it.
 
void shtns_free (void *p)
 free memory allocated with shtns_malloc
 
special values
double sh00_1 (shtns_cfg)
 return the spherical harmonic representation of 1 (l=0,m=0)
 
double sh10_ct (shtns_cfg)
 return the spherical harmonic representation of cos(theta) (l=1,m=0)
 
double sh11_st (shtns_cfg)
 return the spherical harmonic representation of sin(theta)*cos(phi) (l=1,m=1)
 
double shlm_e1 (shtns_cfg, int l, int m)
 return the l,m SH coefficient corresponding to unit energy.
 
int shtns_gauss_wts (shtns_cfg, double *wts)
 fill the given array with Gauss weights. returns the number of weights written (0 if not a Gauss grid).
 
Rotation functions
void SH_Zrotate (shtns_cfg, cplx *Qlm, double alpha, cplx *Rlm)
 Rotate a SH representation Qlm around the z-axis by angle alpha (in radians), which is the same as rotating the reference frame by angle -alpha.
 
void SH_Yrotate (shtns_cfg, cplx *Qlm, double alpha, cplx *Rlm)
 rotate Qlm around Y axis by arbitrary angle, using composition of rotations. Store the result in Rlm.
 
void SH_Yrotate90 (shtns_cfg, cplx *Qlm, cplx *Rlm)
 rotate Qlm by 90 degrees around Y axis and store the result in Rlm.
 
void SH_Xrotate90 (shtns_cfg, cplx *Qlm, cplx *Rlm)
 rotate Qlm by 90 degrees around X axis and store the result in Rlm.
 
shtns_rot shtns_rotation_create (const int lmax, const int mmax, int norm)
 Creation of a spherical harmonic rotation object for degrees and orders up to lmax and mmax respectively.
 
void shtns_rotation_destroy (shtns_rot r)
 Release memory allocated for the rotation object.
 
void shtns_rotation_set_angles_ZYZ (shtns_rot r, double alpha, double beta, double gamma)
 Defines a rotation parameterized by the 3 intrinsic Euler angles with ZYZ convention.
 
void shtns_rotation_set_angles_ZXZ (shtns_rot r, double alpha, double beta, double gamma)
 Defines a rotation parameterized by the 3 intrinsic Euler angles with ZXZ convention.
 
void shtns_rotation_set_angle_axis (shtns_rot r, double theta, double Vx, double Vy, double Vz)
 Defines a rotation of angle theta around axis of cartesian coordinates (Vx,Vy,Vz).
 
int shtns_rotation_wigner_d_matrix (shtns_rot r, const int l, double *mx)
 Generate spherical-harmonic rotation matrix for given degree l and orthonormal convention (Wigner-d matrix)
 
void shtns_rotation_apply_cplx (shtns_rot r, cplx *Zlm, cplx *Rlm)
 apply rotation to the orthonormal spherical harmonic expansion Zlm of a complex-valued field, store the result into Rlm (can be the same as Zlm)
 
void shtns_rotation_apply_real (shtns_rot r, cplx *Qlm, cplx *Rlm)
 apply rotation to the orthonormal spherical harmonic expansion Qlm of a real field, store the result into Rlm (can be the same as Qlm)
 
Generation of Legendre associated functions
int legendre_sphPlm_array (shtns_cfg shtns, const int lmax, const int im, const double x, double *yl)
 Compute values of legendre polynomials noramalized for spherical harmonics, for a range of l=m..lmax, at given m and x, using stable recurrence.
 
int legendre_sphPlm_deriv_array (shtns_cfg shtns, const int lmax, const int im, const double x, const double sint, double *yl, double *dyl)
 
Special operator functions
void mul_ct_matrix (shtns_cfg, double *mx)
 compute the matrix (stored in mx, a double array of size 2*NLM) required to multiply an SH representation by cos(theta) using SH_mul_mx.
 
void st_dt_matrix (shtns_cfg, double *mx)
 compute the matrix (stored in mx, a double array of size 2*NLM) required to apply sin(theta)*d/dtheta to an SH representation using SH_mul_mx.
 
void SH_mul_mx (shtns_cfg, double *mx, cplx *Qlm, cplx *Rlm)
 Apply a matrix involving l+1 and l-1 to an SH representation Qlm. Result stored in Rlm (must be different from Qlm).
 
Scalar transforms
void spat_to_SH (shtns_cfg shtns, double *Vr, cplx *Qlm)
 transform the scalar field Vr into its spherical harmonic representation Qlm.
 
void SH_to_spat (shtns_cfg shtns, cplx *Qlm, double *Vr)
 transform the spherical harmonic coefficients Qlm into its spatial representation Vr.
 
void SH_to_spat_cplx (shtns_cfg shtns, cplx *alm, cplx *z)
 complex scalar synthesis.
 
void spat_cplx_to_SH (shtns_cfg shtns, cplx *z, cplx *alm)
 complex scalar analysis.
 
2D vector transforms
void spat_to_SHsphtor (shtns_cfg, double *Vt, double *Vp, cplx *Slm, cplx *Tlm)
 transform the theta and phi components (Vt,Vp) of a vector into its spheroidal-toroidal spherical harmonic representation (Slm,Tlm).
 
void SHsphtor_to_spat (shtns_cfg, cplx *Slm, cplx *Tlm, double *Vt, double *Vp)
 transform spheroidal-toroidal spherical harmonic coefficients (Slm,Tlm) to the spatial theta and phi components (Vt,Vp).
 
void SHsph_to_spat (shtns_cfg, cplx *Slm, double *Vt, double *Vp)
 transform spheroidal spherical harmonic coefficients Slm to the spatial theta and phi components (Vt,Vp), effectively computing the gradient of S.
 
void SHtor_to_spat (shtns_cfg, cplx *Tlm, double *Vt, double *Vp)
 transform toroidal spherical harmonic coefficients Tlm to the spatial theta and phi components (Vt,Vp).
 
void spat_cplx_to_SHsphtor (shtns_cfg, cplx *Vt, cplx *Vp, cplx *Slm, cplx *Tlm)
 transform the theta and phi components (Vt,Vp) of a complex valued vector into its spheroidal-toroidal spherical harmonic representation (Slm,Tlm).
 
void SHsphtor_to_spat_cplx (shtns_cfg, cplx *Slm, cplx *Tlm, cplx *Vt, cplx *Vp)
 transform spheroidal-toroidal spherical harmonic coefficients (Slm,Tlm) to the spatial theta and phi complex-valued components (Vt,Vp).
 
3D transforms (combine scalar and vector)
void spat_to_SHqst (shtns_cfg, double *Vr, double *Vt, double *Vp, cplx *Qlm, cplx *Slm, cplx *Tlm)
 3D vector transform from spherical coordinates to radial-spheroidal-toroidal spectral components (see Radial - Spheroidal - Toroidal decomposition).
 
void SHqst_to_spat (shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double *Vr, double *Vt, double *Vp)
 3D vector transform from radial-spheroidal-toroidal spectral components (see Radial - Spheroidal - Toroidal decomposition) to spherical coordinates.
 
void spat_cplx_to_SHqst (shtns_cfg, cplx *Vr, cplx *Vt, cplx *Vp, cplx *Qlm, cplx *Slm, cplx *Tlm)
 complex vector transform (3D).
 
void SHqst_to_spat_cplx (shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, cplx *Vr, cplx *Vt, cplx *Vp)
 complex vector transform (3D).
 
Truncated transforms at given degree l

with ltr <= lmax used for setup.

void spat_to_SH_l (shtns_cfg, double *Vr, cplx *Qlm, int ltr)
 
void SH_to_spat_l (shtns_cfg, cplx *Qlm, double *Vr, int ltr)
 
void SHsphtor_to_spat_l (shtns_cfg, cplx *Slm, cplx *Tlm, double *Vt, double *Vp, int ltr)
 
void SHsph_to_spat_l (shtns_cfg, cplx *Slm, double *Vt, double *Vp, int ltr)
 
void SHtor_to_spat_l (shtns_cfg, cplx *Tlm, double *Vt, double *Vp, int ltr)
 
void spat_to_SHsphtor_l (shtns_cfg, double *Vt, double *Vp, cplx *Slm, cplx *Tlm, int ltr)
 
void spat_to_SHqst_l (shtns_cfg, double *Vr, double *Vt, double *Vp, cplx *Qlm, cplx *Slm, cplx *Tlm, int ltr)
 
void SHqst_to_spat_l (shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double *Vr, double *Vt, double *Vp, int ltr)
 
Legendre transform at given m (no fft) and truncated at given degree l <= lmax

The input and output arrays contain only the specified m=im*mres, that is spatial size is nlat and spectral size is lmax+1-m, containing only the coefficients of the given m.

void spat_to_SH_ml (shtns_cfg, int im, cplx *Vr, cplx *Ql, int ltr)
 
void SH_to_spat_ml (shtns_cfg, int im, cplx *Ql, cplx *Vr, int ltr)
 
void spat_to_SHsphtor_ml (shtns_cfg, int im, cplx *Vt, cplx *Vp, cplx *Sl, cplx *Tl, int ltr)
 
void SHsphtor_to_spat_ml (shtns_cfg, int im, cplx *Sl, cplx *Tl, cplx *Vt, cplx *Vp, int ltr)
 
void SHsph_to_spat_ml (shtns_cfg, int im, cplx *Sl, cplx *Vt, cplx *Vp, int ltr)
 
void SHtor_to_spat_ml (shtns_cfg, int im, cplx *Tl, cplx *Vt, cplx *Vp, int ltr)
 
void spat_to_SHqst_ml (shtns_cfg, int im, cplx *Vr, cplx *Vt, cplx *Vp, cplx *Ql, cplx *Sl, cplx *Tl, int ltr)
 
void SHqst_to_spat_ml (shtns_cfg, int im, cplx *Ql, cplx *Sl, cplx *Tl, cplx *Vr, cplx *Vt, cplx *Vp, int ltr)
 
Local and partial evalutions of a SH representation :

Does not require a call to shtns_set_grid_auto

double SH_to_point (shtns_cfg, cplx *Qlm, double cost, double phi)
 Evaluate scalar SH representation Qlm at physical point defined by cost = cos(theta) and phi.
 
void SH_to_point_cplx (shtns_cfg, cplx *alm, double cost, double phi, cplx *z_out)
 writes value to z_out.
 
void SH_to_grad_point (shtns_cfg, cplx *DrSlm, cplx *Slm, double cost, double phi, double *vr, double *vt, double *vp)
 
void SHqst_to_point (shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double cost, double phi, double *vr, double *vt, double *vp)
 Evaluate vector SH representation Qlm at physical point defined by cost = cos(theta) and phi.
 
void SH_to_lat (shtns_cfg shtns, cplx *Qlm, double cost, double *vr, int nphi, int ltr, int mtr)
 synthesis at a given latitude, on nphi equispaced longitude points.
 
void SHqst_to_lat (shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double cost, double *vr, double *vt, double *vp, int nphi, int ltr, int mtr)
 synthesis at a given latitude, on nphi equispaced longitude points.
 

Detailed Description

shtns.h is the definition file for SHTns : include this file in your source code to use SHTns.

Macro Definition Documentation

◆ LM_L_LOOP

#define LM_L_LOOP (   shtns,
  action 
)    { int lm=0; do { int l=shtns->li[lm]; action } while(++lm < shtns->nlm); }

LM_L_LOOP : loop over all (l,im) and perform "action" : l and lm are defined (but NOT m and im).

lm and m must be declared int's. lm is the loop counter and SH array index, while l is the SH degree. more info : Spherical Harmonic coefficients layout

◆ LM_LOOP

#define LM_LOOP (   shtns,
  action 
)    { int lm=0; do { action } while(++lm < shtns->nlm); }

LM_LOOP( shtns, action ) : macro that performs "action" for every (l,m), with lm set, but neither l, m nor im.

lm must be a declared int and is the loop counter and the SH array index. more info : Spherical Harmonic coefficients layout

Examples
SHT_example.c.

◆ NSPAT_ALLOC

#define NSPAT_ALLOC (   shtns)    (shtns->nspat)

total number of 'doubles' required for a spatial field (includes FFTW reserved space).

only the first shtns.nlat*shtns.nphi are real spatial data, the remaining is used by the Fourier Transform. more info : Spatial data layouts and grids used by SHTns

Examples
SHT_example.c.

Enumeration Type Documentation

◆ shtns_norm

enum shtns_norm

Different Spherical Harmonic normalizations. See also section Normalization for details.

Enumerator
sht_orthonormal 

orthonormalized spherical harmonics (default).

sht_fourpi 

Geodesy and spectral analysis : 4.pi normalization.

sht_schmidt 

Schmidt semi-normalized : 4.pi/(2l+1)

◆ shtns_type

enum shtns_type

different SHT types and algorithms

Enumerator
sht_gauss 

use gaussian grid and quadrature. Allows on-the-fly or matrix-based algorithms.

sht_auto 

use fastest algorithm and grid. This currently implies using the Gauss grid and is thus the same as sht_gauss.

sht_reg_fast 

quick initialization of a regular grid. The grid is the same as with sht_reg_dct

sht_reg_dct 

slow initialization of a regular grid (self-tuning). The grid is equispaced and avoids the poles (FĂ©jer quadrature).

sht_quick_init 

gauss grid, with minimum initialization time (useful for pre/post-processing)

sht_reg_poles 

quick initialization of a regular grid including poles (Clenshaw-Curtis quadrature). Useful for vizualisation.

sht_gauss_fly 

legendre polynomials are recomputed on-the-fly for each transform (may be faster on some machines, saves memory and bandwidth).

Function Documentation

◆ legendre_sphPlm_array()

int legendre_sphPlm_array ( shtns_cfg  shtns,
const int  lmax,
const int  im,
const double  x,
double *  yl 
)

Compute values of legendre polynomials noramalized for spherical harmonics, for a range of l=m..lmax, at given m and x, using stable recurrence.

Requires a previous call to shtns_create(). Output compatible with the GSL function gsl_sf_legendre_sphPlm_array(lmax, m, x, yl)

Parameters
[in]lmaxmaximum degree computed
[in]im= m/MRES with m the SH order
[in]xargument, x=cos(theta).
[out]ylis a double array of size (lmax-m+1) filled with the values (of increasing degree l).

◆ sh00_1()

double sh00_1 ( shtns_cfg  shtns)

return the spherical harmonic representation of 1 (l=0,m=0)

return the spherical harmonic representation of 1 (l=0,m=0)

◆ sh10_ct()

double sh10_ct ( shtns_cfg  shtns)

return the spherical harmonic representation of cos(theta) (l=1,m=0)

return the spherical harmonic representation of cos(theta) (l=1,m=0)

◆ sh11_st()

double sh11_st ( shtns_cfg  shtns)

return the spherical harmonic representation of sin(theta)*cos(phi) (l=1,m=1)

return the spherical harmonic representation of sin(theta)*cos(phi) (l=1,m=1)

◆ shlm_e1()

double shlm_e1 ( shtns_cfg  shtns,
int  l,
int  m 
)

return the l,m SH coefficient corresponding to unit energy.

return the l,m SH coefficient corresponding to unit energy.