pyccl.tracers module

Tracers represent projected quantities that can be cross-correlated to get angular power spectra (see angular_cl()). In the most general case, the angular power spectrum between two tracers is given by

\[C^{\alpha\beta}_\ell=\frac{2}{\pi}\int d\chi_1\,d\chi_2\,dk\,k^2 P_{\alpha\beta}(k,\chi_1,\chi_2)\,\Delta^\alpha_\ell(k,\chi_1)\, \Delta^\beta_\ell(k,\chi_2),\]

where \(P_{\alpha\beta}\) is a generalized power spectrum (see Pk2D), and \(\Delta^\alpha_\ell(k,\chi)\) is a sum over different contributions associated to tracer \(\alpha\), where every contribution takes the form:

\[\Delta^\alpha_\ell(k,\chi)=f^\alpha_\ell\,W_\alpha(\chi)\, T_\alpha(k,\chi)\,j^{(n_\alpha)}_\ell(k\chi).\]

Here, \(f^\alpha_\ell\) is an \(\ell\)-dependent prefactor, \(W_\alpha(\chi)\) is the radial kernel, \(T_\alpha(k,\chi)\) is the transfer function, and \(j^{(n)}_\ell(x)\) is a generalized version of the spherical Bessel functions. Descriptions of each of these ingredients, and how to implement generalised Tracers as sub-classes of the Tracer base class can be found below. The documentation of the base Tracer class is a good place to start.

pyccl.tracers.get_density_kernel(cosmo, *, dndz)[source]

This convenience function returns the radial kernel for galaxy-clustering-like tracers. Given an unnormalized redshift distribution, it returns two arrays: \(\chi\), \(W(\chi)\), where \(\chi\) is an array of radial distances in units of Mpc and \(W(\chi) = p(z)\,H(z)\), where \(H(z)\) is the expansion rate in units of \({\rm Mpc}^{-1}\), and \(p(z)\) is the normalized redshift distribution.

Parameters:
  • cosmo (Cosmology) – cosmology object used to transform redshifts into distances.

  • dndz (tuple) – A tuple of arrays (z, N(z)) giving the redshift distribution of the objects. The units are arbitrary; N(z) will be normalized to unity.

pyccl.tracers.get_lensing_kernel(cosmo, *, dndz, mag_bias=None, n_chi=None)[source]

This convenience function returns the radial kernel for weak-lensing-like. Given an unnormalized redshift distribution and an optional magnification bias function, it returns two arrays: \(\chi\), \(W(\chi)\), where \(\chi\) is an array of radial distances in units of Mpc and \(W(\chi)\) is the lensing kernel:

\[W(\chi)=\frac{3H_0^2\Omega_m}{2a}\chi\, \int_{z(\chi)}^\infty dz\,\left(1-\frac{5s(z)}{2}\right)\, \frac{\chi(z)-\chi}{\chi(z)}\]

Note

If using this function to compute the magnification bias kernel, the result must be multiplied by -2.

Parameters:
  • cosmo (Cosmology) – cosmology object used to transform redshifts into distances.

  • dndz (tuple) – A tuple of arrays (z, N(z)) giving the redshift distribution of the objects. The units are arbitrary; N(z) will be normalized to unity.

  • mag_bias (tuple) – A tuple of arrays (z, s(z)) giving the magnification bias as a function of redshift. If None, s=0 will be assumed.

pyccl.tracers.get_kappa_kernel(cosmo, *, z_source, n_samples=100)[source]

This convenience function returns the radial kernel for CMB-lensing-like tracers.

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • z_source (float) – Redshift of source plane for CMB lensing.

  • n_samples (int) – number of samples over which the kernel is desired. These will be equi-spaced in radial distance. The kernel is quite smooth, so usually O(100) samples is enough.

class pyccl.tracers.Tracer(self)[source]

Bases: CCLObject

Tracers contain the information necessary to describe the contribution of a given sky observable to its cross-power spectrum with any other tracer. Tracers are composed of 4 main ingredients:

  • A radial kernel: this expresses the support in redshift/distance over which this tracer extends.

  • A transfer function: this is a function of wavenumber and scale factor that describes the connection between the tracer and the power spectrum on different scales and at different cosmic times.

  • An ell-dependent prefactor: normally associated with angular derivatives of a given fundamental quantity.

  • The order of the derivative of the Bessel functions with which they enter the computation of the angular power spectrum.

A Tracer object will in reality be a list of different such tracers that get combined linearly when computing power spectra.

property chi_min

Returns the minimum comoving distance over which this tracer’s radial kernel is defined, if it exists. For tracers with more than one contribution with an associated chi_min, the lowest value is returned.

property chi_max

Returns the maximum comoving distance over which this tracer’s radial kernel is defined, if it exists. For tracers with more than one contribution with an associated chi_max, the highest value is returned.

get_kernel(chi=None)[source]

Get the radial kernels for all tracers contained in this Tracer.

Parameters:

chi (float or array) – values of the comoving radial distance in increasing order and in Mpc. If None, returns the kernel at the internal spline nodes.

Returns:

list of radial kernels for each tracer. The shape will be (n_tracer, chi.size), where n_tracer is the number of tracers. The last dimension will be squeezed if the input is a scalar.

Return type:

array

get_f_ell(ell)[source]

Get the \(\ell\)-dependent prefactors for all tracers contained in this Tracer.

Parameters:

ell (float or array) – angular multipole values.

Returns:

list of prefactors for each tracer. The shape will be (n_tracer, ell.size), where n_tracer is the number of tracers. The last dimension will be squeezed if the input ell is a scalar.

Return type:

array

get_transfer(lk, a)[source]

Get the transfer functions for all tracers contained in this Tracer.

Parameters:
  • lk (float or array) – values of the natural logarithm of the wave number (in units of inverse Mpc) in increasing order.

  • a (float or array) – values of the scale factor.

Returns:

list of transfer functions for each tracer. The shape will be (n_tracer, lk.size, a.size), where n_tracer is the number of tracers. The other dimensions will be squeezed if the inputs are scalars.

Return type:

array

get_bessel_derivative()[source]

Get list of Bessel function derivative orders for all tracers contained in this Tracer.

Returns:

list of Bessel derivative orders for each tracer.

Return type:

array

get_angles_derivative()[source]

Get list of the \(\ell\)-dependent prefactor order for all tracers contained in this Tracer.

Returns:

list of angular derivative orders for each tracer.

Return type:

array

get_avg_weighted_a()[source]

Get list of kernel-averaged scale factors for all tracers contained in this Tracer.

Returns:

list of kernel-averaged scale factors for each tracer.

Return type:

array

add_tracer(cosmo, *, kernel=None, transfer_ka=None, transfer_k=None, transfer_a=None, der_bessel=0, der_angles=0, is_logt=False, extrap_order_lok=0, extrap_order_hik=2)[source]

Adds one more tracer to the list contained in this Tracer.

Parameters:
  • cosmo (Cosmology) – cosmology object.

  • kernel (tuple) – A tuple of arrays (chi, w_chi) describing the radial kernel of this tracer. chi should contain values of the comoving radial distance in increasing order, and w_chi should contain the values of the kernel at those values of the radial distance. The kernel will be assumed to be zero outside the range of distances covered by chi. If kernel is None a constant kernel \(W(\chi)=1\) will be assumed everywhere.

  • transfer_ka (tuple) – a tuple of arrays (a, lk,t_ka) describing the most general transfer function for a tracer. a should be an array of scale factor values in increasing order. lk should be an array of values of the natural logarithm of the wave number (in units of inverse Mpc) in increasing order. t_ka should be an array of shape (na, nk), where na and nk are the sizes of a and lk respectively.``t_ka`` should hold the values of the transfer function at the corresponding values of a and lk. If your transfer function is factorizable (i.e. \(T(a,k) = A(a)\, K(k))\), it is more efficient to set this to None and use transfer_k and transfer_a to describe \(K\) and \(A\) respectively. The transfer function will be assumed continuous and constant outside the range of scale factors covered by a. It will be extrapolated using polynomials of order extrap_order_lok and extrap_order_hik below and above the range of wavenumbers covered by lk respectively. If this argument is not None, the values of transfer_k and transfer_a will be ignored.

  • transfer_k (tuple) – a tuple of arrays (lk,t_k) describing the scale-dependent part of a factorizable transfer function. lk should be an array of values of the natural logarithm of the wave number (in units of inverse Mpc) in increasing order. t_k should be an array of the same size holding the values of the k-dependent part of the transfer function at those wavenumbers. It will be extrapolated using polynomials of order extrap_order_lok and extrap_order_hik below and above the range of wavenumbers covered by lk respectively. If None, the k-dependent part of the transfer function will be set to 1 everywhere.

  • transfer_a (tuple) – a tuple of arrays (a, t_a)` describing the time-dependent part of a factorizable transfer function. a should be an array of scale factor values in increasing order. t_a should contain the time-dependent part of the transfer function at those values of the scale factor. The time dependence will be assumed continuous and constant outside the range covered by a. If None, the time-dependent part of the transfer function will be set to 1 everywhere.

  • der_bessel (int) – order of the derivative of the Bessel functions with which this tracer enters the calculation of the power spectrum. Allowed values are -1, 0, 1 and 2. 0, 1 and 2 correspond to the raw functions, their first derivatives or their second derivatives. -1 corresponds to the raw functions divided by the square of their argument. We enable this special value because this type of dependence is ubiquitous for many common tracers (lensing, IAs), and makes the corresponding transfer functions more stable for small \(k\) or \(\chi\).

  • der_angles (int) – integer describing the ell-dependent prefactor associated with this tracer. Allowed values are 0, 1 and 2. 0 means no prefactor. 1 means a prefactor \(\ell(\ell+1)\), associated with the angular Laplacian and used e.g. for lensing convergence and magnification. 2 means a prefactor \(\sqrt{(\ell+2)!/(\ell-2)!}\), associated with the angular derivatives of spin-2 fields (e.g. cosmic shear, IAs).

  • is_logt (bool) – if True, transfer_ka, transfer_k and transfer_a will contain the natural logarithm of the transfer function (or their factorizable parts).

  • extrap_order_lok (int) – polynomial order used to extrapolate the transfer functions for low wavenumbers not covered by the input arrays.

  • extrap_order_hik (int) – polynomial order used to extrapolate the transfer functions for high wavenumbers not covered by the input arrays.

classmethod from_z_power(cosmo, *, A, alpha, z_min=0.0, z_max=6.0, n_chi=1024)[source]

Constructor for tracers associated with a radial kernel of the form

\[W(\chi) = \frac{A}{(1+z)^\alpha},\]

where \(A\) is an amplitude and \(\alpha\) is a power law index. The kernel only has support in the redshift range [z_min, z_max].

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • A (float) – amplitude parameter.

  • alpha (float) – power law index.

  • z_min (float) – minimum redshift from to which we define the kernel.

  • z_max (float) – maximum redshift up to which we define the kernel.

  • n_chi (float) – number of intervals in the radial comoving distance on which we sample the kernel.

class pyccl.tracers.NzTracer(self)[source]

Bases: Tracer

Specific base class for tracers with an internal _dndz redshift distribution interpolator. These include NumberCountsTracer() and WeakLensingTracer().

get_dndz(z)[source]

Get the redshift distribution for this tracer.

Parameters:

z (float or array) – redshift values.

Returns:

redshift distribution evaluated at the input values of z.

Return type:

array

pyccl.tracers.NumberCountsTracer(cosmo, *, dndz, bias=None, mag_bias=None, has_rsd, n_samples=256)[source]

Specific Tracer associated to galaxy clustering with linear scale-independent bias, including redshift-space distortions and magnification. The associated contributions are described in detail in Section 2.4.1 of the CCL paper.

Warning

When including redshift-space distortions, the current implementation assumes linear, scale-independent growth Although this should be valid in \(\Lambda\) CDM and on the large scales (especially when considering a broad \(N(z)\)), this approximation should be borne in mind.

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • dndz (tuple) – A tuple of arrays (z, N(z)) giving the redshift distribution of the objects. The units are arbitrary; N(z) will be normalized to unity.

  • bias (tuple) – A tuple of arrays (z, b(z)) giving the galaxy bias. If None, this tracer won’t include a term proportional to the matter density contrast.

  • mag_bias (tuple) – A tuple of arrays (z, s(z)) giving the magnification bias as a function of redshift. If None, the tracer is assumed to not have magnification bias terms.

  • has_rsd (bool) – If True, this tracer will include a redshift-space distortion term.

  • n_samples (int) – number of samples over which the magnification lensing kernel is desired. These will be equi-spaced in radial distance. The kernel is quite smooth, so usually O(100) samples is enough.

pyccl.tracers.WeakLensingTracer(cosmo, *, dndz, has_shear=True, ia_bias=None, use_A_ia=True, n_samples=256)[source]

Specific Tracer associated to galaxy shape distortions including lensing shear and intrinsic alignments within the L-NLA model. The associated contributions are described in detail in Section 2.4.1 of the CCL paper.

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • dndz (tuple) – A tuple of arrays (z, N(z)) giving the redshift distribution of the objects. The units are arbitrary; N(z) will be normalized to unity.

  • has_shear (bool) – set to False if you want to omit the lensing shear contribution from this tracer.

  • ia_bias (tuple) – A tuple of arrays (z, A_IA(z)) giving the intrinsic alignment amplitude A_IA(z). If None, the tracer is assumed to not have intrinsic alignments.

  • use_A_ia (bool) – set to True to use the conventional IA normalization. Set to False to use the raw input amplitude, which will usually be 1 for use with perturbation theory IA modeling.

  • n_samples (int) – number of samples over which the lensing kernel is desired. These will be equi-spaced in radial distance. The kernel is quite smooth, so usually O(100) samples is enough.

pyccl.tracers.CMBLensingTracer(cosmo, *, z_source, n_samples=100)[source]

A Tracer for CMB lensing convergence \(\kappa\). The associated kernel and transfer function are described in Eq. 31 of the CCL paper.

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • z_source (float) – Redshift of source plane for CMB lensing.

  • n_samples (int) – number of samples over which the kernel is desired. These will be equi-spaced in radial distance. The kernel is quite smooth, so usually O(100) samples is enough.

pyccl.tracers.tSZTracer(cosmo, *, z_max=6.0, n_chi=1024)[source]

Specific Tracer associated with the thermal Sunyaev Zel’dovich Compton-y parameter. The radial kernel for this tracer is simply given by

\[W(\chi) = \frac{\sigma_T}{m_ec^2} \frac{1}{1+z},\]

where \(\sigma_T\) is the Thomson scattering cross section and \(m_e\) is the electron mass.

Any angular power spectra computed with this tracer, should use a three-dimensional power spectrum involving the electron pressure in physical (non-comoving) units of \(eV\,{\rm cm}^{-3}\).

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • z_max (float) – maximum redshift up to which we define the kernel.

  • n_chi (float) – number of intervals in the radial comoving distance on which we sample the kernel.

pyccl.tracers.CIBTracer(cosmo, *, z_min=0.0, z_max=6.0, n_chi=1024)[source]

Specific Tracer associated with the cosmic infrared background (CIB). The radial kernel for this tracer is simply

\[W(\chi) = \frac{1}{1+z}.\]

Any angular power spectra computed with this tracer, should use a three-dimensional power spectrum involving the CIB emissivity density in units of \({\rm Jy}\,{\rm Mpc}^{-1}\,{\rm srad}^{-1}\) (or multiples thereof – see e.g. HaloProfileCIBShang12).

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • z_min (float) – minimum redshift down to which we define the kernel.

  • z_max (float) – maximum redshift up to which we define the kernel.

  • n_chi (float) – number of intervals in the radial comoving distance on which we sample the kernel.

pyccl.tracers.ISWTracer(cosmo, *, z_max=6.0, n_chi=1024)[source]

Specific Tracer associated with the integrated Sachs-Wolfe effect (ISW). Useful when cross-correlating any low-redshift probe with the primary CMB anisotropies. The ISW contribution to the temperature fluctuations is:

\[\Delta T_{\rm CMB} = 2T_{\rm CMB} \int_0^{\chi_{LSS}}d\chi a\,\dot{\phi}\]

Any angular power spectra computed with this tracer, should use a three-dimensional power spectrum involving the matter power spectrum.

Warning

The current implementation of this tracer assumes a standard Poisson equation relating \(\phi\) and \(\delta\), and linear, scale-independent structure growth. Although this should be valid in \(\Lambda\) CDM and on the large scales the ISW is sensitive to, these approximations must be borne in mind.

Parameters:
  • cosmo (Cosmology) – Cosmology object.

  • z_max (float) – maximum redshift up to which we define the kernel.

  • n_chi (float) – number of intervals in the radial comoving distance on which we sample the kernel.