s2wav is currently in an open beta, please provide feedback on GitHub

Filter functions#

s2wav.filters.filters_axisym(L: int, J_min: int = 0, lam: float = 2.0) Tuple[ndarray, ndarray]#
Computes wavelet kernels \(\Psi^j_{\ell m}\) and scaling kernel

\(\Phi_{\ell m}\) in harmonic space.

Specifically, these kernels are derived in [1], where the wavelet kernels are defined (15) for scale \(j\) to be

\[\Psi^j_{\ell m} \equiv \sqrt{\frac{2\ell+1}{4\pi}} \kappa_{\lambda}(\frac{\ell}{\lambda^j})\delta_{m0},\]

where \(\kappa_{\lambda} = \sqrt{k_{\lambda}(t/\lambda) - k_{\lambda}(t)}\) for \(k_{\lambda}\) given in k_lam(). Similarly, the scaling kernel is defined (16) as

\[\Phi_{\ell m} \equiv \sqrt{\frac{2\ell+1}{4\pi}} \nu_{\lambda} (\frac{\ell}{\lambda^{J_0}})\delta_{m0},\]

where \(\nu_{\lambda} = \sqrt{k_{\lambda}(t)}\) for \(k_{\lambda}\) given in k_lam(). Notice that \(\delta_{m0}\) enforces that these kernels are axisymmetric, i.e. coefficients for \(m \not = \ell\) are zero. In this implementation the normalisation constant has been omitted as it is nulled in subsequent functions.

Parameters:
  • L (int) – Harmonic band-limit.

  • J_min (int, optional) – Lowest frequency wavelet scale to be used. Defaults to 0.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

Raises:

ValueError – J_min is negative or greater than J.

Returns:

Unnormalised wavelet kernels \(\Psi^j_{\ell m}\) with shape \([(J+1)L]\), and scaling kernel \(\Phi_{\el m}\) with shape \([L]\) in harmonic space.

Return type:

Tuple[np.ndarray, np.ndarray]

Note

[1] B. Leidstedt et. al., “S2LET: A code to perform fast wavelet analysis on the sphere”, A&A, vol. 558, p. A128, 2013.

s2wav.filters.filters_axisym_jax(L: int, J_min: int = 0, lam: float = 2.0) Tuple[Array, Array]#

JAX version of filters_axisym_vectorised().

Parameters:
  • L (int) – Harmonic band-limit.

  • J_min (int, optional) – Lowest frequency wavelet scale to be used. Defaults to 0.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

Raises:

ValueError – J_min is negative or greater than J.

Returns:

Unnormalised wavelet kernels \(\Psi^j_{\ell m}\) with shape \([(J+1)L], and scaling kernel :math:\)Phi_{ell m}` with shape \([L]\) in harmonic space.

Return type:

Tuple[np.ndarray, np.ndarray]

s2wav.filters.filters_axisym_vectorised(L: int, J_min: int = 0, lam: float = 2.0) Tuple[ndarray, ndarray]#

Vectorised version of filters_axisym().

Parameters:
  • L (int) – Harmonic band-limit.

  • J_min (int, optional) – Lowest frequency wavelet scale to be used. Defaults to 0.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

Raises:

ValueError – J_min is negative or greater than J.

Returns:

Unnormalised wavelet kernels \(\Psi^j_{\ell m}\) with shape \([(J+1)L], and scaling kernel :math:\)Phi_{ell m}` with shape \([L]\) in harmonic space.

Return type:

Tuple[np.ndarray, np.ndarray]

s2wav.filters.filters_directional(L: int, N: int = 1, J_min: int = 0, lam: float = 2.0, spin: int = 0, spin0: int = 0, using_torch: bool = False) Tuple[ndarray, ndarray]#

Generates the harmonic coefficients for the directional tiling wavelets.

This implementation is based on equation 36 in the wavelet computation paper [1].

Parameters:
  • L (int) – Harmonic band-limit.

  • N (int, optional) – Upper azimuthal band-limit. Defaults to 1.

  • J_min (int, optional) – Lowest frequency wavelet scale to be used. Defaults to 0.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

  • spin (int, optional) – Spin (integer) to perform the transform. Defaults to 0.

  • spin0 (int, optional) – Spin number the wavelet was lowered from. Defaults to 0.

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

Returns:

Tuple of wavelet and scaling kernels

(\(\Psi^j_{\ell n}\), \(\Phi_{\ell m}\))

Return type:

Tuple[np.ndarray, np.ndarray]

Notes

[1] J. McEwen et. al., “Directional spin wavelets on the sphere”, arXiv preprint arXiv:1509.06749 (2015).

s2wav.filters.filters_directional_jax(L: int, N: int = 1, J_min: int = 0, lam: float = 2.0, spin: int = 0, spin0: int = 0) Tuple[Array, Array]#

JAX version of filters_directional().

Parameters:
  • L (int) – Harmonic band-limit.

  • N (int, optional) – Upper azimuthal band-limit. Defaults to 1.

  • J_min (int, optional) – Lowest frequency wavelet scale to be used. Defaults to 0.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

  • spin (int, optional) – Spin (integer) to perform the transform. Defaults to 0.

  • spin0 (int, optional) – Spin number the wavelet was lowered from. Defaults to 0.

Returns:

Tuple of wavelet and scaling kernels

(\(\Psi^j_{\ell n}\), \(\Phi_{\ell m}\)).

Return type:

Tuple[np.ndarray, np.ndarray]

s2wav.filters.filters_directional_vectorised(L: int, N: int = 1, J_min: int = 0, lam: float = 2.0, spin: int = 0, spin0: int = 0, using_torch: bool = False) Tuple[ndarray, ndarray]#

Vectorised version of filters_directional().

Parameters:
  • L (int) – Harmonic band-limit.

  • N (int, optional) – Upper azimuthal band-limit. Defaults to 1.

  • J_min (int, optional) – Lowest frequency wavelet scale to be used. Defaults to 0.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

  • spin (int, optional) – Spin (integer) to perform the transform. Defaults to 0.

  • spin0 (int, optional) – Spin number the wavelet was lowered from. Defaults to 0.

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

Returns:

Tuple of wavelet and scaling kernels

(\(\Psi^j_{\ell n}\), \(\Phi_{\ell m}\)).

Return type:

Tuple[np.ndarray, np.ndarray]

s2wav.filters.k_lam(L: int, lam: float = 2.0, quad_iters: int = 300) float#

Compute function \(k_{\lambda}\) used as a wavelet generating function.

Specifically, this function is derived in [1] and is given by

\[k_{\lambda} \equiv \frac{ \int_t^1 \frac{\text{d}t^{\prime}}{t^{\prime}} s_{\lambda}^2(t^{\prime})}{ \int_{\frac{1}{\lambda}}^1 \frac{\text{d}t^{\prime}}{t^{\prime}} s_{\lambda}^2(t^{\prime})},\]

where the integrand is defined to be

\[s_{\lambda} \equiv s \Big ( \frac{2\lambda}{\lambda - 1}(t-\frac{1}{\lambda}) - 1 \Big ),\]

for infinitely differentiable Cauchy-Schwartz function \(s(t) \in C^{\infty}\).

Parameters:
  • L (int) – Harmonic band-limit.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

  • quad_iters (int, optional) – Total number of iterations for quadrature integration. Defaults to 300.

Returns:

Value of \(k_{\lambda}\) computed for values between

\(\frac{1}{\lambda}\) and 1, parametrised by \(\ell\) as required to compute the axisymmetric filters in tiling_axisym().

Return type:

(np.ndarray)

Note

[1] B. Leidstedt et. al., “S2LET: A code to perform fast wavelet analysis on the

sphere”, A&A, vol. 558, p. A128, 2013.

s2wav.filters.k_lam_jax(L: int, lam: float = 2.0, quad_iters: int = 300) float#
JAX version of k_lam. Compute function \(k_{\lambda}\) used as a wavelet

generating function.

Specifically, this function is derived in [1] and is given by

\[k_{\lambda} \equiv \frac{ \int_t^1 \frac{\text{d}t^{\prime}}{t^{\prime}} s_{\lambda}^2(t^{\prime})}{ \int_{\frac{1}{\lambda}}^1 \frac{\text{d}t^{\prime}}{t^{\prime}} s_{\lambda}^2(t^{\prime})},\]

where the integrand is defined to be

\[s_{\lambda} \equiv s \Big ( \frac{2\lambda}{\lambda - 1}(t-\frac{1}{\lambda}) - 1 \Big ),\]

for infinitely differentiable Cauchy-Schwartz function \(s(t) \in C^{\infty}\).

Parameters:
  • L (int) – Harmonic band-limit.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales. Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

  • quad_iters (int, optional) – Total number of iterations for quadrature integration. Defaults to 300.

Returns:

Value of \(k_{\lambda}\) computed for values between

\(\frac{1}{\lambda}\) and 1, parametrised by \(\ell\) as required to compute the axisymmetric filters in tiling_axisym().

Return type:

(np.ndarray)

Note

[1] B. Leidstedt et. al., “S2LET: A code to perform fast wavelet analysis on the

sphere”, A&A, vol. 558, p. A128, 2013.

s2wav.filters.part_scaling_fn(a: float, b: float, n: int, lam: float = 2.0) float#

Computes integral used to calculate smoothly decreasing function \(k_{\lambda}\).

Intermediate step used to compute the wavelet and scaling function generating functions. Uses the trapezium method to integrate tiling_integrand() in the limits from \(a \rightarrow b\) with scaling parameter \(\lambda\). One of the basic mathematical functions needed to carry out the tiling of the harmonic space.

Parameters:
  • a (float) – Lower limit of the numerical integration.

  • b (float) – Upper limit of the numerical integration.

  • n (int) – Number of steps to be performed during integration.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales.Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

Returns:

Integral of the tiling integrand from \(a \rightarrow b\).

Return type:

float

s2wav.filters.tiling_direction(L: int, N: int = 1) ndarray#
Generates the harmonic coefficients for the directionality component of the

tiling functions.

Formally, this function implements the follow equation

\[_{s}\eta_{\el m} = \nu \vu \sqrt{\frac{1}{2^{\gamma}} \big ( \binom{\gamma}{ (\gamma - m)/2} \big )}\]

which was first derived in [1].

Parameters:
  • L (int) – Harmonic band-limit.

  • N (int, optional) – Upper orientational band-limit. Defaults to 1.

Returns:

Harmonic coefficients of directionality components

\(_{s}\eta_{\el m}\).

Return type:

np.ndarray

Notes

[1] J. McEwen et. al., “Directional spin wavelets on the sphere”, arXiv preprint

arXiv:1509.06749 (2015).

s2wav.filters.tiling_direction_jax(L: int, N: int = 1) ndarray#
JAX version of tiling_direction. Generates the harmonic coefficients for the

directionality component of the tiling functions.

Formally, this function implements the follow equation

\[_{s}\eta_{\ell m} = \nu \vu \sqrt{\frac{1}{2^{\gamma}} \big ( \binom{\gamma}{ (\gamma - m)/2} \big )}\]

which was first derived in [1].

Parameters:
  • L (int) – Harmonic band-limit.

  • N (int, optional) – Upper orientational band-limit. Defaults to 1.

Returns:

Harmonic coefficients of directionality components

\(_{s}\eta_{\ell m}\).

Return type:

np.ndarray

Notes

[1] J. McEwen et. al., “Directional spin wavelets on the sphere”, arXiv preprint

arXiv:1509.06749 (2015).

s2wav.filters.tiling_integrand(t: float, lam: float = 2.0) float#

Tiling integrand for scale-discretised wavelets [1].

Intermediate step used to compute the wavelet and scaling function generating functions. One of the basic mathematical functions needed to carry out the tiling of the harmonic space.

Parameters:
  • t (float) – Real argument over which we integrate.

  • lam (float, optional) – Wavelet parameter which determines the scale factor between consecutive wavelet scales.Note that \(\lambda = 2\) indicates dyadic wavelets. Defaults to 2.

Returns:

Value of tiling integrand for given \(t\) and scaling factor.

Return type:

float

Note

[1] B. Leidstedt et. al., “S2LET: A code to perform fast wavelet analysis on

the sphere”, A&A, vol. 558, p. A128, 2013.