My Project
Functions
ssht_core.h File Reference
#include <complex.h>

Go to the source code of this file.

Functions

void ssht_core_mw_inverse_sov_sym (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_inverse_sov_sym (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_real (double *f, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_inverse_sov_sym_real (double *f, const SSHT_COMPLEX(double) *flm, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_forward_sov_conv_sym (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_real (SSHT_COMPLEX(double) *flm, const double *f, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_forward_sov_conv_sym_real (SSHT_COMPLEX(double) *flm, const double *f, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_pole (SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *f_sp, double *phi_sp, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_real_pole (double *f, double *f_sp, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_pole (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) f_sp, double phi_sp, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_real_pole (SSHT_COMPLEX(double) *flm, const double *f, double f_sp, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mwdirect_inverse (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
 
void ssht_core_mwdirect_inverse_sov (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_ss (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_inverse_sov_sym_ss (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_ss_real (double *f, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_inverse_sov_sym_ss_real (double *f, const SSHT_COMPLEX(double) *flm, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_ss (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_forward_sov_conv_sym_ss (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_ss_real (SSHT_COMPLEX(double) *flm, const double *f, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_lb_forward_sov_conv_sym_ss_real (SSHT_COMPLEX(double) *flm, const double *f, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_ss_pole (SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *f_np, double *phi_np, SSHT_COMPLEX(double) *f_sp, double *phi_sp, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_inverse_sov_sym_ss_real_pole (double *f, double *f_np, double *f_sp, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_ss_pole (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) f_np, double phi_np, SSHT_COMPLEX(double) f_sp, double phi_sp, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mw_forward_sov_conv_sym_ss_real_pole (SSHT_COMPLEX(double) *flm, const double *f, double f_np, double f_sp, int L, ssht_dl_method_t dl_method, int verbosity)
 
void ssht_core_mwdirect_inverse_ss (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
 
void ssht_core_gl_inverse_sov (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
 
void ssht_core_gl_inverse_sov_real (double *f, const SSHT_COMPLEX(double) *flm, int L, int verbosity)
 
void ssht_core_gl_forward_sov (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, int verbosity)
 
void ssht_core_gl_forward_sov_real (SSHT_COMPLEX(double) *flm, const double *f, int L, int verbosity)
 
void ssht_core_dh_inverse_sov (SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
 
void ssht_core_dh_inverse_sov_real (double *f, const SSHT_COMPLEX(double) *flm, int L, int verbosity)
 
void ssht_core_dh_forward_sov (SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, int verbosity)
 
void ssht_core_dh_forward_sov_real (SSHT_COMPLEX(double) *flm, const double *f, int L, int verbosity)
 

Function Documentation

◆ ssht_core_dh_forward_sov()

void ssht_core_dh_forward_sov ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
int  L,
int  spin,
int  verbosity 
)

Compute forward transform using Driscoll and Healy quadrature with separation of variables.

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 4236 of file ssht_core.c.

References SSHT_COMPLEX, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, and SSHT_PROMPT.

Referenced by main().

◆ ssht_core_dh_forward_sov_real()

void ssht_core_dh_forward_sov_real ( SSHT_COMPLEX(double) *  flm,
const double *  f,
int  L,
int  verbosity 
)

Compute forward transform of real scalar signal using Driscoll and Healy quadrature with separation of variables (symmetries for real signals are exploited).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 4373 of file ssht_core.c.

References SSHT_COMPLEX, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, and SSHT_PROMPT.

Referenced by main().

◆ ssht_core_dh_inverse_sov()

void ssht_core_dh_inverse_sov ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
int  verbosity 
)

Compute inverse transform using direct method with separation of variables for DH sampling.

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3996 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_beta_kostelec_line_table(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_dh_t2theta().

Referenced by main().

◆ ssht_core_dh_inverse_sov_real()

void ssht_core_dh_inverse_sov_real ( double *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  verbosity 
)

Compute inverse transform of real scalar signal using direct method with separation of variables for DH sampling (symmetries for real signals are exploited).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 4115 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_beta_kostelec_halfline_table(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_dh_t2theta().

Referenced by main().

◆ ssht_core_gl_forward_sov()

void ssht_core_gl_forward_sov ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
int  L,
int  spin,
int  verbosity 
)

Compute forward transform using Gauss-Legendgre quadrature with separation of variables.

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3685 of file ssht_core.c.

References SSHT_COMPLEX, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_gl_thetas_weights().

Referenced by main().

◆ ssht_core_gl_forward_sov_real()

void ssht_core_gl_forward_sov_real ( SSHT_COMPLEX(double) *  flm,
const double *  f,
int  L,
int  verbosity 
)

Compute forward transform of real scalar signal using Gauss-Legendgre quadrature with separation of variables (symmetries for real signals are exploited).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3832 of file ssht_core.c.

References SSHT_COMPLEX, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_gl_thetas_weights().

Referenced by main().

◆ ssht_core_gl_inverse_sov()

void ssht_core_gl_inverse_sov ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
int  verbosity 
)

Compute inverse transform using direct method with separation of variables for GL sampling.

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3421 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_beta_kostelec_line_table(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_gl_thetas_weights().

Referenced by main().

◆ ssht_core_gl_inverse_sov_real()

void ssht_core_gl_inverse_sov_real ( double *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  verbosity 
)

Compute inverse transform of real scalar signal using direct method with separation of variables for GL sampling (symmetries for real signals are exploited).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3552 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_beta_kostelec_halfline_table(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_gl_thetas_weights().

Referenced by main().

◆ ssht_core_mw_forward_sov_conv_sym()

void ssht_core_mw_forward_sov_conv_sym ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (for complex spin signal).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 810 of file ssht_core.c.

References ssht_core_mw_lb_forward_sov_conv_sym().

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_pole().

◆ ssht_core_mw_forward_sov_conv_sym_pole()

void ssht_core_mw_forward_sov_conv_sym_pole ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
SSHT_COMPLEX(double)  f_sp,
double  phi_sp,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

South pole wrapper for forward transform for MW method. The South pole is defined by a single sample and its corresponding phi angle, rather than specifying samples for all phi at the South pole (which are simply related by the rotation of a spin function in its tangent plane).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere (excluding South pole).
[in]f_spFunction sample on South pole.
[in]phi_spPhi angle corresponding to quoted sample at South pole.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1702 of file ssht_core.c.

References SSHT_COMPLEX, ssht_core_mw_forward_sov_conv_sym(), SSHT_ERROR_MEM_ALLOC_CHECK, and ssht_sampling_mw_p2phi().

Referenced by main().

◆ ssht_core_mw_forward_sov_conv_sym_real()

void ssht_core_mw_forward_sov_conv_sym_real ( SSHT_COMPLEX(double) *  flm,
const double *  f,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method of real scalar signal using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1203 of file ssht_core.c.

References ssht_core_mw_lb_forward_sov_conv_sym_real().

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_real_pole().

◆ ssht_core_mw_forward_sov_conv_sym_real_pole()

void ssht_core_mw_forward_sov_conv_sym_real_pole ( SSHT_COMPLEX(double) *  flm,
const double *  f,
double  f_sp,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

South pole wrapper for forward transform of real scalar function for MW method. The South pole is defined by a single sample, rather than specifying samples for all phi at the South pole (which for a scalar function are identical).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere (excluding South pole).
[in]f_spFunction sample on South pole.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1750 of file ssht_core.c.

References ssht_core_mw_forward_sov_conv_sym_real(), and SSHT_ERROR_MEM_ALLOC_CHECK.

Referenced by main().

◆ ssht_core_mw_forward_sov_conv_sym_ss()

void ssht_core_mw_forward_sov_conv_sym_ss ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method with symmetric sampling using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (for complex spin signal).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2425 of file ssht_core.c.

References ssht_core_mw_lb_forward_sov_conv_sym_ss().

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_ss_pole().

◆ ssht_core_mw_forward_sov_conv_sym_ss_pole()

void ssht_core_mw_forward_sov_conv_sym_ss_pole ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
SSHT_COMPLEX(double)  f_np,
double  phi_np,
SSHT_COMPLEX(double)  f_sp,
double  phi_sp,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

North-South pole wrapper for forward transform for MW method with symmetric sampling. The poles are defined by single samples and their corresponding phi angle, rather than specifying samples for all phi at the poles (which are simply related by the rotation of a spin function in its tangent plane).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere (excluding poles).
[in]f_npFunction sample on North pole.
[in]phi_npPhi angle corresponding to quoted sample at North pole.
[in]f_spFunction sample on South pole.
[in]phi_spPhi angle corresponding to quoted sample at South pole.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3316 of file ssht_core.c.

References SSHT_COMPLEX, ssht_core_mw_forward_sov_conv_sym_ss(), SSHT_ERROR_MEM_ALLOC_CHECK, and ssht_sampling_mw_ss_p2phi().

Referenced by main().

◆ ssht_core_mw_forward_sov_conv_sym_ss_real()

void ssht_core_mw_forward_sov_conv_sym_ss_real ( SSHT_COMPLEX(double) *  flm,
const double *  f,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method using symmetric sampling of real scalar signal using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]LHarmonic band-limit.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2820 of file ssht_core.c.

References ssht_core_mw_lb_forward_sov_conv_sym_ss_real().

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_ss_real_pole().

◆ ssht_core_mw_forward_sov_conv_sym_ss_real_pole()

void ssht_core_mw_forward_sov_conv_sym_ss_real_pole ( SSHT_COMPLEX(double) *  flm,
const double *  f,
double  f_np,
double  f_sp,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

North-South pole wrapper for forward transform of real scalar function for MW method with symmetric sampling. The poles are defined by single samples, rather than specifying samples for all phi at the poles (which for a scalar function are identical).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere (excluding poles).
[in]f_npFunction sample on North pole.
[in]f_spFunction sample on South pole.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3369 of file ssht_core.c.

References ssht_core_mw_forward_sov_conv_sym_ss_real(), and SSHT_ERROR_MEM_ALLOC_CHECK.

Referenced by main().

◆ ssht_core_mw_inverse_sov_sym()

void ssht_core_mw_inverse_sov_sym ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method using separation of variables, fast Fourier transforms and exploiting all symmetries (for complex spin signal).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 51 of file ssht_core.c.

References ssht_core_mw_lb_inverse_sov_sym().

Referenced by main(), and ssht_core_mw_inverse_sov_sym_pole().

◆ ssht_core_mw_inverse_sov_sym_pole()

void ssht_core_mw_inverse_sov_sym_pole ( SSHT_COMPLEX(double) *  f,
SSHT_COMPLEX(double) *  f_sp,
double *  phi_sp,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

South pole wrapper for inverse transform for MW method. The South pole is defined by a single sample and its corresponding phi angle, rather than specifying samples for all phi at the South pole (which are simply related by the rotation of a spin function in its tangent plane).

Parameters
[out]fFunction on sphere (excluding South pole).
[out]f_spFunction sample on South pole.
[out]phi_spPhi angle corresponding to quoted sample at South pole.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1608 of file ssht_core.c.

References SSHT_COMPLEX, ssht_core_mw_inverse_sov_sym(), SSHT_ERROR_MEM_ALLOC_CHECK, and ssht_sampling_mw_p2phi().

Referenced by main().

◆ ssht_core_mw_inverse_sov_sym_real()

void ssht_core_mw_inverse_sov_sym_real ( double *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method of real scalar signal using separation of variables, fast Fourier transforms and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 322 of file ssht_core.c.

References ssht_core_mw_lb_inverse_sov_sym_real().

Referenced by main(), and ssht_core_mw_inverse_sov_sym_real_pole().

◆ ssht_core_mw_inverse_sov_sym_real_pole()

void ssht_core_mw_inverse_sov_sym_real_pole ( double *  f,
double *  f_sp,
const SSHT_COMPLEX(double) *  flm,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

South pole wrapper for inverse transform of real scalar function for MW method. The South pole is defined by a single sample, rather than specifying samples for all phi at the South pole (which for a scalar function are identical).

Parameters
[out]fFunction on sphere (excluding South pole).
[out]f_spFunction sample on South pole.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1654 of file ssht_core.c.

References ssht_core_mw_inverse_sov_sym_real(), and SSHT_ERROR_MEM_ALLOC_CHECK.

Referenced by main().

◆ ssht_core_mw_inverse_sov_sym_ss()

void ssht_core_mw_inverse_sov_sym_ss ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method with symmetric sampling using separation of variables, fast Fourier transforms and exploiting all symmetries (for complex spin signal).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1798 of file ssht_core.c.

References ssht_core_mw_lb_inverse_sov_sym_ss().

Referenced by main(), and ssht_core_mw_inverse_sov_sym_ss_pole().

◆ ssht_core_mw_inverse_sov_sym_ss_pole()

void ssht_core_mw_inverse_sov_sym_ss_pole ( SSHT_COMPLEX(double) *  f,
SSHT_COMPLEX(double) *  f_np,
double *  phi_np,
SSHT_COMPLEX(double) *  f_sp,
double *  phi_sp,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

North-South pole wrapper for inverse transform for MW method with symmetric sampling. The poles are defined by single samples and their corresponding phi angle, rather than specifying samples for all phi at the poles (which are simply related by the rotation of a spin function in its tangent plane).

Parameters
[out]fFunction on sphere (excluding poles).
[out]f_npFunction sample on North pole.
[out]phi_npPhi angle corresponding to quoted sample at North pole.
[out]f_spFunction sample on South pole.
[out]phi_spPhi angle corresponding to quoted sample at South pole.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3209 of file ssht_core.c.

References SSHT_COMPLEX, ssht_core_mw_inverse_sov_sym_ss(), SSHT_ERROR_MEM_ALLOC_CHECK, and ssht_sampling_mw_ss_p2phi().

Referenced by main().

◆ ssht_core_mw_inverse_sov_sym_ss_real()

void ssht_core_mw_inverse_sov_sym_ss_real ( double *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method with symmetric sampling of real scalar signal using separation of variables, fast Fourier transforms and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2068 of file ssht_core.c.

References ssht_core_mw_lb_inverse_sov_sym_ss_real().

Referenced by main(), and ssht_core_mw_inverse_sov_sym_ss_real_pole().

◆ ssht_core_mw_inverse_sov_sym_ss_real_pole()

void ssht_core_mw_inverse_sov_sym_ss_real_pole ( double *  f,
double *  f_np,
double *  f_sp,
const SSHT_COMPLEX(double) *  flm,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

North-South pole wrapper for inverse transform of real scalar function for MW method with symmetric sampling. The poles are defined by single samples, rather than specifying samples for all phi at the poles (which for a scalar function are identical).

Parameters
[out]fFunction on sphere (excluding poles).
[out]f_spFunction sample on South pole.
[out]f_npFunction sample on North pole.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 3261 of file ssht_core.c.

References ssht_core_mw_inverse_sov_sym_ss_real(), and SSHT_ERROR_MEM_ALLOC_CHECK.

Referenced by main().

◆ ssht_core_mw_lb_forward_sov_conv_sym()

void ssht_core_mw_lb_forward_sov_conv_sym ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
int  L0,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (for complex spin signal).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 837 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym().

◆ ssht_core_mw_lb_forward_sov_conv_sym_real()

void ssht_core_mw_lb_forward_sov_conv_sym_real ( SSHT_COMPLEX(double) *  flm,
const double *  f,
int  L0,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method of real scalar signal using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1230 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_real().

◆ ssht_core_mw_lb_forward_sov_conv_sym_ss()

void ssht_core_mw_lb_forward_sov_conv_sym_ss ( SSHT_COMPLEX(double) *  flm,
const SSHT_COMPLEX(double) *  f,
int  L0,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method with symmetric sampling using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (for complex spin signal).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2452 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_ss().

◆ ssht_core_mw_lb_forward_sov_conv_sym_ss_real()

void ssht_core_mw_lb_forward_sov_conv_sym_ss_real ( SSHT_COMPLEX(double) *  flm,
const double *  f,
int  L0,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute forward transform for MW method using symmetric sampling of real scalar signal using separation of variables, fast Fourier transforms, performing convolution with weights as product in transformed space and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]flmHarmonic coefficients.
[in]fFunction on sphere.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2847 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_forward_sov_conv_sym_ss_real().

◆ ssht_core_mw_lb_inverse_sov_sym()

void ssht_core_mw_lb_inverse_sov_sym ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L0,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method using separation of variables, fast Fourier transforms and exploiting all symmetries (for complex spin signal).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 77 of file ssht_core.c.

References MAX, SSHT_COMPLEX, ssht_dl_beta_risbo_eighth_table(), ssht_dl_beta_risbo_fill_eighth2quarter_table(), ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), ssht_dl_halfpi_trapani_eighth_table(), ssht_dl_halfpi_trapani_fill_eighth2quarter_table(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_DL_TRAPANI, SSHT_ERROR_GENERIC, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_inverse_sov_sym().

◆ ssht_core_mw_lb_inverse_sov_sym_real()

void ssht_core_mw_lb_inverse_sov_sym_real ( double *  f,
const SSHT_COMPLEX(double) *  flm,
int  L0,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method of real scalar signal using separation of variables, fast Fourier transforms and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 348 of file ssht_core.c.

References MAX, SSHT_COMPLEX, ssht_dl_beta_risbo_eighth_table(), ssht_dl_beta_risbo_fill_eighth2quarter_table(), ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), ssht_dl_halfpi_trapani_eighth_table(), ssht_dl_halfpi_trapani_fill_eighth2quarter_table(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_DL_TRAPANI, SSHT_ERROR_GENERIC, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_inverse_sov_sym_real().

◆ ssht_core_mw_lb_inverse_sov_sym_ss()

void ssht_core_mw_lb_inverse_sov_sym_ss ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L0,
int  L,
int  spin,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method with symmetric sampling using separation of variables, fast Fourier transforms and exploiting all symmetries (for complex spin signal).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 1824 of file ssht_core.c.

References MAX, SSHT_COMPLEX, ssht_dl_beta_risbo_eighth_table(), ssht_dl_beta_risbo_fill_eighth2quarter_table(), ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), ssht_dl_halfpi_trapani_eighth_table(), ssht_dl_halfpi_trapani_fill_eighth2quarter_table(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_DL_TRAPANI, SSHT_ERROR_GENERIC, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_inverse_sov_sym_ss().

◆ ssht_core_mw_lb_inverse_sov_sym_ss_real()

void ssht_core_mw_lb_inverse_sov_sym_ss_real ( double *  f,
const SSHT_COMPLEX(double) *  flm,
int  L0,
int  L,
ssht_dl_method_t  dl_method,
int  verbosity 
)

Compute inverse transform for MW method with symmetric sampling of real scalar signal using separation of variables, fast Fourier transforms and exploiting all symmetries (including additional symmetries for real signals).

Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]L0Lower harmonic band-limit.
[in]LUpper harmonic band-limit.
[in]spinSpin number.
[in]dl_methodMethod to use when compute Wigner functions.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2095 of file ssht_core.c.

References MAX, SSHT_COMPLEX, ssht_dl_beta_risbo_eighth_table(), ssht_dl_beta_risbo_fill_eighth2quarter_table(), ssht_dl_calloc(), ssht_dl_get_offset(), ssht_dl_get_stride(), ssht_dl_halfpi_trapani_eighth_table(), ssht_dl_halfpi_trapani_fill_eighth2quarter_table(), SSHT_DL_QUARTER, SSHT_DL_QUARTER_EXTENDED, SSHT_DL_RISBO, SSHT_DL_TRAPANI, SSHT_ERROR_GENERIC, SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PION2, and SSHT_PROMPT.

Referenced by main(), and ssht_core_mw_inverse_sov_sym_ss_real().

◆ ssht_core_mwdirect_inverse()

void ssht_core_mwdirect_inverse ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
int  verbosity 
)

Compute inverse transform using direct method for MW sampling.

Warning
This algorithm is very slow and is included for verification purposes only.
Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 590 of file ssht_core.c.

References ssht_dl_beta_risbo_full_table(), ssht_dl_calloc(), SSHT_DL_FULL, ssht_dl_get_offset(), ssht_dl_get_stride(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_mw_t2theta().

◆ ssht_core_mwdirect_inverse_sov()

void ssht_core_mwdirect_inverse_sov ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
int  verbosity 
)

Compute inverse transform using direct method with separation of variables for MW sampling.

Warning
This algorithm is included for verification purposes only.
Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 689 of file ssht_core.c.

References SSHT_COMPLEX, ssht_dl_beta_kostelec_line_table(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_mw_t2theta().

◆ ssht_core_mwdirect_inverse_ss()

void ssht_core_mwdirect_inverse_ss ( SSHT_COMPLEX(double) *  f,
const SSHT_COMPLEX(double) *  flm,
int  L,
int  spin,
int  verbosity 
)

Compute inverse transform using direct method for MW symmetric sampling.

Warning
This algorithm is very slow and is included for verification purposes only.
Parameters
[out]fFunction on sphere.
[in]flmHarmonic coefficients.
[in]LHarmonic band-limit.
[in]spinSpin number.
[in]verbosityVerbosity flag in range [0,5].
Return values
none
Author
Jason McEwen

Definition at line 2326 of file ssht_core.c.

References ssht_dl_beta_risbo_full_table(), ssht_dl_calloc(), SSHT_DL_FULL, ssht_dl_get_offset(), ssht_dl_get_stride(), SSHT_ERROR_MEM_ALLOC_CHECK, SSHT_PI, SSHT_PROMPT, and ssht_sampling_mw_ss_t2theta().