SHTns 3.6.1
Macros | Functions | Variables
Fortran API.

Call from fortran without the trailing '_'. More...

Macros

#define MAX_CONFIGS   (4)
 

Functions

void shtns_verbose_ (int *v)
 Set verbosity level.
 
void shtns_use_threads_ (int *num_threads)
 Enable threads.
 
void shtns_print_cfg_ ()
 Print info.
 
void shtns_init_sh_gauss_ (int *layout, int *lmax, int *mmax, int *mres, int *nlat, int *nphi)
 Initializes spherical harmonic transforms of given size using Gauss algorithm with default polar optimization.
 
void shtns_init_sh_auto_ (int *layout, int *lmax, int *mmax, int *mres, int *nlat, int *nphi)
 Initializes spherical harmonic transforms of given size using Fastest available algorithm and polar optimization.
 
void shtns_init_sh_reg_fast_ (int *layout, int *lmax, int *mmax, int *mres, int *nlat, int *nphi)
 Initializes spherical harmonic transforms of given size using a regular grid and agressive optimizations.
 
void shtns_init_sh_poles_ (int *layout, int *lmax, int *mmax, int *mres, int *nlat, int *nphi)
 Initializes spherical harmonic transform SYNTHESIS ONLY of given size using a regular grid including poles.
 
void shtns_set_size_ (int *lmax, int *mmax, int *mres, int *norm)
 Defines the size and convention of the transform.
 
void shtns_precompute_ (int *type, int *layout, double *eps, int *nlat, int *nphi)
 Precompute matrices for synthesis and analysis.
 
void shtns_precompute_auto_ (int *type, int *layout, double *eps, int *nl_order, int *nlat, int *nphi)
 Same as shtns_precompute_ but choose optimal nlat and/or nphi.
 
void shtns_reset_ ()
 Clear everything.
 
void shtns_save_cfg_ (unsigned *n)
 Saves the current shtns configuration with tag n.
 
void shtns_load_cfg_ (unsigned *n)
 Loads the previously saved configuration tagged n.
 
void shtns_calc_nlm_ (int *nlm, const int *const lmax, const int *const mmax, const int *const mres)
 returns nlm, the number of complex*16 elements in an SH array.
 
int shtns_lmidx_fortran (shtns_cfg shtns, const int *const l, const int *const m)
 
void shtns_lmidx_ (int *lm, const int *const l, const int *const m)
 returns lm, the index in an SH array of mode (l,m).
 
int shtns_lm2l_fortran (shtns_cfg shtns, const int *lm)
 
int shtns_lm2m_fortran (shtns_cfg shtns, const int *lm)
 
void shtns_l_m_ (int *l, int *m, const int *const lm)
 returns l and m, degree and order of an index in SH array lm.
 
void shtns_cos_array_ (double *costh)
 fills the given array with the cosine of the co-latitude angle (NLAT real*8) if no grid has been set, the first element will be set to zero.
 
void shtns_gauss_wts_ (double *wts)
 fills the given array with the gaussian quadrature weights ((NLAT+1)/2 real*8).
 
void shtns_sh_zrotate_ (cplx *Qlm, double *alpha, cplx *Rlm)
 
void shtns_sh_yrotate_ (cplx *Qlm, double *alpha, cplx *Rlm)
 
void shtns_sh_xrotate90_ (cplx *Qlm, cplx *Rlm)
 
void shtns_sh_yrotate90_ (cplx *Qlm, cplx *Rlm)
 
void shtns_sh_cplx_zrotate_ (cplx *Qlm, double *alpha, cplx *Rlm)
 
void shtns_sh_cplx_yrotate_ (cplx *Qlm, double *alpha, cplx *Rlm)
 
void shtns_sh_cplx_xrotate90_ (cplx *Qlm, cplx *Rlm)
 
void shtns_sh_cplx_yrotate90_ (cplx *Qlm, cplx *Rlm)
 

Variables

shtns_cfg shtns_configs [(4)]
 

Point evaluation of Spherical Harmonics

Evaluate at a given point ( $cos(\theta)$ and $\phi$) a spherical harmonic representation.

void shtns_sh_to_point_ (double *spat, cplx *Qlm, double *cost, double *phi)
 
void shtns_qst_to_point_ (double *vr, double *vt, double *vp, cplx *Qlm, cplx *Slm, cplx *Tlm, double *cost, double *phi)
 

Detailed Description

Call from fortran without the trailing '_'.

see the Fortran example for a simple usage of SHTns from Fortran language.

Function Documentation

◆ shtns_calc_nlm_()

void shtns_calc_nlm_ ( int *  nlm,
const int *const  lmax,
const int *const  mmax,
const int *const  mres 
)

returns nlm, the number of complex*16 elements in an SH array.

call from fortran using

call shtns_calc_nlm(nlm, lmax, mmax, mres)

◆ shtns_gauss_wts_()

void shtns_gauss_wts_ ( double *  wts)

fills the given array with the gaussian quadrature weights ((NLAT+1)/2 real*8).

when there is no gaussian grid, the first element is set to zero.

◆ shtns_l_m_()

void shtns_l_m_ ( int *  l,
int *  m,
const int *const  lm 
)

returns l and m, degree and order of an index in SH array lm.

call from fortran using

call shtns_l_m(l, m, lm)

◆ shtns_lmidx_()

void shtns_lmidx_ ( int *  lm,
const int *const  l,
const int *const  m 
)

returns lm, the index in an SH array of mode (l,m).

call from fortran using

call shtns_lmidx(lm, l, m)

◆ shtns_precompute_()

void shtns_precompute_ ( int *  type,
int *  layout,
double *  eps,
int *  nlat,
int *  nphi 
)

Precompute matrices for synthesis and analysis.

Allow to choose polar optimization threshold and algorithm type.

See also
shtns_set_grid

◆ shtns_precompute_auto_()

void shtns_precompute_auto_ ( int *  type,
int *  layout,
double *  eps,
int *  nl_order,
int *  nlat,
int *  nphi 
)

Same as shtns_precompute_ but choose optimal nlat and/or nphi.

See also
shtns_set_grid_auto

◆ shtns_qst_to_point_()

void shtns_qst_to_point_ ( double *  vr,
double *  vt,
double *  vp,
cplx Qlm,
cplx Slm,
cplx Tlm,
double *  cost,
double *  phi 
)
See also
SHqst_to_point for argument description

◆ shtns_set_size_()

void shtns_set_size_ ( int *  lmax,
int *  mmax,
int *  mres,
int *  norm 
)

Defines the size and convention of the transform.

Allow to choose the normalization and whether or not to include the Condon-Shortley phase.

See also
shtns_create

◆ shtns_sh_to_point_()

void shtns_sh_to_point_ ( double *  spat,
cplx Qlm,
double *  cost,
double *  phi 
)
See also
SH_to_point for argument description