SHTns  2.6.5
Functions
Fortran API.

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

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. More...
 
void shtns_precompute_ (int *type, int *layout, double *eps, int *nlat, int *nphi)
 Precompute matrices for synthesis and analysis. More...
 
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. More...
 
void shtns_reset_ ()
 Clear everything.
 
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. More...
 
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). More...
 
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. More...
 
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. More...
 
void shtns_gauss_wts_ (double *wts)
 fills the given array with the gaussian quadrature weights ((NLAT+1)/2 real*8). More...
 

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

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

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

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)
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)
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
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
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
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
void shtns_sh_to_point_ ( double *  spat,
cplx *  Qlm,
double *  cost,
double *  phi 
)
See Also
SH_to_point for argument description