SHTns  2.6.5
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 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 NSPAT_ALLOC(shtns)   (shtns->nspat)
 total number of 'doubles' required for a spatial field (includes FFTW reserved space). More...
 
#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)   (M_PI/(shtns->nphi*shtns->mres))*(2*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. Alias for SHsph_to_spat.
 
#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. Alias for SHsph_to_spat_l.
 
#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. Alias for SHsph_to_spat_l.
 
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)   ( shtns->lmidx[im] + (l) )
 LiM(l,im) : macro returning array index for given l and im.
 
#define LM(shtns, l, m)   ( shtns->lmidx[((unsigned)(m))/shtns->mres] + (l) )
 LM(l,m) : macro returning array index for given l and m.
 
#define LM_LOOP(shtns, action)   { int lm=0; do { action } while(++lm < shtns->nlm); }
 LM_LOOP( action ) : macro that performs "action" for every (l,m), with lm set, but neither l, m nor im. More...
 
#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). More...
 

Typedefs

typedef complex double cplx
 
typedef struct shtns_infoshtns_cfg
 pointer to data structure describing an SHT, returned by shtns_init() or shtns_create().
 

Enumerations

enum  shtns_norm { sht_orthonormal, sht_fourpi, sht_schmidt }
 different Spherical Harmonic normalizations. 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) for given size parameters. Does not require any previous setup. More...
 
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.
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
void shtns_reset (void)
 destroy all configs, free memory, and go back to initial state. More...
 
void shtns_destroy (shtns_cfg)
 free memory of given config, which cannot be used afterwards. More...
 
void shtns_unset_grid (shtns_cfg)
 unset the grid. More...
 
special values
double sh00_1 (shtns_cfg)
 return the spherical harmonic representation of 1 (l=0,m=0) More...
 
double sh10_ct (shtns_cfg)
 return the spherical harmonic representation of cos(theta) (l=1,m=0) More...
 
double sh11_st (shtns_cfg)
 return the spherical harmonic representation of sin(theta)*cos(phi) (l=1,m=1) More...
 
double shlm_e1 (shtns_cfg, int l, int m)
 return the l,m SH coefficient corresponding to unit energy. More...
 
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). More...
 
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. More...
 
void SH_Yrotate (shtns_cfg, cplx *Qlm, double alpha, cplx *Rlm)
 Rotate SH representation around Y axis by alpha (in radians). More...
 
void SH_Yrotate90 (shtns_cfg, cplx *Qlm, cplx *Rlm)
 Rotate SH representation around Y axis by 90 degrees. More...
 
void SH_Xrotate90 (shtns_cfg, cplx *Qlm, cplx *Rlm)
 Rotate SH representation around X axis by 90 degrees. More...
 
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. More...
 
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. More...
 
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). More...
 
Scalar transforms
void spat_to_SH (shtns_cfg shtns, double *Vr, cplx *Qlm)
 transform the scalar field Vr into its spherical harmonic representation Qlm. More...
 
void SH_to_spat (shtns_cfg shtns, cplx *Qlm, double *Vr)
 transform the spherical harmonic coefficients Qlm into its spatial representation Vr. More...
 
void SH_to_spat_cplx (shtns_cfg shtns, cplx *alm, cplx *z)
 complex scalar synthesis. More...
 
void spat_cplx_to_SH (shtns_cfg shtns, cplx *z, cplx *alm)
 complex scalar analysis. More...
 
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). More...
 
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). More...
 
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. More...
 
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). More...
 
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). More...
 
void SHqst_to_spat (shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double *Vr, double *Vt, double *Vp)
 3D vector transform from spherical coordinates to radial-spheroidal-toroidal spectral components (see Radial - Spheroidal - Toroidal decomposition). More...
 
Truncated transforms at given degree l

wiht l <= 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
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_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. More...
 
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. More...
 

Detailed Description

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

Macro Definition Documentation

#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

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

LM_LOOP( 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.
#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

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)

enum shtns_type

different SHT types and algorithms

Enumerator
sht_gauss 

use gaussian grid and quadrature.

sht_auto 

use a regular grid if dct is faster, otherwise defaults to gauss.

sht_reg_fast 

use fastest algorithm, on a regular grid, mixing dct and regular quadrature.

sht_reg_dct 

use pure dct algorithm, on a regular grid.

sht_quick_init 

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

sht_reg_poles 

use a synthesis only algo including poles, not suitable for computations. 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

long nlm_calc ( long  lmax,
long  mmax,
long  mres 
)

compute number of spherical harmonics modes (l,m) for given size parameters. Does not require any previous setup.

return (mmax+1)*(lmax+1) - mres*(mmax*(mmax+1))/2;
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)

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)

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)

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.