SHTns 3.6.1
Functions
Pseudo-spectral rotations of Spherical Harmonic fields (deprecated).

deprecated use Rotations of Spherical Harmonic fields. instead. More...

Functions

void SH_Zrotate (shtns_cfg shtns, cplx *Qlm, double alpha, cplx *Rlm)
 Rotate a SH representation Qlm around the z-axis by angle alpha (in radians), which is the same as rotating the reference frame by angle -alpha.
 
void SH_Xrotate90 (shtns_cfg shtns, cplx *Qlm, cplx *Rlm)
 rotate Qlm by 90 degrees around X axis and store the result in Rlm.
 
void SH_Yrotate90 (shtns_cfg shtns, cplx *Qlm, cplx *Rlm)
 rotate Qlm by 90 degrees around Y axis and store the result in Rlm.
 
void SH_Yrotate (shtns_cfg shtns, cplx *Qlm, double alpha, cplx *Rlm)
 rotate Qlm around Y axis by arbitrary angle, using composition of rotations. Store the result in Rlm.
 

Rotation functions

shtns_rot shtns_rotation_create (const int lmax, const int mmax, int norm)
 Creation of a spherical harmonic rotation object for degrees and orders up to lmax and mmax respectively.
 
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)
 Defines a rotation parameterized by the 3 intrinsic Euler angles with ZYZ convention.
 
void shtns_rotation_set_angles_ZXZ (shtns_rot r, double alpha, double beta, double gamma)
 Defines a rotation parameterized by the 3 intrinsic Euler angles with ZXZ convention.
 
void shtns_rotation_set_angle_axis (shtns_rot r, double theta, double Vx, double Vy, double Vz)
 Defines a rotation of angle theta around axis of cartesian coordinates (Vx,Vy,Vz).
 
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_cplx (shtns_rot r, cplx *Zlm, cplx *Rlm)
 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)
 
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)
 

Detailed Description

deprecated use Rotations of Spherical Harmonic fields. instead.

Rotation around axis other than Z should be considered of beta quality (they have been tested but may still contain bugs). They also require mmax = lmax. They use an Algorithm inspired by the pseudospectral rotation described in Gimbutas Z. and Greengard L. 2009 "A fast and stable method for rotating spherical harmonic expansions" Journal of Computational Physics. doi:10.1016/j.jcp.2009.05.014

These functions do only require a call to shtns_create, but not to shtns_set_grid.

Function Documentation

◆ SH_Xrotate90()

void SH_Xrotate90 ( shtns_cfg  shtns,
cplx Qlm,
cplx Rlm 
)

rotate Qlm by 90 degrees around X axis and store the result in Rlm.

shtns->mres MUST be 1, and lmax=mmax.

◆ SH_Yrotate()

void SH_Yrotate ( shtns_cfg  shtns,
cplx Qlm,
double  alpha,
cplx Rlm 
)

rotate Qlm around Y axis by arbitrary angle, using composition of rotations. Store the result in Rlm.

Deprecated:
Rotate SH representation around Y axis by alpha (in radians).

◆ SH_Yrotate90()

void SH_Yrotate90 ( shtns_cfg  shtns,
cplx Qlm,
cplx Rlm 
)

rotate Qlm by 90 degrees around Y axis and store the result in Rlm.

shtns->mres MUST be 1, and lmax=mmax.

◆ SH_Zrotate()

void SH_Zrotate ( shtns_cfg  shtns,
cplx Qlm,
double  alpha,
cplx Rlm 
)

Rotate a SH representation Qlm around the z-axis by angle alpha (in radians), which is the same as rotating the reference frame by angle -alpha.

Result is stored in Rlm (which can be the same array as Qlm).

◆ shtns_rotation_apply_cplx()

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

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)

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 
)

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

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 
)

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

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 
)

Defines a rotation parameterized by the 3 intrinsic Euler angles with 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 
)

Defines a rotation parameterized by the 3 intrinsic Euler angles with 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.
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)