Representation#
- s2scat.representation.scatter(flm: Array, L: int, N: int, J_min: int = 0, reality: bool = True, config: List[Array] = None, norm: List[Array] = None, recursive: bool = False, isotropic: bool = False, delta_j: int = None) List[Array] #
Compute directional scattering covariances on the sphere.
- Parameters:
flm (jnp.ndarray) – Spherical harmonic coefficients.
L (int) – Spherical harmonic bandlimit.
N (int) – Azimuthal bandlimit (directionality).
J_min (int, optional) – Minimum dyadic wavelet scale to consider. Defaults to 0.
reality (bool, optional) – Whether \(f \in \mathbb{R}\), if True exploits hermitian symmetry of harmonic coefficients. Defaults to True.
config (List[jnp.ndarray], optional) – All necessary precomputed arrays. Defaults to None.
norm (List[jnp.ndarray], optional) – Covariance normalisation values. Defaults to None.
recursive (bool, optional) – Whether to perform a memory efficient recursive transform, or a faster but less memory efficient fully precompute transform. Defaults to False.
isotropic (bool, optional) – Whether to return isotropic coefficients, i.e. average over directionality. Defaults to False.
delta_j (int, optional) – Range of wavelet scales over which to compute covariances. If None, covariances between all scales will be considered. Defaults to None.
- Raises:
ValueError – If one does not pass configuration arrays.
- Returns:
Directional scattering covariance statistics.
- Return type:
Tuple[jnp.ndarray]
Notes
The recursive transform, outlined in Price & McEwen (2023), requires \(\mathcal{O}(NL^2)\) memory overhead and can scale to high bandlimits \(L\). Conversely, the fully precompute transform requires \(\mathcal{O}(NL^3)\) memory overhead which can be large. However, the transform will be much faster. For applications at \(L \leq 512\) the precompute approach is a better choice, beyond which we recommend the users switch to recursive transforms or the C backend functionality.
If isotropic is true, the statistics will be contracted across \(n\). This will dramatically compress the covariance representation, but will be somewhat less sensitive to directional structure.
- s2scat.representation.scatter_c(flm: Array, L: int, N: int, J_min: int = 0, reality: bool = False, config: List[Array] | None = None, norm: List[Array] | None = None, isotropic: bool = False, delta_j: int | None = None) List[Array] #
Compute directional scattering covariances on the sphere using a custom C backend.
- Parameters:
flm (jnp.ndarray) – Spherical harmonic coefficients.
L (int) – Spherical harmonic bandlimit.
N (int) – Azimuthal bandlimit (directionality).
J_min (int, optional) – Minimum dyadic wavelet scale to consider. Defaults to 0.
reality (bool, optional) – Whether \(f \in \mathbb{R}\), if True exploits hermitian symmetry of harmonic coefficients. Defaults to False.
config (List[jnp.ndarray], optional) – All necessary precomputed arrays. Defaults to None.
norm (List[jnp.ndarray], optional) – Covariance normalisation values. Defaults to None.
isotropic (bool, optional) – Whether to return isotropic coefficients, i.e. average over directionality. Defaults to False.
delta_j (int, optional) – Range of wavelet scales over which to compute covariances. If None, covariances between all scales will be considered. Defaults to None.
- Raises:
ValueError – If one does not pass an array of wavelet filters.
- Returns:
Directional scattering covariance statistics.
- Return type:
Tuple[jnp.ndarray]
Notes
This variant of the directional scattering covariance transform leverages the JAX frontend for highly optimised C spherical harmonic libraries provided by S2FFT. As such, it is currently limited to CPU compute and cannot be JIT compiled. However, this approach can still be very fast as the underlying spherical harmonic libraries are extremely optimised. Reverse mode gradient functionality is supported, peak memory overhead is \(\mathcal{O}(NL^2)\), and this variant can scale to very high \(L \geq 4096\).
If isotropic is true, the statistics will be contracted across \(n\). This will dramatically compress the covariance representation, but will be somewhat less sensitive to directional structure.