| 1 | /*
|
|---|
| 2 | * Copyright (c) 2010-2015 Centre National de la Recherche Scientifique.
|
|---|
| 3 | * written by Nathanael Schaeffer (CNRS, ISTerre, Grenoble, France).
|
|---|
| 4 | *
|
|---|
| 5 | * nathanael.schaeffer@ujf-grenoble.fr
|
|---|
| 6 | *
|
|---|
| 7 | * This software is governed by the CeCILL license under French law and
|
|---|
| 8 | * abiding by the rules of distribution of free software. You can use,
|
|---|
| 9 | * modify and/or redistribute the software under the terms of the CeCILL
|
|---|
| 10 | * license as circulated by CEA, CNRS and INRIA at the following URL
|
|---|
| 11 | * "http://www.cecill.info".
|
|---|
| 12 | *
|
|---|
| 13 | * The fact that you are presently reading this means that you have had
|
|---|
| 14 | * knowledge of the CeCILL license and that you accept its terms.
|
|---|
| 15 | *
|
|---|
| 16 | */
|
|---|
| 17 |
|
|---|
| 18 | /** \file shtns.h
|
|---|
| 19 | \brief shtns.h is the definition file for SHTns : include this file in your source code to use SHTns.
|
|---|
| 20 | **/
|
|---|
| 21 |
|
|---|
| 22 | #ifdef _COMPLEX_H
|
|---|
| 23 | /// double precision complex number data type
|
|---|
| 24 | typedef complex double cplx;
|
|---|
| 25 | #else
|
|---|
| 26 | #ifdef __cplusplus
|
|---|
| 27 | #include <complex>
|
|---|
| 28 | typedef std::complex<double> cplx;
|
|---|
| 29 | #else
|
|---|
| 30 | #include <complex.h>
|
|---|
| 31 | typedef complex double cplx;
|
|---|
| 32 | #endif
|
|---|
| 33 | #endif
|
|---|
| 34 |
|
|---|
| 35 |
|
|---|
| 36 | #ifdef __cplusplus
|
|---|
| 37 | extern "C" {
|
|---|
| 38 | #endif /* __cplusplus */
|
|---|
| 39 |
|
|---|
| 40 | /// pointer to data structure describing an SHT, returned by shtns_init() or shtns_create().
|
|---|
| 41 | typedef struct shtns_info* shtns_cfg;
|
|---|
| 42 |
|
|---|
| 43 | /// different Spherical Harmonic normalizations.
|
|---|
| 44 | /// see also section \ref norm for details.
|
|---|
| 45 | enum shtns_norm {
|
|---|
| 46 | sht_orthonormal, ///< orthonormalized spherical harmonics (default).
|
|---|
| 47 | sht_fourpi, ///< Geodesy and spectral analysis : 4.pi normalization.
|
|---|
| 48 | sht_schmidt ///< Schmidt semi-normalized : 4.pi/(2l+1)
|
|---|
| 49 | };
|
|---|
| 50 | #define SHT_NO_CS_PHASE (256*4) ///< don't include Condon-Shortley phase (add to last argument of \ref shtns_create)
|
|---|
| 51 | #define SHT_REAL_NORM (256*8) ///< use a "real" normalization. (add to last argument of \ref shtns_create)
|
|---|
| 52 |
|
|---|
| 53 | /// different SHT types and algorithms
|
|---|
| 54 | enum shtns_type {
|
|---|
| 55 | sht_gauss, ///< use <b>gaussian grid</b> and quadrature.
|
|---|
| 56 | sht_auto, ///< use a regular grid if dct is faster, otherwise defaults to gauss.
|
|---|
| 57 | sht_reg_fast, ///< use fastest algorithm, on a <b>regular grid</b>, mixing dct and regular quadrature.
|
|---|
| 58 | sht_reg_dct, ///< use pure dct algorithm, on a <b>regular grid</b>.
|
|---|
| 59 | sht_quick_init, ///< gauss grid, with minimum initialization time (useful for pre/post-processing)
|
|---|
| 60 | sht_reg_poles, ///< use a <b>synthesis only</b> algo <b>including poles</b>, not suitable for computations. Useful for vizualisation.
|
|---|
| 61 | sht_gauss_fly ///< legendre polynomials are recomputed on-the-fly for each transform (may be faster on some machines, saves memory and bandwidth).
|
|---|
| 62 | };
|
|---|
| 63 | #define SHT_NATIVE_LAYOUT 0 ///< Tells shtns_init to use \ref native
|
|---|
| 64 | #define SHT_THETA_CONTIGUOUS 256 ///< use \ref theta_fast
|
|---|
| 65 | #define SHT_PHI_CONTIGUOUS (256*2) ///< use \ref phi_fast
|
|---|
| 66 | #define SHT_SOUTH_POLE_FIRST (256*32) ///< latitudinal data are stored starting from south pole.
|
|---|
| 67 |
|
|---|
| 68 | #define SHT_SCALAR_ONLY (256*16) ///< don't compute vector matrices. (add to flags in shtns_set_grid)
|
|---|
| 69 | #define SHT_LOAD_SAVE_CFG (256*64) ///< try to load and save the config. (add to flags in shtns_set_grid)
|
|---|
| 70 |
|
|---|
| 71 |
|
|---|
| 72 | #ifndef SHTNS_PRIVATE
|
|---|
| 73 | struct shtns_info { // allow read-only access to some data (useful for optimization and helper macros)
|
|---|
| 74 | const unsigned int nlm; ///< total number of (l,m) spherical harmonics components.
|
|---|
| 75 | const unsigned short lmax; ///< maximum degree (lmax) of spherical harmonics.
|
|---|
| 76 | const unsigned short mmax; ///< maximum order (mmax*mres) of spherical harmonics.
|
|---|
| 77 | const unsigned short mres; ///< the periodicity along the phi axis.
|
|---|
| 78 | const unsigned short nphi; ///< number of spatial points in Phi direction (longitude)
|
|---|
| 79 | const unsigned short nlat; ///< number of spatial points in Theta direction (latitude) ...
|
|---|
| 80 | const unsigned short nlat_2; ///< ...and half of it (using (shtns.nlat+1)/2 allows odd shtns.nlat.)
|
|---|
| 81 | const int *const lmidx; ///< (virtual) index in SH array of given im (size mmax+1) : LiM(l,im) = lmidx[im] + l
|
|---|
| 82 | const unsigned short *const li; ///< degree l for given mode index (size nlm) : li[lm]
|
|---|
| 83 | const unsigned short *const mi; ///< order m for given mode index (size nlm) : mi[lm]
|
|---|
| 84 | const double *const ct; ///< cos(theta) array (size nlat)
|
|---|
| 85 | const double *const st; ///< sin(theta) array (size nlat)
|
|---|
| 86 | const unsigned int nspat; ///< number of real numbers that must be allocated in a spatial field.
|
|---|
| 87 | };
|
|---|
| 88 | #endif
|
|---|
| 89 |
|
|---|
| 90 | // MACROS //
|
|---|
| 91 |
|
|---|
| 92 | /*! \name Access to spherical harmonic components
|
|---|
| 93 | * The following macros give access to single spherical harmonic coefficient or perform loops spanning all of them.
|
|---|
| 94 | **///@{
|
|---|
| 95 | ///LiM(l,im) : macro returning array index for given l and im.
|
|---|
| 96 | #define LiM(shtns, l,im) ( shtns->lmidx[im] + (l) )
|
|---|
| 97 | //#define LiM(shtns, l,im) ( ((im)*(2*shtns->lmax + 2 - ((im)+1)*shtns->mres))>>1 + (l) )
|
|---|
| 98 | /// LM(l,m) : macro returning array index for given l and m.
|
|---|
| 99 | #define LM(shtns, l,m) ( shtns->lmidx[((unsigned)(m))/shtns->mres] + (l) )
|
|---|
| 100 | //#define LM(shtns, l,m) ( ((((unsigned)(m))/shtns->mres)*(2*shtns->lmax + 2 - ((m)+shtns->mres)))>>1 + (l) )
|
|---|
| 101 | /// LM_LOOP( action ) : macro that performs "action" for every (l,m), with lm set, but neither l, m nor im.
|
|---|
| 102 | /// \c lm must be a declared int and is the loop counter and the SH array index. more info : \ref spec_data
|
|---|
| 103 | #define LM_LOOP( shtns, action ) { int lm=0; do { action } while(++lm < shtns->nlm); }
|
|---|
| 104 | /// LM_L_LOOP : loop over all (l,im) and perform "action" : l and lm are defined (but NOT m and im).
|
|---|
| 105 | /// \c lm and \c m must be declared int's. \c lm is the loop counter and SH array index, while \c l is the SH degree. more info : \ref spec_data
|
|---|
| 106 | #define LM_L_LOOP( shtns, action ) { int lm=0; do { int l=shtns->li[lm]; action } while(++lm < shtns->nlm); }
|
|---|
| 107 | //@}
|
|---|
| 108 |
|
|---|
| 109 |
|
|---|
| 110 | /// total number of 'doubles' required for a spatial field (includes FFTW reserved space).
|
|---|
| 111 | /// only the first shtns.nlat*shtns.nphi are real spatial data, the remaining is used by the Fourier Transform. more info : \ref spat
|
|---|
| 112 | #define NSPAT_ALLOC(shtns) (shtns->nspat)
|
|---|
| 113 |
|
|---|
| 114 | // HELPER MACROS //
|
|---|
| 115 |
|
|---|
| 116 | /// phi angle value in degrees for given index ip.
|
|---|
| 117 | #define PHI_DEG(shtns, ip) (360./(shtns->nphi*shtns->mres))*(ip)
|
|---|
| 118 | /// phi angle value in radians for given index ip.
|
|---|
| 119 | #define PHI_RAD(shtns, ip) (M_PI/(shtns->nphi*shtns->mres))*(2*ip)
|
|---|
| 120 |
|
|---|
| 121 |
|
|---|
| 122 | // FUNCTIONS //
|
|---|
| 123 |
|
|---|
| 124 | /// compute number of spherical harmonics modes (l,m) for given size parameters. Does not require any previous setup.
|
|---|
| 125 | long nlm_calc(long lmax, long mmax, long mres);
|
|---|
| 126 |
|
|---|
| 127 | void shtns_verbose(int); ///< controls output during initialization: 0=no output (default), 1=some output, 2=degug (if compiled in)
|
|---|
| 128 | void shtns_print_version(void); ///< print version information to stdout.
|
|---|
| 129 |
|
|---|
| 130 | #ifndef SWIG
|
|---|
| 131 |
|
|---|
| 132 | void shtns_print_cfg(shtns_cfg); ///< print information about given config to stdout.
|
|---|
| 133 |
|
|---|
| 134 | /// \name initialization
|
|---|
| 135 | //@{
|
|---|
| 136 | /// Simple initialization of the spherical harmonic transforms of given size. Calls \ref shtns_create and \ref shtns_set_grid_auto.
|
|---|
| 137 | shtns_cfg shtns_init(enum shtns_type flags, int lmax, int mmax, int mres, int nlat, int nphi);
|
|---|
| 138 | /// Defines the sizes of the spectral description. Use for advanced initialization.
|
|---|
| 139 | shtns_cfg shtns_create(int lmax, int mmax, int mres, enum shtns_norm norm);
|
|---|
| 140 | /// Precompute everything for a given spatial grid. Use for advanced initialization, after \ref shtns_create.
|
|---|
| 141 | int shtns_set_grid(shtns_cfg, enum shtns_type flags, double eps, int nlat, int nphi);
|
|---|
| 142 | /// Precompute everything and choose the optimal nlat and nphi for a given non-linear order.
|
|---|
| 143 | int shtns_set_grid_auto(shtns_cfg, enum shtns_type flags, double eps, int nl_order, int *nlat, int *nphi);
|
|---|
| 144 | /// Copy a given config but allow a different (smaller) mmax and the possibility to enable/disable fft.
|
|---|
| 145 | shtns_cfg shtns_create_with_grid(shtns_cfg, int mmax, int nofft);
|
|---|
| 146 | /// Enables multi-thread transform using OpenMP with num_threads (if available). Returns number of threads that will be used.
|
|---|
| 147 | int shtns_use_threads(int num_threads);
|
|---|
| 148 |
|
|---|
| 149 | void shtns_reset(void); ///< destroy all configs, free memory, and go back to initial state.
|
|---|
| 150 | void shtns_destroy(shtns_cfg); ///< free memory of given config, which cannot be used afterwards.
|
|---|
| 151 | void shtns_unset_grid(shtns_cfg); ///< unset the grid.
|
|---|
| 152 |
|
|---|
| 153 | //@}
|
|---|
| 154 |
|
|---|
| 155 | /// \name special values
|
|---|
| 156 | //@{
|
|---|
| 157 | double sh00_1(shtns_cfg); ///< return the spherical harmonic representation of 1 (l=0,m=0)
|
|---|
| 158 | double sh10_ct(shtns_cfg); ///< return the spherical harmonic representation of cos(theta) (l=1,m=0)
|
|---|
| 159 | double sh11_st(shtns_cfg); ///< return the spherical harmonic representation of sin(theta)*cos(phi) (l=1,m=1)
|
|---|
| 160 | double shlm_e1(shtns_cfg, int l, int m); ///< return the l,m SH coefficient corresponding to unit energy.
|
|---|
| 161 | /// fill the given array with Gauss weights. returns the number of weights written (0 if not a Gauss grid).
|
|---|
| 162 | int shtns_gauss_wts(shtns_cfg, double *wts);
|
|---|
| 163 |
|
|---|
| 164 | //@}
|
|---|
| 165 |
|
|---|
| 166 | /// \name Rotation functions
|
|---|
| 167 | //@{
|
|---|
| 168 | /// Rotate a SH representation Qlm around the z-axis by angle alpha (in radians),
|
|---|
| 169 | /// which is the same as rotating the reference frame by angle -alpha.
|
|---|
| 170 | /// Result is stored in Rlm (which can be the same array as Qlm).
|
|---|
| 171 | void SH_Zrotate(shtns_cfg, cplx *Qlm, double alpha, cplx *Rlm);
|
|---|
| 172 | /// Rotate SH representation around Y axis by alpha (in radians).
|
|---|
| 173 | void SH_Yrotate(shtns_cfg, cplx *Qlm, double alpha, cplx *Rlm);
|
|---|
| 174 | /// Rotate SH representation around Y axis by 90 degrees.
|
|---|
| 175 | void SH_Yrotate90(shtns_cfg, cplx *Qlm, cplx *Rlm);
|
|---|
| 176 | /// Rotate SH representation around X axis by 90 degrees.
|
|---|
| 177 | void SH_Xrotate90(shtns_cfg, cplx *Qlm, cplx *Rlm);
|
|---|
| 178 | //@}
|
|---|
| 179 |
|
|---|
| 180 | /// \name Special operator functions
|
|---|
| 181 | //@{
|
|---|
| 182 | /// compute the matrix (stored in mx, a double array of size 2*NLM) required
|
|---|
| 183 | /// to multiply an SH representation by cos(theta) using \ref SH_mul_mx.
|
|---|
| 184 | void mul_ct_matrix(shtns_cfg, double* mx);
|
|---|
| 185 | /// compute the matrix (stored in mx, a double array of size 2*NLM) required
|
|---|
| 186 | /// to apply sin(theta)*d/dtheta to an SH representation using \ref SH_mul_mx.
|
|---|
| 187 | void st_dt_matrix(shtns_cfg, double* mx);
|
|---|
| 188 | /// Apply a matrix involving l+1 and l-1 to an SH representation Qlm. Result stored in Rlm (must be different from Qlm).
|
|---|
| 189 | void SH_mul_mx(shtns_cfg, double* mx, cplx *Qlm, cplx *Rlm);
|
|---|
| 190 | //@}
|
|---|
| 191 |
|
|---|
| 192 | /** \addtogroup sht Spherical Harmonic transform functions.
|
|---|
| 193 | * All these function perform a global spherical harmonic transform.
|
|---|
| 194 | * Their first argument is a shtns_cfg variable (which is a pointer to a \ref shtns_info struct)
|
|---|
| 195 | * obtained by a previous call to \ref shtns_create or \ref shtns_init.
|
|---|
| 196 | * \see \ref spat \see \ref spec \see \ref vsh
|
|---|
| 197 | */
|
|---|
| 198 | //@{
|
|---|
| 199 |
|
|---|
| 200 | /// \name Scalar transforms
|
|---|
| 201 | //@{
|
|---|
| 202 | /// transform the scalar field Vr into its spherical harmonic representation Qlm.
|
|---|
| 203 | /// \param[in] shtns = a configuration created by \ref shtns_create with a grid set by \ref shtns_set_grid or \ref shtns_set_grid_auto
|
|---|
| 204 | /// \param[in] Vr = spatial scalar field : double array of size shtns->nspat.
|
|---|
| 205 | /// \param[out] Qlm = spherical harmonics coefficients : cplx array of size shtns->nlm.
|
|---|
| 206 | void spat_to_SH(shtns_cfg shtns, double *Vr, cplx *Qlm);
|
|---|
| 207 | /// transform the spherical harmonic coefficients Qlm into its spatial representation Vr.
|
|---|
| 208 | /// \param[in] shtns = a configuration created by \ref shtns_create with a grid set by \ref shtns_set_grid or \ref shtns_set_grid_auto
|
|---|
| 209 | /// \param[in] Qlm = spherical harmonics coefficients : cplx array of size shtns->nlm.
|
|---|
| 210 | /// \param[out] Vr = spatial scalar field : double array of size shtns->nspat.
|
|---|
| 211 | void SH_to_spat(shtns_cfg shtns, cplx *Qlm, double *Vr);
|
|---|
| 212 | /// complex scalar synthesis.
|
|---|
| 213 | /// \param[in] shtns = a configuration created by \ref shtns_create with a grid set by \ref shtns_set_grid or \ref shtns_set_grid_auto
|
|---|
| 214 | /// \param[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]
|
|---|
| 215 | /// \param[out] z = complex spatial field
|
|---|
| 216 | void SH_to_spat_cplx(shtns_cfg shtns, cplx *alm, cplx *z);
|
|---|
| 217 | /// complex scalar analysis.
|
|---|
| 218 | /// \param[in] shtns = a configuration created by \ref shtns_create with a grid set by \ref shtns_set_grid or \ref shtns_set_grid_auto
|
|---|
| 219 | /// \param[in] z = complex spatial field
|
|---|
| 220 | /// \param[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]
|
|---|
| 221 | void spat_cplx_to_SH(shtns_cfg shtns, cplx *z, cplx *alm);
|
|---|
| 222 | //@}
|
|---|
| 223 |
|
|---|
| 224 | /// \name 2D vector transforms
|
|---|
| 225 | //@{
|
|---|
| 226 | /// transform the theta and phi components (Vt,Vp) of a vector into its spheroidal-toroidal spherical harmonic representation (Slm,Tlm). \see \ref vsh
|
|---|
| 227 | void spat_to_SHsphtor(shtns_cfg, double *Vt, double *Vp, cplx *Slm, cplx *Tlm);
|
|---|
| 228 | /// transform spheroidal-toroidal spherical harmonic coefficients (Slm,Tlm) to the spatial theta and phi components (Vt,Vp). \see \ref vsh
|
|---|
| 229 | void SHsphtor_to_spat(shtns_cfg, cplx *Slm, cplx *Tlm, double *Vt, double *Vp);
|
|---|
| 230 | /// transform spheroidal spherical harmonic coefficients Slm to the spatial theta and phi components (Vt,Vp), effectively computing the gradient of S. \see \ref vsh
|
|---|
| 231 | void SHsph_to_spat(shtns_cfg, cplx *Slm, double *Vt, double *Vp);
|
|---|
| 232 | /// transform toroidal spherical harmonic coefficients Tlm to the spatial theta and phi components (Vt,Vp). \see \ref vsh
|
|---|
| 233 | void SHtor_to_spat(shtns_cfg, cplx *Tlm, double *Vt, double *Vp);
|
|---|
| 234 | //@}
|
|---|
| 235 | /// Compute the spatial representation of the gradient of a scalar SH field. Alias for \ref SHsph_to_spat
|
|---|
| 236 | #define SH_to_grad_spat(shtns, S,Gt,Gp) SHsph_to_spat(shtns, S, Gt, Gp)
|
|---|
| 237 |
|
|---|
| 238 | /// \name 3D transforms (combine scalar and vector)
|
|---|
| 239 | //@{
|
|---|
| 240 | /// 3D vector transform from spherical coordinates to radial-spheroidal-toroidal spectral components (see \ref vsh_def).
|
|---|
| 241 | /// They should be prefered over separate calls to scalar and 2D vector transforms as they are often significantly faster.
|
|---|
| 242 | void spat_to_SHqst(shtns_cfg, double *Vr, double *Vt, double *Vp, cplx *Qlm, cplx *Slm, cplx *Tlm);
|
|---|
| 243 | /// 3D vector transform from spherical coordinates to radial-spheroidal-toroidal spectral components (see \ref vsh_def).
|
|---|
| 244 | /// They should be prefered over separate calls to scalar and 2D vector transforms as they are often significantly faster.
|
|---|
| 245 | void SHqst_to_spat(shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double *Vr, double *Vt, double *Vp);
|
|---|
| 246 | //@}
|
|---|
| 247 |
|
|---|
| 248 | /// \name Truncated transforms at given degree l
|
|---|
| 249 | /// wiht l <= lmax used for setup.
|
|---|
| 250 | //@{
|
|---|
| 251 | void spat_to_SH_l(shtns_cfg, double *Vr, cplx *Qlm, int ltr);
|
|---|
| 252 | void SH_to_spat_l(shtns_cfg, cplx *Qlm, double *Vr, int ltr);
|
|---|
| 253 |
|
|---|
| 254 | void SHsphtor_to_spat_l(shtns_cfg, cplx *Slm, cplx *Tlm, double *Vt, double *Vp, int ltr);
|
|---|
| 255 | void SHsph_to_spat_l(shtns_cfg, cplx *Slm, double *Vt, double *Vp, int ltr);
|
|---|
| 256 | void SHtor_to_spat_l(shtns_cfg, cplx *Tlm, double *Vt, double *Vp, int ltr);
|
|---|
| 257 | void spat_to_SHsphtor_l(shtns_cfg, double *Vt, double *Vp, cplx *Slm, cplx *Tlm, int ltr);
|
|---|
| 258 |
|
|---|
| 259 | void spat_to_SHqst_l(shtns_cfg, double *Vr, double *Vt, double *Vp, cplx *Qlm, cplx *Slm, cplx *Tlm, int ltr);
|
|---|
| 260 | void SHqst_to_spat_l(shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double *Vr, double *Vt, double *Vp, int ltr);
|
|---|
| 261 | //@}
|
|---|
| 262 | /// Compute the spatial representation of the gradient of a scalar SH field. Alias for \ref SHsph_to_spat_l
|
|---|
| 263 | #define SH_to_grad_spat_l(shtns, S,Gt,Gp,ltr) SHsph_to_spat_l(shtns, S, Gt, Gp, ltr)
|
|---|
| 264 |
|
|---|
| 265 | /// \name Legendre transform at given m (no fft) and truncated at given degree l <= lmax
|
|---|
| 266 | //@{
|
|---|
| 267 | void spat_to_SH_ml(shtns_cfg, int im, cplx *Vr, cplx *Ql, int ltr);
|
|---|
| 268 | void SH_to_spat_ml(shtns_cfg, int im, cplx *Ql, cplx *Vr, int ltr);
|
|---|
| 269 |
|
|---|
| 270 | void spat_to_SHsphtor_ml(shtns_cfg, int im, cplx *Vt, cplx *Vp, cplx *Sl, cplx *Tl, int ltr);
|
|---|
| 271 | void SHsphtor_to_spat_ml(shtns_cfg, int im, cplx *Sl, cplx *Tl, cplx *Vt, cplx *Vp, int ltr);
|
|---|
| 272 | void SHsph_to_spat_ml(shtns_cfg, int im, cplx *Sl, cplx *Vt, cplx *Vp, int ltr);
|
|---|
| 273 | void SHtor_to_spat_ml(shtns_cfg, int im, cplx *Tl, cplx *Vt, cplx *Vp, int ltr);
|
|---|
| 274 |
|
|---|
| 275 | void spat_to_SHqst_ml(shtns_cfg, int im, cplx *Vr, cplx *Vt, cplx *Vp, cplx *Ql, cplx *Sl, cplx *Tl, int ltr);
|
|---|
| 276 | void SHqst_to_spat_ml(shtns_cfg, int im, cplx *Ql, cplx *Sl, cplx *Tl, cplx *Vr, cplx *Vt, cplx *Vp, int ltr);
|
|---|
| 277 | //@}
|
|---|
| 278 | /// Compute the spatial representation of the gradient of a scalar SH field. Alias for \ref SHsph_to_spat_l
|
|---|
| 279 | #define SH_to_grad_spat_ml(shtns, im, S,Gt,Gp,ltr) SHsph_to_spat_ml(shtns, im, S, Gt, Gp, ltr)
|
|---|
| 280 |
|
|---|
| 281 | //@}
|
|---|
| 282 |
|
|---|
| 283 | /// \name Local and partial evalutions of a SH representation :
|
|---|
| 284 | /// Does not require a call to \ref shtns_set_grid_auto
|
|---|
| 285 | //@{
|
|---|
| 286 | double SH_to_point(shtns_cfg, cplx *Qlm, double cost, double phi);
|
|---|
| 287 | void SH_to_grad_point(shtns_cfg, cplx *DrSlm, cplx *Slm,
|
|---|
| 288 | double cost, double phi, double *vr, double *vt, double *vp);
|
|---|
| 289 | void SHqst_to_point(shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm,
|
|---|
| 290 | double cost, double phi, double *vr, double *vt, double *vp);
|
|---|
| 291 |
|
|---|
| 292 | void SH_to_lat(shtns_cfg shtns, cplx *Qlm, double cost,
|
|---|
| 293 | double *vr, int nphi, int ltr, int mtr);
|
|---|
| 294 | void SHqst_to_lat(shtns_cfg, cplx *Qlm, cplx *Slm, cplx *Tlm, double cost,
|
|---|
| 295 | double *vr, double *vt, double *vp, int nphi, int ltr, int mtr);
|
|---|
| 296 | //@}
|
|---|
| 297 |
|
|---|
| 298 | #endif
|
|---|
| 299 |
|
|---|
| 300 | #ifdef __cplusplus
|
|---|
| 301 | }
|
|---|
| 302 | #endif /* __cplusplus */
|
|---|