epsproc.sphCalc module

ePSproc spherical function calculations.

Collection of functions for calculating Ylm, wignerD etc.

For spherical harmonics, currently using scipy.special.sph_harm

For other functions, using Moble’s spherical_functions package https://github.com/moble/spherical_functions

See tests/Spherical function testing Aug 2019.ipynb

27/08/19 Added wDcalc for Wigner D functions. 14/08/19 v1 Implmented sphCalc

epsproc.sphCalc.sphCalc(Lmax, Lmin=0, res=None, angs=None, XFlag=True)[source]

Calculate set of spherical harmonics Ylm(theta,phi) on a grid.

Parameters:
  • Lmax (int) – Maximum L for the set. Ylm calculated for Lmin:Lmax, all m.
  • Lmin (int, optional, default 0) – Min L for the set. Ylm calculated for Lmin:Lmax, all m.
  • res (int, optional, default None) – (Theta, Phi) grid resolution, outputs will be of dim [res,res].
  • angs (list of 2D np.arrays, [thetea, phi], optional, default None) – If passed, use these grids for calculation
  • XFlag (bool, optional, default True) – Flag for output. If true, output is Xarray. If false, np.arrays
  • that either res OR angs needs to be passed. (Note) –
  • Outputs
  • -------
  • if XFlag - (-) –
  • YlmX – 3D Xarray, dims (lm,theta,phi)
  • else - (-) –
  • lm (Ylm,) – 3D np.array of values, dims (lm,theta,phi), plus list of lm pairs
Currently set for scipy.special.sph_harm as calculation routine.

Example

>>> YlmX = sphCalc(2, res = 50)
epsproc.sphCalc.wDcalc(Lrange=[0, 1], Nangs=None, eAngs=None, R=None, XFlag=True)[source]

Calculate set of Wigner D functions D(l,m,mp,R) on a grid.

Parameters:
  • Lrange (list, optional, default [0, 1]) – Range of L to calculate parameters for. If len(Lrange) == 2 assumed to be of form [Lmin, Lmax], otherwise list is used directly. For a given l, all (m, mp) combinations are calculated.
  • for setting angles (use one only) (Options) –
  • Nangs (int, optional, default None) – If passed, use this to define Euler angles sampled. Ranges will be set as (theta, phi, chi) = (0:pi, 0:pi/2, 0:pi) in Nangs steps.
  • eAngs (np.array, optional, default None) – If passed, use this to define Euler angles sampled. Array of angles, [theta,phi,chi], in radians
  • R (np.array, optional, default None) – If passed, use this to define Euler angles sampled. Array of quaternions, as given by quaternion.from_euler_angles(eAngs).
XFlag : bool, optional, default True
Flag for output. If true, output is Xarray. If false, np.arrays
  • if XFlag -
wDX
Xarray, dims (lmmp,Euler)
  • else -
wD, R, lmmp
np.arrays of values, dims (lmmp,Euler), plus list of angles and lmmp sets.
Uses Moble's spherical_functions package for wigner D function.
https://github.com/moble/spherical_functions
Moble's quaternion package for angles and conversions.
https://github.com/moble/quaternion

Examples

>>> wDX1 = wDcalc(eAngs = np.array([0,0,0]))
>>> wDX2 = wDcalc(Nangs = 10)