signal generator#

s2fft.utils.signal_generator.complex_el_and_m_indices(L: int, min_el: int) tuple[ndarray, ndarray]#

Generate pairs of el, m indices for accessing complex harmonic coefficients.

Equivalent to nested list-comprehension based implementation

``` el_indices, m_indices = np.array(

[(el, m) for el in range(min_el, L) for m in range(1, el + 1)]

).T#

For L, min_el = 1024, 0, this implementation is around 80x quicker in benchmarks compared to list-comprehension implementation.

param L:

Harmonic band-limit.

param min_el:

Inclusive lower-bound for el indices.

returns:

Tuple (el_indices, m_indices) with both entries 1D integer-valued NumPy arrays of same size, with values of corresponding entries corresponding to pairs of el and m indices.

s2fft.utils.signal_generator.complex_normal(rng: np.random.Generator, size: int | tuple[int], var: float) np.ndarray#

Generate array of samples from zero-mean complex normal distribution.

For z ~ ComplexNormal(0, var) we have that imag(z) ~ Normal(0, var/2) and real(z) ~ Normal(0, var/2) where Normal(μ, σ²) is the (real-valued) normal distribution with mean parameter μ and variance parameter σ².

Parameters:
  • rng – Numpy random generator object to generate samples using.

  • size – Output shape of array to generate.

  • var – Variance of complex normal distribution to generate samples from.

Returns:

Complex-valued array of shape size contained generated samples.

s2fft.utils.signal_generator.generate_flm(rng: np.random.Generator, L: int, L_lower: int = 0, spin: int = 0, reality: bool = False, using_torch: bool = False) np.ndarray | torch.Tensor#

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.

  • using_torch (bool, optional) – Desired frontend functionality. Defaults to False.

Returns:

Random set of spherical harmonic coefficients.

Return type:

np.ndarray

s2fft.utils.signal_generator.generate_flmn(rng: np.random.Generator, L: int, N: int = 1, L_lower: int = 0, reality: bool = False, using_torch: bool = False) np.ndarray | torch.Tensor#

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, optional) – Number of Fourier coefficients for tangent plane rotations (i.e. directionality). Defaults to 1.

  • L_lower (int, optional) – Harmonic lower bound. Defaults to 0.

  • reality (bool, optional) – Reality of signal. Defaults to False.

  • using_torch (bool, optional) – Desired frontend functionality. Defaults to False.

Returns:

Random set of Wigner coefficients.

Return type:

np.ndarray