tjpcov.wigner_transform#

Module Contents#

Classes#

WignerTransform

Class to compute curved sky Hankel transforms with wigner-d matrices.

Functions#

wigner_d(s1, s2, theta, ell[, l_use_bessel])

Function to compute the wigner-d matrices.

wigner_d_parallel(s1, s2, theta, ell[, ncpu, l_use_bessel])

Compute the wigner-d matrix in parallel using multiprocessing Pool.

bin_cov(r, cov, r_bins)

Function to apply the binning operator.

class tjpcov.wigner_transform.WignerTransform(theta, ell, s1_s2, ncpu=None)[source]#

Class to compute curved sky Hankel transforms with wigner-d matrices.

Initialize the class for the given angles, scales and spins.

Parameters:
  • theta – Values of angular separation, theta, at which the Hankel transform is done. Should be in radians.

  • ell – ell values at which the Hankel transform is done. Should be integers

  • s1_s2

    List of spin pairs of the tracers. Each spin pair should be a tuple. e.g. for 3X2 analysis, pass [(0,0),(0,2),(2,2),(2,-2)].

    (0,0): (galaxy,galaxy) (0,2): (galaxy,shear). (2,0) is equivalent. (2,2): (shear,shear), xi+ (2,-2): (shear,shear), xi-

  • ncpu – Number of python processes to use when computing wigner-d matrices.

cl_grid(ell_cl, cl, taper=False, **taper_kwargs)[source]#

Interpolate the input C_ell.

This is done in case the ell values of C_ell are different from the grid on which wigner-d matrices were computed during intialization.

Parameters:
  • ell_cl (array) – ell at which the input C_ell is computed.

  • cl (array) – input C_ell

  • taper (bool) – if True apply the tapering to the input C_ell. Tapering can help in reducing ringing.

Returns:

C_ell evaluated at the initialization ells.

Return type:

array

cl_cov_grid(ell_cl, cl_cov, taper=False, **taper_kwargs)[source]#

Interpolate the input 2D covariances.

This is done in case in case the ell values are different from the grid on which wigner-d matrices were computed during intialization.

Parameters:
  • ell_cl (array) – ell at which the input covariance was computed.

  • cl (array) – input covariance

  • taper (bool) – if True apply the tapering to the input C_ell. Tapering can help in reducing ringing.

  • **taper_kwargs – Arguments to pass to the tapering method

Returns:

covariance evaluated at the initialization ells.

Return type:

array

projected_covariance(ell_cl, cl_cov, s1_s2, s1_s2_cross=None, taper=False, **kwargs)[source]#

Convert C_ell covariance to correlation function.

This function assumes that cl_cov is 2D matrix.

Parameters:
  • cl_cov – C_ell covariance matrix.

  • ell_cl – ell values at which input C_ell is computed.

  • s1_s2 – Tuple of the spin factors of the first set of tracers. Used to identify the correct wigner-d matrix to use.

  • s1_s2_cross – Tuple of the spin factors of the second set of tracers, if different from s1_s2. Used to identify the correct wigner-d matrix to use.

  • taper (bool) – if True apply the tapering to the input C_ell. Tapering can help in reducing ringing.

  • **kwargs – Arguments to pass to the tapering method

Returns:

  • theta (array): angles given at initialization

  • cov (array): real space covariance at the given angles

Return type:

tuple

abstract taper(ell, large_k_lower=10, large_k_upper=100, low_k_lower=0, low_k_upper=1e-05)[source]#

Apply tapering to input C_ell.

Tapering is useful to reduce the ringing. This function uses the cosine function to apply the tapering. See eq. 71 in https://arxiv.org/pdf/2105.04548.pdf for the function and meaning of input parameters.

Parameters:
  • ell – ell values at which input C_ell is computed.

  • large_k_lower

  • large_k_upper

  • low_k_lower

  • low_k_upper

diagonal_err(cov)[source]#

Returns the diagonal error from the covariance.

Useful for errorbar plots.

Parameters:

cov (array) – Covariance

Returns:

Diagonal errors (i.e. sqrt(diag(cov)))

Return type:

array

tjpcov.wigner_transform.wigner_d(s1, s2, theta, ell, l_use_bessel=10000.0)[source]#

Function to compute the wigner-d matrices.

Parameters:
  • s1 – Spin factors for the wigner-d matrix.

  • s2 – Spin factors for the wigner-d matrix.

  • theta – Angular separation for which to compute the wigner-d matrix. The matrix depends on cos(theta).

  • ell – The spherical harmonics mode ell for which to compute the matrix.

  • l_use_bessel – Due to numerical issues, we need to switch from wigner-d matrix to bessel functions at high ell (see the note below). This defines the scale at which the switch happens.

Returns:

Wigner-d matrix

Return type:

array

tjpcov.wigner_transform.wigner_d_parallel(s1, s2, theta, ell, ncpu=None, l_use_bessel=10000.0)[source]#

Compute the wigner-d matrix in parallel using multiprocessing Pool.

This function calls the wigner-d function defined above.

Parameters:
  • s1 – Spin factors for the wigner-d matrix.

  • s2 – Spin factors for the wigner-d matrix.

  • theta – Angular separation for which to compute the wigner-d matrix. The matrix depends on cos(theta).

  • ell – The spherical harmonics mode ell for which to compute the matrix.

  • ncpu – number of processes to use for computing the matrix.

  • l_use_bessel – Due to numerical issues, we need to switch from wigner-d matrix to bessel functions at high ell (see the note below). This defines the scale at which the switch happens.

Returns:

Wigner-d matrix

Return type:

array

tjpcov.wigner_transform.bin_cov(r, cov, r_bins)[source]#

Function to apply the binning operator.

This function works on both one dimensional vectors and two dimensional covariance covrices.

Parameters:
  • r – theta or ell values at which the un-binned vector is computed.

  • cov – Unbinned covariance. It also works for a vector of C_ell or xi

  • r_bins – theta or ell bins to which the values should be binned.

Returns:

Binned covariance or vector of C_ell or xi

Return type:

array_like