SHTns 3.6.1
Functions
Rotations of Spherical Harmonic fields.

Rotation of spherical harmonics, using an on-the-fly algorithm (does not store the rotation matrix) inspired by the GUMEROV's algorithm to generate the Wigner-d matrices describing rotation of Spherical Harmonics. More...

Functions

shtns_rot shtns_rotation_create (const int lmax, const int mmax, int norm)
 Allocate memory and precompute some recurrence coefficients for rotation (independent of angle).
 
void shtns_rotation_destroy (shtns_rot r)
 Release memory allocated for the rotation object.
 
void shtns_rotation_set_angles_ZYZ (shtns_rot r, double alpha, double beta, double gamma)
 Set the rotation angles, and compute associated Legendre functions, given the 3 intrinsic Euler angles in ZYZ convention.
 
void shtns_rotation_set_angles_ZXZ (shtns_rot r, double alpha, double beta, double gamma)
 Set the rotation angles, and compute associated Legendre functions, given the 3 intrinsic Euler angles in ZXZ convention.
 
void shtns_rotation_set_angle_axis (shtns_rot r, double theta, double Vx, double Vy, double Vz)
 Sets a rotation by angle theta around axis of cartesian coordinates (Vx,Vy,Vz), and compute associated Legendre functions.
 
int quarter_wigner_d_matrix (shtns_rot r, const int l, double *mx, const int compressed)
 
int shtns_rotation_wigner_d_matrix (shtns_rot r, const int l, double *mx)
 Generate spherical-harmonic rotation matrix for given degree l and orthonormal convention (Wigner-d matrix)
 
void shtns_rotation_apply_real (shtns_rot r, cplx *Qlm, cplx *Rlm)
 apply rotation to the orthonormal spherical harmonic expansion Qlm of a real field, store the result into Rlm (can be the same as Qlm)
 
void shtns_rotation_apply_cplx (shtns_rot r, cplx *Zlm, cplx *Rlm)
 rotate Zlm, store the result in Rlm, without ever storing the wigner-d matrices (on-the-fly operation)
 

Detailed Description

Rotation of spherical harmonics, using an on-the-fly algorithm (does not store the rotation matrix) inspired by the GUMEROV's algorithm to generate the Wigner-d matrices describing rotation of Spherical Harmonics.

See https://arxiv.org/abs/1403.7698 or https://doi.org/10.1007/978-3-319-13230-3_5 Thanks to Alex J. Yuffa, ayuff.nosp@m.a@gm.nosp@m.ail.c.nosp@m.om for his suggestions and help.

These functions are more accurate and faster than the pseudo-spectral rotations of Gimbutas and do not require mmax=lmax. Note that if mmax<lmax, the rotations are apprixmate only, with the quality of the approximation decreasing with increasing 'beta' angle (rotation angle around Y-axis).

Function Documentation

◆ shtns_rotation_apply_cplx()

void shtns_rotation_apply_cplx ( shtns_rot  r,
cplx Zlm,
cplx Rlm 
)

rotate Zlm, store the result in Rlm, without ever storing the wigner-d matrices (on-the-fly operation)

apply rotation to the orthonormal spherical harmonic expansion Zlm of a complex-valued field, store the result into Rlm (can be the same as Zlm)

◆ shtns_rotation_create()

shtns_rot shtns_rotation_create ( const int  lmax,
const int  mmax,
int  norm 
)

Allocate memory and precompute some recurrence coefficients for rotation (independent of angle).

Creation of a spherical harmonic rotation object for degrees and orders up to lmax and mmax respectively.

Setting mmax < lmax will result in approximate rotations if not aligned with the z-axis, while mmax=lmax leads to exact rotations.

◆ shtns_rotation_set_angle_axis()

void shtns_rotation_set_angle_axis ( shtns_rot  r,
double  theta,
double  Vx,
double  Vy,
double  Vz 
)

Sets a rotation by angle theta around axis of cartesian coordinates (Vx,Vy,Vz), and compute associated Legendre functions.

Defines a rotation of angle theta around axis of cartesian coordinates (Vx,Vy,Vz).

◆ shtns_rotation_set_angles_ZXZ()

void shtns_rotation_set_angles_ZXZ ( shtns_rot  r,
double  alpha,
double  beta,
double  gamma 
)

Set the rotation angles, and compute associated Legendre functions, given the 3 intrinsic Euler angles in ZXZ convention.

Defines a rotation parameterized by the 3 intrinsic Euler angles with ZXZ convention.

◆ shtns_rotation_set_angles_ZYZ()

void shtns_rotation_set_angles_ZYZ ( shtns_rot  r,
double  alpha,
double  beta,
double  gamma 
)

Set the rotation angles, and compute associated Legendre functions, given the 3 intrinsic Euler angles in ZYZ convention.

Defines a rotation parameterized by the 3 intrinsic Euler angles with ZYZ convention.

◆ shtns_rotation_wigner_d_matrix()

int shtns_rotation_wigner_d_matrix ( shtns_rot  r,
const int  l,
double *  mx 
)

Generate spherical-harmonic rotation matrix for given degree l and orthonormal convention (Wigner-d matrix)

Parameters
[out]mxis an (2*l+1)*(2*l+1) array that will be filled with the Wigner-d matrix elements (rotation matrix along Y-axis in orthonormal spherical harmonic space).
Returns
0 if error, or 2*l+1 (size of the square matrix) otherwise.
Warning
The returned rotation matrix only applies to orthonormal convention (see Normalization)