SHTns 3.6.1
Macros
Spherical Harmonic transform functions.

All these function perform a global spherical harmonic transform. More...

Macros

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

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)
 

Detailed Description

All these function perform a global spherical harmonic transform.

Their first argument is a shtns_cfg variable (which is a pointer to a shtns_info struct) obtained by a previous call to shtns_create or shtns_init.

See also
Spatial data layouts and grids used by SHTns
Spherical Harmonics storage and normalization
Vector Spherical Harmonics as implemented in SHTns

Macro Definition Documentation

◆ SH_to_grad_spat

#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

◆ SH_to_grad_spat_l

#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

◆ SH_to_grad_spat_ml

#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

Function Documentation

◆ SH_to_spat()

void SH_to_spat ( shtns_cfg  shtns,
cplx Qlm,
double *  Vr 
)

transform the spherical harmonic coefficients Qlm into its spatial representation Vr.

Parameters
[in]shtns= a configuration created by shtns_create with a grid set by shtns_set_grid or shtns_set_grid_auto
[in]Qlm= spherical harmonics coefficients : cplx array of size shtns->nlm.
[out]Vr= spatial scalar field : double array of size shtns->nspat.

transform the spherical harmonic coefficients Qlm into its spatial representation Vr.

Examples
SHT_example.c.

◆ SH_to_spat_cplx()

void SH_to_spat_cplx ( shtns_cfg  shtns,
cplx alm,
cplx z 
)

complex scalar synthesis.

Parameters
[in]shtns= a configuration created by shtns_create with a grid set by shtns_set_grid or shtns_set_grid_auto
[in]alm[l*(l+1)+m]is the SH coefficient of order l and degree m (with -l <= m <= l) [total of (LMAX+1)^2 coefficients]
[out]z= complex spatial field

complex scalar synthesis.

in: alm[LM_cplx(shtns,l,m)] is the SH coefficients of order l and degree m (with -l <= m <= l) for a total of nlm_cplx_calc(lmax,mmax,mres) coefficients. out: complex spatial field z.

◆ SHqst_to_spat()

void SHqst_to_spat ( shtns_cfg  shtns,
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.

They should be prefered over separate calls to scalar and 2D vector transforms as they can be significantly faster.

◆ SHqst_to_spat_cplx()

void SHqst_to_spat_cplx ( shtns_cfg  shtns,
cplx qlm,
cplx slm,
cplx tlm,
cplx zr,
cplx zt,
cplx zp 
)

complex vector transform (3D).

in: {qlm,slm,tlm}[LM_cplx(l,m)] are the SH coefficients of order l and degree m (with -l <= m <= l) out: zr,zt,zp: r,theta,phi components of the complex spatial vector field.

◆ SHsph_to_spat()

void SHsph_to_spat ( shtns_cfg  shtns,
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.

See also
Vector Spherical Harmonics as implemented in SHTns

◆ SHsphtor_to_spat()

void SHsphtor_to_spat ( shtns_cfg  shtns,
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).

See also
Vector Spherical Harmonics as implemented in SHTns
Examples
SHT_example.c.

◆ SHsphtor_to_spat_cplx()

void SHsphtor_to_spat_cplx ( shtns_cfg  shtns,
cplx slm,
cplx tlm,
cplx zt,
cplx zp 
)

transform spheroidal-toroidal spherical harmonic coefficients (Slm,Tlm) to the spatial theta and phi complex-valued components (Vt,Vp).

See also
Vector Spherical Harmonics as implemented in SHTns

transform spheroidal-toroidal spherical harmonic coefficients (Slm,Tlm) to the spatial theta and phi complex-valued components (Vt,Vp).

in: slm, tlm are the spheroidal/toroidal SH coefficients of order l and degree m (with -l <= m <= l) out: zt, zp are respectively the theta and phi components of the complex spatial vector field.

◆ SHtor_to_spat()

void SHtor_to_spat ( shtns_cfg  shtns,
cplx Tlm,
double *  Vt,
double *  Vp 
)

transform toroidal spherical harmonic coefficients Tlm to the spatial theta and phi components (Vt,Vp).

See also
Vector Spherical Harmonics as implemented in SHTns
Examples
SHT_example.c.

◆ spat_cplx_to_SH()

void spat_cplx_to_SH ( shtns_cfg  shtns,
cplx z,
cplx alm 
)

complex scalar analysis.

Parameters
[in]shtns= a configuration created by shtns_create with a grid set by shtns_set_grid or shtns_set_grid_auto
[in]z= complex spatial field
[out]alm[l*(l+1)+m]is the SH coefficient of order l and degree m (with -l <= m <= l) [total of (LMAX+1)^2 coefficients]

complex scalar analysis.

in: complex spatial field z. out: alm[LM(shtns,l,m)] is the SH coefficients of order l and degree m (with -l <= m <= l) for a total of nlm_cplx_calc(lmax,mmax,mres) coefficients.

◆ spat_cplx_to_SHqst()

void spat_cplx_to_SHqst ( shtns_cfg  shtns,
cplx zr,
cplx zt,
cplx zp,
cplx qlm,
cplx slm,
cplx tlm 
)

complex vector transform (3D).

in: zr,zt,zp are the r,theta,phi components of the complex spatial vector field. out: {qlm,slm,tlm}[LM_cplx(l,m)] are the SH coefficients of order l and degree m (with -l <= m <= l)

◆ spat_cplx_to_SHsphtor()

void spat_cplx_to_SHsphtor ( shtns_cfg  shtns,
cplx zt,
cplx zp,
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).

See also
Vector Spherical Harmonics as implemented in SHTns

transform the theta and phi components (Vt,Vp) of a complex valued vector into its spheroidal-toroidal spherical harmonic representation (Slm,Tlm).

zt,zp: theta,phi components of the complex spatial vector field. out: slm[LM_cplx(l,m)] and tlm[LM_cplx(l,m)] are the SH coefficients of order l and degree m (with -l <= m <= l) for a total of shtns->nlm_cplx = nlm_cplx_calc(lmax, mmax, mres) coefficients.

◆ spat_to_SH()

void spat_to_SH ( shtns_cfg  shtns,
double *  Vr,
cplx Qlm 
)

transform the scalar field Vr into its spherical harmonic representation Qlm.

Parameters
[in]shtns= a configuration created by shtns_create with a grid set by shtns_set_grid or shtns_set_grid_auto
[in]Vr= spatial scalar field : double array of size shtns->nspat; NOT GUARANTEED TO BE PRESERVED.
[out]Qlm= spherical harmonics coefficients : cplx array of size shtns->nlm.
Examples
SHT_example.c.

◆ spat_to_SHqst()

void spat_to_SHqst ( shtns_cfg  shtns,
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).

They should be prefered over separate calls to scalar and 2D vector transforms as they can be significantly faster.

◆ spat_to_SHsphtor()

void spat_to_SHsphtor ( shtns_cfg  shtns,
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).

See also
Vector Spherical Harmonics as implemented in SHTns