All these function perform a global spherical harmonic transform.
More...
|
#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.
|
|
|
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...
|
|
|
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...
|
|
|
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) |
|
|
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) |
|
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
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. |
- Examples:
- SHT_example.c.
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[l*(l+1)+m] is the SH coefficients of order l and degree m (with -l <= m <= l) for a total of (LMAX+1)^2 coefficients. out: complex spatial field.
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).
They should be prefered over separate calls to scalar and 2D vector transforms as they are often significantly faster.
void SHsph_to_spat |
( |
shtns_cfg |
, |
|
|
cplx * |
Slm, |
|
|
double * |
Vt, |
|
|
double * |
Vp |
|
) |
| |
void SHsphtor_to_spat |
( |
shtns_cfg |
, |
|
|
cplx * |
Slm, |
|
|
cplx * |
Tlm, |
|
|
double * |
Vt, |
|
|
double * |
Vp |
|
) |
| |
void SHtor_to_spat |
( |
shtns_cfg |
, |
|
|
cplx * |
Tlm, |
|
|
double * |
Vt, |
|
|
double * |
Vp |
|
) |
| |
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. out: alm[l*(l+1)+m] is the SH coefficients of order l and degree m (with -l <= m <= l) for a total of (LMAX+1)^2 coefficients.
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. |
[out] | Qlm | = spherical harmonics coefficients : cplx array of size shtns->nlm. |
- Examples:
- SHT_example.c.
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).
They should be prefered over separate calls to scalar and 2D vector transforms as they are often significantly faster.
void spat_to_SHsphtor |
( |
shtns_cfg |
, |
|
|
double * |
Vt, |
|
|
double * |
Vp, |
|
|
cplx * |
Slm, |
|
|
cplx * |
Tlm |
|
) |
| |