|
Ssht
1.3.3
Fast and exact spin spherical harmonic transforms
|
#include "ssht/ssht_types.h"
Macros | |
| #define | ABS(x) ((x) >= 0 ? (x) : -(x)) |
Functions | |
| void | gauleg (double x1, double x2, double *x, double *w, int n) |
| ssht_complex_double | ssht_sampling_weight_mw (int p) |
| double | ssht_sampling_weight_dh (double theta_t, int L) |
| void | ssht_sampling_gl_thetas_weights (double *thetas, double *weights, int L) |
| int | ssht_sampling_mw_ntheta (int L) |
| double | ssht_sampling_mw_t2theta (int t, int L) |
| int | ssht_sampling_mw_nphi (int L) |
| double | ssht_sampling_mw_p2phi (int p, int L) |
| int | ssht_sampling_mw_n (int L) |
| double | ssht_sampling_mw_ss_t2theta (int t, int L) |
| int | ssht_sampling_mw_ss_ntheta (int L) |
| double | ssht_sampling_mw_ss_p2phi (int p, int L) |
| int | ssht_sampling_mw_ss_nphi (int L) |
| int | ssht_sampling_mw_ss_n (int L) |
| double | ssht_sampling_dh_t2theta (int t, int L) |
| int | ssht_sampling_dh_ntheta (int L) |
| double | ssht_sampling_dh_p2phi (int p, int L) |
| int | ssht_sampling_dh_nphi (int L) |
| int | ssht_sampling_dh_n (int L) |
| int | ssht_sampling_gl_ntheta (int L) |
| double | ssht_sampling_gl_p2phi (int p, int L) |
| int | ssht_sampling_gl_nphi (int L) |
| int | ssht_sampling_gl_n (int L) |
Functionality to define sample positions for various algorithms, to compute weights and to convert 1D and 2D harmonic indices.
| #define ABS | ( | x | ) | ((x) >= 0 ? (x) : -(x)) |
| void gauleg | ( | double | x1, |
| double | x2, | ||
| double * | x, | ||
| double * | w, | ||
| int | n | ||
| ) |
Given the lower and upper limits of integration x1 and x2, this routine returns arrays x[1..n] and w[1..n] of length n, containing the abscissas and weights of the Gauss-Legendre n-point quadrature formula.
| [in] | x1 | Lower bound of range. |
| [in] | x2 | Upper bound of range. |
| [out] | x | Node positions (i.e. roots of Legendre polynomials). |
| [out] | w | Corresponding weights. |
| [in] | n | Number of points (note memory must already be allocated for x and w to store n terms). |
| int ssht_sampling_dh_n | ( | int | L | ) |
Compute total number of samples for Driscoll and Healy sampling.
| [in] | L | Harmonic band-limit. |
| n | Number of samples. |
| int ssht_sampling_dh_nphi | ( | int | L | ) |
Compute number of phi samples for Driscoll and Healy sampling.
| [in] | L | Harmonic band-limit. |
| nphi | Number of phi samples. |
| int ssht_sampling_dh_ntheta | ( | int | L | ) |
Compute number of theta samples for Driscoll and Healy sampling.
| [in] | L | Harmonic band-limit. |
| ntheta | Number of theta samples. |
| double ssht_sampling_dh_p2phi | ( | int | p, |
| int | L | ||
| ) |
Convert phi index to angle for Driscoll and Healy sampling.
| [in] | p | Phi index. |
| [in] | L | Harmonic band-limit. |
| phi | Phi angle. |
| double ssht_sampling_dh_t2theta | ( | int | t, |
| int | L | ||
| ) |
Convert theta index to angle for Driscoll and Healy sampling.
| [in] | t | Theta index. |
| [in] | L | Harmonic band-limit. |
| theta | Theta angle. |
| int ssht_sampling_gl_n | ( | int | L | ) |
Compute total number of samples for Gauss-Legendre sampling.
| [in] | L | Harmonic band-limit. |
| n | Number of samples. |
| int ssht_sampling_gl_nphi | ( | int | L | ) |
Compute number of phi samples for Gauss-Legendre sampling.
| [in] | L | Harmonic band-limit. |
| nphi | Number of phi samples. |
| int ssht_sampling_gl_ntheta | ( | int | L | ) |
Compute number of theta samples for Gauss-Legendre sampling.
| [in] | L | Harmonic band-limit. |
| ntheta | Number of theta samples. |
| double ssht_sampling_gl_p2phi | ( | int | p, |
| int | L | ||
| ) |
Convert phi index to angle for Gauss-Legendre sampling.
| [in] | p | Phi index. |
| [in] | L | Harmonic band-limit. |
| phi | Phi angle. |
| void ssht_sampling_gl_thetas_weights | ( | double * | thetas, |
| double * | weights, | ||
| int | L | ||
| ) |
Compute Gauss-Legendre theta positions (roots of Legendre polynomials) and corresponding weights.
| [out] | thetas | L theta positions (memory must already be allocated to store the L theta positions). |
| [out] | weights | Corresponding weights (memory must already be allocated to store the L weights corresponding to each theta position). |
| [in] | L | Harmonic band-limit. |
| none |
| int ssht_sampling_mw_n | ( | int | L | ) |
Compute total number of samples for McEwen and Wiaux sampling.
/note Computes number of samples on sphere, not over extended domain.
| [in] | L | Harmonic band-limit. |
| n | Number of samples. |
| int ssht_sampling_mw_nphi | ( | int | L | ) |
Compute number of phi samples for McEwen and Wiaux sampling.
| [in] | L | Harmonic band-limit. |
| nphi | Number of phi samples. |
| int ssht_sampling_mw_ntheta | ( | int | L | ) |
Compute number of theta samples for McEwen and Wiaux sampling.
/note Computes number of samples in (0,pi], not over extended domain.
| [in] | L | Harmonic band-limit. |
| ntheta | Number of theta samples. |
| double ssht_sampling_mw_p2phi | ( | int | p, |
| int | L | ||
| ) |
Convert phi index to angle for McEwen and Wiaux sampling.
| [in] | p | Phi index. |
| [in] | L | Harmonic band-limit. |
| phi | Phi angle. |
| int ssht_sampling_mw_ss_n | ( | int | L | ) |
Compute total number of samples for McEwen and Wiaux symmetric sampling.
/note Computes number of samples on sphere, not over extended domain.
| [in] | L | Harmonic band-limit. |
| n | Number of samples. |
| int ssht_sampling_mw_ss_nphi | ( | int | L | ) |
Compute number of phi samples for McEwen and Wiaux symmetric sampling.
| [in] | L | Harmonic band-limit. |
| nphi | Number of phi samples. |
| int ssht_sampling_mw_ss_ntheta | ( | int | L | ) |
Compute number of theta samples for McEwen and Wiaux symmetric sampling.
/note Computes number of samples in [0,pi], not over extended domain.
| [in] | L | Harmonic band-limit. |
| ntheta | Number of theta samples. |
| double ssht_sampling_mw_ss_p2phi | ( | int | p, |
| int | L | ||
| ) |
Convert phi index to angle for McEwen and Wiaux symmetric sampling.
| [in] | p | Phi index. |
| [in] | L | Harmonic band-limit. |
| phi | Phi angle. |
| double ssht_sampling_mw_ss_t2theta | ( | int | t, |
| int | L | ||
| ) |
Convert theta index to angle for McEwen and Wiaux symmetric sampling.
| [in] | t | Theta index. |
| [in] | L | Harmonic band-limit. |
| theta | Theta angle. |
| double ssht_sampling_mw_t2theta | ( | int | t, |
| int | L | ||
| ) |
Convert theta index to angle for McEwen and Wiaux sampling.
| [in] | t | Theta index. |
| [in] | L | Harmonic band-limit. |
| theta | Theta angle. |
| double ssht_sampling_weight_dh | ( | double | theta_t, |
| int | L | ||
| ) |
Compute Driscoll and Healy weights.
| [in] | theta_t | Theta value to compute weight for. |
| [in] | L | Harmonic band-limit. return w Weight value computed. |
| ssht_complex_double ssht_sampling_weight_mw | ( | int | p | ) |
Compute weights for toroidal extension.
| [in] | p | Integer index to compute weight for. |
| [out] | w | Corresponding weight. |
| Weight | value computed. |