utils#
- s2ball.utils.elm2ind(el: int, m: int) int #
Convert from spherical harmonic 2D indexing of \((\ell,m)\) to 1D index.
1D index is defined by el**2 + el + m.
Warning
Note that 1D storage of spherical harmonic coefficients is not the default.
- Parameters:
el (int) – Harmonic degree \(\ell\).
m (int) – Harmonic order \(m\).
- Returns:
Corresponding 1D index value.
- Return type:
int
- s2ball.utils.elmn2ind(el: int, m: int, n: int, L: int, N: int) int #
Convert from Wigner space 3D indexing of \((\ell,m, n)\) to 1D index. :param el: Harmonic degree \(\ell\). :type el: int :param m: Harmonic order \(m\). :type m: int :param n: Directional order \(n\). :type n: int :param L: Harmonic band-limit. :type L: int :param N: Directional band-limit. Defaults to 1. :type N: int, optional
- Returns:
Corresponding 1D index in Wigner space.
- Return type:
int
- s2ball.utils.flm_1d_to_2d(flm_1d: ndarray, L: int) ndarray #
Convert from 1D indexed harmnonic coefficients to 2D indexed coefficients.
Note
Storage conventions for harmonic coefficients \(flm_{(\ell,m)}\), for e.g. \(L = 3\), are as follows.
\[\begin{split}\text{ 2D data format}: \begin{bmatrix} 0 & 0 & flm_{(0,0)} & 0 & 0 \\ 0 & flm_{(1,-1)} & flm_{(1,0)} & flm_{(1,1)} & 0 \\ flm_{(2,-2)} & flm_{(2,-1)} & flm_{(2,0)} & flm_{(2,1)} & flm_{(2,2)} \end{bmatrix}\end{split}\]\[\text{1D data format}: [flm_{0,0}, flm_{1,-1}, flm_{1,0}, flm_{1,1}, \dots]\]- Parameters:
flm_1d (np.ndarray) – 1D indexed harmonic coefficients.
L (int) – Harmonic band-limit.
- Returns:
2D indexed harmonic coefficients.
- Return type:
np.ndarray
- s2ball.utils.flm_2d_to_1d(flm_2d: ndarray, L: int) ndarray #
Convert from 2D indexed harmonic coefficients to 1D indexed coefficients.
Note
Storage conventions for harmonic coefficients \(flm_{(\ell,m)}\), for e.g. \(L = 3\), are as follows.
\[\begin{split}\text{ 2D data format}: \begin{bmatrix} 0 & 0 & flm_{(0,0)} & 0 & 0 \\ 0 & flm_{(1,-1)} & flm_{(1,0)} & flm_{(1,1)} & 0 \\ flm_{(2,-2)} & flm_{(2,-1)} & flm_{(2,0)} & flm_{(2,1)} & flm_{(2,2)} \end{bmatrix}\end{split}\]\[\text{1D data format}: [flm_{0,0}, flm_{1,-1}, flm_{1,0}, flm_{1,1}, \dots]\]- Parameters:
flm_2d (np.ndarray) – 2D indexed harmonic coefficients.
L (int) – Harmonic band-limit.
- Returns:
1D indexed harmonic coefficients.
- Return type:
np.ndarray
- s2ball.utils.flmn_1d_to_3d(flmn_1d: ndarray, L: int, N: int) ndarray #
Convert from 1D indexed Wigner coefficients to 3D indexed coefficients. :param flm_1d: 1D indexed Wigner coefficients, C flatten index priority
\(n, \ell, m\).
- Parameters:
L (int) – Harmonic band-limit.
N (int) – Directional band-limit.
- Raises:
ValueError – flmn is already 3D indexed.
ValueError – flmn is not 1D.
- Returns:
3D indexed Wigner coefficients, index order \([n, \ell, m]\).
- Return type:
np.ndarray
- s2ball.utils.flmn_3d_to_1d(flmn_3d: ndarray, L: int, N: int) ndarray #
Convert from 3D indexed Wigner coefficients to 1D indexed coefficients. :param flm_3d: 3D indexed Wigner coefficients, index order
\([n, \ell, m]\).
- Parameters:
L (int) – Harmonic band-limit.
N (int, optional) – Directional band-limit.
- Raises:
ValueError – flmn is already 1D indexed.
ValueError – flmn is not 3D.
- Returns:
1D indexed Wigner coefficients, C flatten index priority \(n, \ell, m\).
- Return type:
np.ndarray
- s2ball.utils.generate_flm(rng: Generator, L: int, L_lower: int = 0, spin: int = 0, reality: bool = False) ndarray #
Generate a 2D set of random harmonic coefficients.
Note
Real signals are explicitly produced from conjugate symmetry.
- Parameters:
rng (Generator) – Random number generator.
L (int) – Harmonic band-limit.
L_lower (int, optional) – Harmonic lower bound. Defaults to 0.
spin (int, optional) – Harmonic spin. Defaults to 0.
reality (bool, optional) – Reality of signal. Defaults to False.
- Returns:
Random set of spherical harmonic coefficients.
- Return type:
np.ndarray
- s2ball.utils.generate_flmn(rng: Generator, L: int, N: int = 1, L_lower: int = 0, reality: bool = False) ndarray #
Generate a 3D set of random Wigner coefficients.
Note
Real signals are explicitly produced from conjugate symmetry.
- Parameters:
rng (Generator) – Random number generator.
L (int) – Harmonic band-limit.
N (int) – Directional band-limit. Must be < L.
L_lower (int, optional) – Harmonic lower bound. Defaults to 0.
reality (bool, optional) – Reality of signal. Defaults to False.
- Returns:
Random set of Wigner coefficients.
- Return type:
np.ndarray
- s2ball.utils.generate_flmnp(rng: Generator, L: int, N: int, P: int) ndarray #
Generate a 4D set of random Wigner-Laguerre coefficients.
- Parameters:
rng (Generator) – Random number generator.
L (int) – Harmonic band-limit.
N (int) – Directional band-limit. Must be < L.
P (int) – Radial band-limit.
- Returns:
Random set of Spherical-Laguerre coefficients.
- Return type:
np.ndarray
- s2ball.utils.generate_flmp(rng: Generator, L: int, P: int) ndarray #
Generate a 3D set of random Spherical-Laguerre coefficients.
- Parameters:
rng (Generator) – Random number generator.
L (int) – Harmonic band-limit.
P (int) – Radial band-limit.
- Returns:
Random set of Spherical-Laguerre coefficients.
- Return type:
np.ndarray
- s2ball.utils.ind2elm(ind: int) tuple #
Convert from 1D spherical harmonic index to 2D index of \((\ell,m)\).
Warning
Note that 1D storage of spherical harmonic coefficients is not the default.
- Parameters:
ind (int) – 1D spherical harmonic index.
- Returns:
(el,m) defining spherical harmonic degree and order.
- Return type:
tuple
- s2ball.utils.ncoeff(L: int) int #
Number of spherical harmonic coefficients for given band-limit L.
- Parameters:
L (int, optional) – Harmonic band-limit.
- Returns:
Number of spherical harmonic coefficients.
- Return type:
int