pyccl.cosmology module

The core functionality of CCL, including the core data types lives in this module. Its focus is the Cosmology class, which plays a central role, carrying the information on cosmological parameters and derived quantities needed in most of the calculations carried out by CCL.

Note

All of the standalone functions in other modules, which take cosmo as their first argument, are methods of Cosmology.

Some important CCL parameters, governing for example the precision and speed of some calculations, are independent of the Cosmology objects, and instead can be accessed at a global level. You can do so as e.g. pyccl.gsl_params['ODE_GROWTH_EPSREL'], pyccl.spline_params['K_MIN'], pyccl.physical_constants['CLIGHT'], or pyccl.gsl_params.ODE_GROWTH_EPSREL, pyccl.spline_params.K_MIN, pyccl.physical_constants.CLIGHT.

class pyccl.cosmology.TransferFunctions(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

BBKS = 'bbks'
EISENSTEIN_HU = 'eisenstein_hu'
EISENSTEIN_HU_NOWIGGLES = 'eisenstein_hu_nowiggles'
BOLTZMANN_CLASS = 'boltzmann_class'
BOLTZMANN_CAMB = 'boltzmann_camb'
BOLTZMANN_ISITGR = 'boltzmann_isitgr'
CALCULATOR = 'calculator'
EMULATOR_LINPK = 'emulator'
class pyccl.cosmology.MatterPowerSpectra(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]

Bases: Enum

LINEAR = 'linear'
HALOFIT = 'halofit'
CAMB = 'camb'
CALCULATOR = 'calculator'
EMULATOR_NLPK = 'emulator'
class pyccl.cosmology.Cosmology(self, *, Omega_c=None, Omega_b=None, h=None, n_s=None, sigma8=None, A_s=None, Omega_k=0.0, Omega_g=None, Neff=None, m_nu=0.0, mass_split='normal', w0=-1.0, wa=0.0, T_CMB=2.7255, T_ncdm=0.71611, transfer_function='boltzmann_camb', matter_power_spectrum='halofit', baryonic_effects=None, mg_parametrization=None, extra_parameters=None)[source]

Bases: CCLObject

Stores information about cosmological parameters and associated data (e.g. distances, power spectra).

The values of cosmological parameters may be looked up by name (e.g. cosmo["sigma8"]). Note that some of the parameters accessible this way are not contained in the signature of Cosmology, but are derived during initialization.

Note

Although some arguments default to None, they will raise a ValueError inside this function if not specified, so they are not optional.

Note

The parameter Omega_g can be used to set the radiation density (not including relativistic neutrinos) to zero. Doing this will give you a model that is physically inconsistent since the temperature of the CMB will still be non-zero.

Parameters:
  • Omega_c (float) – Cold dark matter density fraction.

  • Omega_b (float) – Baryonic matter density fraction.

  • h (float) – Hubble constant divided by 100 km/s/Mpc; unitless.

  • A_s (float) – Power spectrum normalization. Exactly one of A_s and sigma_8 is required.

  • sigma8 (float) – Variance of matter density perturbations at an 8 Mpc/h scale. Exactly one of A_s and sigma_8 is required. Note that, if a value of sigma8 is passed, CCL will enforce the linear matter power spectrum to be correctly normalised to this value of \(\sigma_8\), even in the presence of other parameters (e.g. modified gravity parameters) that might affect the overall power spectrum normalization.

  • n_s (float) – Primordial scalar perturbation spectral index.

  • Omega_k (float) – Curvature density fraction. Defaults to 0.

  • Omega_g (float) – Density in relativistic species except massless neutrinos. The default of None corresponds to setting this from the CMB temperature. Note that if a non-None value is given, this may result in a physically inconsistent model because the CMB temperature will still be non-zero in the parameters.

  • Neff (float) – Effective number of massless neutrinos present. Defaults to 3.044.

  • m_nu (float or array) – Mass in eV of the massive neutrinos present. Defaults to 0. If a sequence is passed, it is assumed that the elements of the sequence represent the individual neutrino masses.

  • mass_split (str) – Type of massive neutrinos. Should be one of ‘single’, ‘equal’, ‘normal’, ‘inverted’. ‘single’ treats the mass as being held by one massive neutrino. The other options split the mass into 3 massive neutrinos. Ignored if a sequence is passed in m_nu. Default is ‘normal’.

  • w0 (float) – First order term of dark energy equation of state. Defaults to -1.

  • wa (float) – Second order term of dark energy equation of state. Defaults to 0.

  • T_CMB (float) – The CMB temperature today. The default value is 2.7255.

  • T_ncdm (float) – Non-CDM temperature in units of photon temperature. The default is 0.71611.

  • transfer_function (str or EmulatorPk) – The transfer function to use. Defaults to ‘boltzmann_camb’.

  • matter_power_spectrum (str or EmulatorPk) – The matter power spectrum to use. Defaults to ‘halofit’.

  • baryonic_effects (Baryons or None) – The baryonic effects model to use. Options are None (no baryonic effects), or a Baryons object.

  • mg_parametrization (ModifiedGravity or None) – The modified gravity parametrization to use. Options are None (no MG), or a ModifiedGravity object. Currently, only MuSigmaMG is supported.

  • extra_parameters (dict) – Dictionary holding extra parameters. Currently supports extra parameters for CAMB. Details described below. Defaults to None.

Currently supported extra parameters for CAMB are:

  • halofit_version

  • HMCode_A_baryon

  • HMCode_eta_baryon

  • HMCode_logT_AGN

  • kmax

  • lmax

  • dark_energy_model

Consult the CAMB documentation for their usage. These parameters are passed in a dict to extra_parameters as:

extra_parameters = {"camb": {"halofit_version": "mead2020_feedback",
                             "HMCode_logT_AGN": 7.8}}

Note

If using camb to compute the non-linear power spectrum with HMCode to include baryonic effects, you should not include any extra baryonic effects (i.e. set baryonic_effects=None).

write_yaml(filename, *, sort_keys=False)[source]

Write a YAML representation of the parameters to file.

Parameters:

filename (str) – file name, file pointer, or stream to write parameters to.

classmethod read_yaml(filename, **kwargs)[source]

Read the parameters from a YAML file.

Parameters:
  • filename (str) – file name, file pointer, or stream to read parameters from.

  • **kwargs (dict) – additional keywords that supersede file contents

compute_growth()[source]

Compute the growth function.

compute_linear_power()[source]

Compute the linear power spectrum.

compute_nonlin_power()[source]

Compute the non-linear power spectrum.

compute_sigma()[source]

Compute the sigma(M) spline.

get_linear_power(name='delta_matter:delta_matter')[source]

Get the Pk2D object associated with the linear power spectrum with name name.

Parameters:

name (str or None) – name of the power spectrum to return.

Returns:

Pk2D object containing the linear power spectrum with name name.

get_nonlin_power(name='delta_matter:delta_matter')[source]

Get the Pk2D object associated with the non-linear power spectrum with name name.

Parameters:

name (str or None) – name of the power spectrum to return.

Returns:

Pk2D object containing the non-linear power spectrum with name name.

property has_distances

Checks if the distances have been precomputed.

property has_growth

Checks if the growth function has been precomputed.

property has_linear_power

Checks if the linear power spectra have been precomputed.

property has_nonlin_power

Checks if the non-linear power spectra have been precomputed.

property has_sigma

Checks if sigma(M) is precomputed.

CIBTracer(*, z_min=0.0, z_max=6.0, n_chi=1024)

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.

CMBLensingTracer(*, z_source, n_samples=100)

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.

ISWTracer(*, z_max=6.0, n_chi=1024)

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.

NumberCountsTracer(*, dndz, bias=None, mag_bias=None, has_rsd, n_samples=256)

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.

WeakLensingTracer(*, dndz, has_shear=True, ia_bias=None, use_A_ia=True, n_samples=256)

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.

age_of_universe(a)

Age of the Universe at some scale factor, in \(\rm Gyr\).

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or (na,) array_like) – Scale factor(s) normalized to 1 today.

Returns:

t_age – Age of the Universe at a. nan if a is out of bounds of the spline parametets stored in cosmo.

Return type:

float or (na,) ndarray

angular_cl(tracer1, tracer2, ell, *, p_of_k_a='delta_matter:delta_matter', l_limber=-1, limber_max_error=0.01, limber_integration_method='qag_quad', non_limber_integration_method='FKEM', p_of_k_a_lin='delta_matter:delta_matter', return_meta=False)

Calculate the angular (cross-)power spectrum for a pair of tracers.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • tracer1 (Tracer) – a Tracer object, of any kind.

  • tracer2 (Tracer) – a second Tracer object.

  • ell (float or array) – Angular multipole(s) at which to evaluate the angular power spectrum.

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to project. If a string, it must correspond to one of the non-linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

  • l_limber (int, float or 'auto') – Angular wavenumber beyond which Limber’s approximation will be used. Defaults to -1. If ‘auto’, then the non-limber integrator will be used to compute the right transition point given the value of limber_max_error.

  • limber_max_error (float) – Maximum fractional error for Limber integration.

  • limber_integration_method (string) – integration method to be used for the Limber integrals. Possibilities: ‘qag_quad’ (GSL’s qag method backed up by quad when it fails) and ‘spline’ (the integrand is splined and then integrated numerically).

  • non_limber_integration_method (string) – integration method to be used for the non-Limber integrals. Currently the only method implemented is 'FKEM' (see the N5K paper for details).

  • p_of_k_a_lin (Pk2D, str or None) – 3D linear Power spectrum to project, for special use in PT calculations using the FKEM non-limber integration technique. If a string, it must correspond to one of the linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

  • return_meta (bool) – if True, also return a dictionary with various metadata about the calculation, such as l_limber as calculated by the non-limber integrator.

Returns:

Angular (cross-)power spectrum values, \(C_\ell\), for the pair of tracers, as a function of \(\ell\).

Return type:

float or array

angular_cl_cov_SSC(tracer1, tracer2, *, ell, t_of_kk_a, tracer3=None, tracer4=None, ell2=None, sigma2_B=None, fsky=1.0, integration_method='qag_quad')

Calculate the super-sample contribution to the connected non-Gaussian covariance for a pair of power spectra \(C_{\ell_1}^{ab}\) and \(C_{\ell_2}^{cd}\), between two pairs of tracers (\((a,b)\) and \((c,d)\)).

Specifically, it computes:

\[{\rm Cov}_{\rm SSC}(\ell_1,\ell_2)= \int \frac{d\chi}{\chi^4} \tilde{\Delta}^a_{\ell_1}(\chi) \tilde{\Delta}^b_{\ell_1}(\chi) \tilde{\Delta}^c_{\ell_2}(\chi) \tilde{\Delta}^d_{\ell_2}(\chi)\, \bar{T}_{abcd}\left[\frac{\ell_1+1/2}{\chi}, \frac{\ell_2+1/2}{\chi}, a(\chi)\right]\]

where \(\Delta^x_\ell(\chi)\) is the transfer function for tracer \(x\) (see Eq. 39 in the CCL note), and \(\bar{T}_{abcd}(k_1,k_2,a)\) is the isotropized connected trispectrum of the four tracers (see the documentation of the Tk3D class for details).

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • tracer1 (Tracer) – a Tracer object.

  • tracer2 (Tracer) – a second Tracer object.

  • ell (float or array) – Angular wavenumber(s) at which to evaluate the first dimension of the angular power spectrum covariance.

  • t_of_kk_a (Tk3D) – 3D connected trispectrum.

  • tracer3 (Tracer) – a Tracer object. If None, tracer1 will be used instead.

  • tracer4 (Tracer) – a Tracer object. If None, tracer2 will be used instead.

  • ell2 (float or array) – Angular wavenumber(s) at which to evaluate the second dimension of the angular power spectrum covariance. If None, ell will be used instead.

  • sigma2_B (tuple or None) – A tuple of arrays (a, sigma2_B(a)) containing the variance of the projected matter overdensity over the footprint as a function of the scale factor. If None, a compact circular footprint will be assumed covering a sky fraction fsky.

  • fsky (float) – sky fraction.

  • integration_method (str) – integration method to be used for the Limber integrals. Possibilities: 'qag_quad' (GSL’s qag method backed up by quad when it fails) and 'spline' (the integrand is splined and then integrated analytically).

Returns:

2D array containing the super-sample Angular power spectrum covariance \({\rm Cov}_{\rm SSC}(\ell_1,\ell_2)\), for the four input tracers, as a function of \(\ell_1\) and \(\ell_2\). The ordering is such that out[i2, i1] = Cov(ell2[i2], ell[i1]).

Return type:

(float or array)

angular_cl_cov_cNG(tracer1, tracer2, *, ell, t_of_kk_a, tracer3=None, tracer4=None, ell2=None, fsky=1.0, integration_method='qag_quad')

Calculate the connected non-Gaussian covariance for a pair of angular power spectra \(C_{\ell_1}^{ab}\) and \(C_{\ell_2}^{cd}\), between two pairs of tracers (\((a,b)\) and \((c,d)\)).

Specifically, it computes:

\[{\rm Cov}_{\rm cNG}(\ell_1,\ell_2)=\frac{1}{4\pi f_{\rm sky}} \int \frac{d\chi}{\chi^6} \tilde{\Delta}^a_{\ell_1}(\chi) \tilde{\Delta}^b_{\ell_1}(\chi) \tilde{\Delta}^c_{\ell_2}(\chi) \tilde{\Delta}^d_{\ell_2}(\chi)\, \bar{T}_{abcd}\left[\frac{\ell_1+1/2}{\chi}, \frac{\ell_2+1/2}{\chi}, a(\chi)\right]\]

where \(\Delta^x_\ell(\chi)\) is the transfer function for tracer \(x\) (see Eq. 39 in the CCL note), and \(\bar{T}_{abcd}(k_1,k_2,a)\) is the isotropized connected trispectrum of the four tracers (see the documentation of the Tk3D class for details).

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • tracer1 (Tracer) – a Tracer object.

  • tracer2 (Tracer) – a second Tracer object.

  • ell (float or array) – Angular wavenumber(s) at which to evaluate the first dimension of the angular power spectrum covariance.

  • t_of_kk_a (Tk3D) – 3D connected trispectrum.

  • tracer3 (Tracer) – a Tracer object. If None, tracer1 will be used instead.

  • tracer4 (Tracer) – a Tracer object. If None, tracer2 will be used instead.

  • ell2 (float or array) – Angular wavenumber(s) at which to evaluate the second dimension of the angular power spectrum covariance. If None, ell will be used instead.

  • fsky (float) – sky fraction.

  • integration_method (str) – integration method to be used for the Limber integrals. Possibilities: 'qag_quad' (GSL’s qag method backed up by quad when it fails) and 'spline' (the integrand is splined and then integrated analytically).

Returns:

2D array containing the connected non-Gaussian Angular power spectrum covariance \({\rm Cov}_{\rm cNG}(\ell_1,\ell_2)\), for the four input tracers, as a function of \(\ell_1\) and \(\ell_2\). The ordering is such that out[i2, i1] = Cov(ell2[i2], ell[i1]).

Return type:

(float or array)

angular_diameter_distance(a1, a2=None)

Angular diameter distance.

The angular diameter distance in Mpc from scale factor a1 to scale factor a2. If a2 is not provided, it is assumed that the distance will be calculated between 1 and a1.

Note

a2 has to be smaller than a1 (i.e. a source at a2 is behind one at a1). You can compute the distance between a single lens at a1 and multiple sources at a2 by passing a scalar a1.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a1 (float or array) – Scale factor(s), normalized to 1 today.

  • a2 (float or array) – Scale factor(s), normalized to 1 today, optional.

Returns:

angular diameter distance; Mpc.

Return type:

(float or array)

comoving_angular_distance(a)

Comoving angular distance.

Note

This quantity is otherwise known as the transverse comoving distance, and is NOT angular diameter distance or angular separation. The comoving angular distance is defined such that the comoving distance between two objects at a fixed scale factor separated by an angle \(\theta\) is \(\theta r_{A}(a)\) where \(r_{A}(a)\) is the comoving angular distance.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Comoving angular distance; Mpc.

Return type:

(float or array)

comoving_radial_distance(a)

Comoving radial distance.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Comoving radial distance; Mpc.

Return type:

(float or array)

comoving_volume(a, *, solid_angle=12.566370614359172)

Comoving volume, in \(\rm Mpc^3\).

\[V_{\rm C} = \int_{\Omega} \mathrm{{d}}\Omega \int_z \mathrm{d}z D_{\rm H} \frac{(1+z)^2 D_{\mathrm{A}}^2}{E(z)}\]

See Eq. 29 in Hogg 2000.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or (na,) array_like) – Scale factor(s) normalized to 1 today.

  • solid_angle (float) – Solid angle subtended in the sky for which the comoving volume is calculated.

Returns:

V_C – Comoving volume at a.

Return type:

float or (na,) ndarray

See also

comoving_volume_element

comoving volume element

comoving_volume_element(a)

Comoving volume element in \(\rm Mpc^3 \, sr^{-1}\).

\[\frac{\mathrm{d}V}{\mathrm{d}a \, \mathrm{d} \Omega}\]
Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or (na,) array_like) – Scale factor(s) normalized to 1 today.

Returns:

dV – Comoving volume per unit scale factor per unit solid angle.

Return type:

float or (na,) numpy.ndarray

See also

comoving_volume

integral of the comoving volume element

compute_distances()

Compute the distance splines.

convert_concentration(*, c_old, Delta_old, Delta_new)

Computes the concentration parameter for a different mass definition. This is done assuming an NFW profile. The output concentration c_new is found by solving the equation:

\[f(c_{\rm old}) \Delta_{\rm old} = f(c_{\rm new}) \Delta_{\rm new}\]

where

\[f(x) = \frac{x^3}{\log(1+x) - x/(1+x)}.\]
Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • c_old (float or array) – concentration to translate from.

  • Delta_old (float) – Delta parameter associated to the input concentration. See description of the MassDef class.

  • Delta_new (float) – Delta parameter associated to the output concentration.

Returns:

concentration parameter for the new mass definition.

Return type:

(float or array)

correlation(*, ell, C_ell, theta, type='NN', method='fftlog')

Compute the angular correlation function.

\[\xi^{ab}_\pm(\theta) = \sum_\ell\frac{2\ell+1}{4\pi}\,(\pm1)^{s_b}\, C^{ab\pm}_\ell\,d^\ell_{s_a,\pm s_b}(\theta)\]

where \(\theta\) is the angle between the two fields \(a\) and \(b\) with spins \(s_a\) and \(s_b\) after alignement of their tangential coordinate. \(d^\ell_{mm'}\) are the Wigner-d matrices and we have defined the power spectra

\[C^{ab\pm}_\ell \equiv (C^{a_Eb_E}_\ell \pm C^{a_Bb_B}_\ell)+i (C^{a_Bb_E}_\ell \mp C^{a_Eb_B}_\ell),\]

which reduces to the \(EE\) power spectrum when all \(B\)-modes are 0.

The different spin combinations are:

  • \(s_a=s_b=0\) e.g. galaxy-galaxy, galaxy-\(\kappa\) and \(\kappa\)-\(\kappa\)

  • \(s_a=2\), \(s_b=0\) e.g. galaxy-shear, and \(\kappa\)-shear

  • \(s_a=s_b=2\) e.g. shear-shear.

Note

For scales smaller than \(\sim 0.1^{\circ}\), the input power spectrum should be sampled to sufficienly high \(\ell\) to ensure the Hankel transform is well-behaved. The following spline parameters, related to FFTLog-sampling may also be modified for accuracy:

  • ccl.spline_params.ELL_MIN_CORR

  • ccl.spline_params.ELL_MAX_CORR

  • ccl.spline_params.N_ELL_CORR.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • ell (array) – Multipoles corresponding to the input angular power spectrum.

  • C_ell (array) – Input angular power spectrum.

  • theta (float or array) – Angular separation(s) at which to calculate the angular correlation function (in degrees).

  • type (str) – Type of correlation function. Choices: 'NN' (0x0), 'NG' (0x2), 'GG+' (2x2, \(\xi_+\)), 'GG-' (2x2, \(\xi_-\)), where numbers refer to the spins of the two quantities being cross-correlated (see Section 2.4.2 of the CCL paper). The naming system roughly follows the nomenclature used in TreeCorr.

  • method (str) – Method to compute the correlation function. Choices: 'Bessel' (direct integration over Bessel function), 'FFTLog' (fast integration with FFTLog), 'Legendre' (brute-force sum over Legendre polynomials).

Returns:

Value(s) of the correlation function at the input angular separations.

Return type:

(float or array)

correlation_3d(*, r, a, p_of_k_a='delta_matter:delta_matter')

Compute the 3D correlation function:

\[\xi(r)\equiv\frac{1}{2\pi^2}\int dk\,k^2\,P(k)\,j_0(kr).\]
Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • r (float or array) – distance(s) at which to calculate the 3D correlation function (in Mpc).

  • a (float) – scale factor.

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to integrate. If a string, it must correspond to one of the non-linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

Returns:

Value(s) of the correlation function at the input distance(s).

correlation_3dRsd(*, r, a, mu, beta, p_of_k_a='delta_matter:delta_matter', use_spline=True)

Compute the 3D correlation function with linear RSDs using multipoles:

\[\xi(r,\mu) = \sum_{\ell\in\{0,2,4\}}\xi_\ell(r)\,P_\ell(\mu)\]

where \(P_\ell(\mu)\) are the Legendre polynomials, and \(\xi_\ell(r)\) are the correlation function multipoles.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • r (float or array) – distance(s) at which to calculate the 3D correlation function (in Mpc).

  • a (float) – scale factor.

  • mu (float) – cosine of the angle at which to calculate the 3D correlation function.

  • beta (float) – growth rate divided by galaxy bias.

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to integrate. If a string, it must correspond to one of the non-linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

  • use_spline (bool) – switch that determines whether the RSD correlation function is calculated using global splines of multipoles.

Returns:

Value(s) of the correlation function at the input distance(s) & angle.

correlation_3dRsd_avgmu(*, r, a, beta, p_of_k_a='delta_matter:delta_matter')

Compute the 3D correlation function averaged over angles with RSDs. Equivalent to calling correlation_multipole() with ell=0.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • r (float or array) – distance(s) at which to calculate the 3D correlation function (in Mpc).

  • a (float) – scale factor.

  • beta (float) – growth rate divided by galaxy bias.

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to integrate. If a string, it must correspond to one of the non-linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

Returns:

Value(s) of the correlation function at the input distance(s) & angle.

correlation_multipole(*, r, a, beta, ell, p_of_k_a='delta_matter:delta_matter')

Compute the correlation function multipoles:

\[\xi_\ell(r)\equiv\frac{i^\ell}{2\pi^2}\int dk\,k^2\,P(k)\,j_\ell(kr).\]
Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • r (float or array) – distance(s) at which to calculate the 3D correlation function (in Mpc).

  • a (float) – scale factor.

  • beta (float) – growth rate divided by galaxy bias.

  • ell (int) – the desired multipole

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to integrate. If a string, it must correspond to one of the non-linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

Returns:

Value(s) of the correlation function at the input distance(s).

correlation_pi_sigma(*, pi, sigma, a, beta, use_spline=True, p_of_k_a='delta_matter:delta_matter')

Compute the 3D correlation in \((\pi,\sigma)\) space. This is just

\[\xi(\pi,\sigma) = \xi(r=\sqrt{\pi^2+\sigma^2},\mu=\pi/r).\]
Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • pi (float) – distance times cosine of the angle (in Mpc).

  • sigma (float or array) – distance(s) times sine of the angle (in Mpc).

  • a (float) – scale factor.

  • beta (float) – growth rate divided by galaxy bias.

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to integrate. If a string, it must correspond to one of the non-linear power spectra stored in cosmo (e.g. ‘delta_matter:delta_matter’).

  • use_spline (bool) – switch that determines whether the RSD correlation function is calculated using global splines of multipoles.

Returns:

Value(s) of the correlation function at the input pi and sigma.

distance_modulus(a)

Distance Modulus, defined as

\[\mu = 5\,\log_{10}(d_L/10\,{\rm pc})\]

where \(d_L\) is the luminosity distance.

Note

The distance modulus can be used to convert between apparent and absolute magnitudes via \(m = M + \mu\), where \(m\) is the apparent magnitude and \(M\) is the absolute magnitude.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Distance modulus at a.

Return type:

(float or array)

get_camb_pk_lin(*, nonlin=False)

Run CAMB and return the linear power spectrum.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters. The cosmological parameters with which to run CAMB.

  • nonlin (bool) – Whether to compute and return the non-linear power spectrum as well.

Returns:

Power spectrum object. The linear power spectrum. If nonlin=True, returns a tuple (pk_lin, pk_nonlin).

Return type:

Pk2D

get_class_pk_lin()

Run CLASS and return the linear power spectrum.

Parameters:

cosmo (Cosmology) – Cosmological parameters. The cosmological parameters with which to run CLASS.

Returns:

Power spectrum object. The linear power spectrum.

Return type:

Pk2D

get_delta_c(a, kind='EdS')

Returns the linear collapse threshold.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • a (float or array) – scale factor.

  • kind (str) –

    prescription to use. Should be one of

    • ’EdS’: the SC prediction in Einstein de-Sitter, \(\delta_c=(3/20)(12\pi)^{2/3}\).

    • ’EdS_approx’: a common approximation to the EdS result \(\delta_c=1.686\).

    • ’NakamuraSuto97’: the prescription from Nakamura & Suto 1997.

    • ’Mead16’: the prescription from Mead et al. 2016.

Returns:

linear collapse threshold.

Return type:

(float or array)

get_density_kernel(*, dndz)

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.

get_isitgr_pk_lin()

Run ISiTGR-CAMB and return the linear power spectrum.

Parameters:

cosmo (Cosmology) – Cosmological parameters. The cosmological parameters with which to run ISiTGR-CAMB.

Returns:

Power spectrum object. The linear power spectrum.

Return type:

Pk2D

get_kappa_kernel(*, z_source, n_samples=100)

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.

get_lensing_kernel(*, dndz, mag_bias=None, n_chi=None)

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.

growth_factor(a)

Growth factor.

Warning

CCL is not able to compute the scale-dependent growth factor for cosmologies with massive neutrinos.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Growth factor, normalized to unity today.

Return type:

(float or array)

growth_factor_unnorm(a)

Unnormalized growth factor.

Warning

CCL is not able to compute the scale-dependent growth factor for cosmologies with massive neutrinos.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Unnormalized growth factor, normalized to the scale factor at early times.

Return type:

(float or array)

growth_rate(a)

Growth rate defined as the logarithmic derivative of the growth factor, \(f\equiv d\log D/d\log a\).

Warning

CCL is not able to compute the scale-dependent growth rate for cosmologies with massive neutrinos.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Growth rate.

Return type:

(float or array)

h_over_h0(a)

Ratio of Hubble constant at a over Hubble constant today.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

H(a)/H0.

Return type:

(float or array)

halomod_Pk2D(hmc, prof, *, prof2=None, prof_2pt=None, p_of_k_a=None, get_1h=True, get_2h=True, lk_arr=None, a_arr=None, extrap_order_lok=1, extrap_order_hik=2, smooth_transition=None, suppress_1h=None, extrap_pk=False)

Returns a Pk2D object containing the halo-model power spectrum for two quantities defined by their respective halo profiles. See halomod_power_spectrum() for more details about the actual calculation.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile.

  • prof2 (HaloProfile) – a second halo profile. If None, prof will be used as prof2.

  • prof_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of the two profiles being correlated. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • get_1h (bool) – if False, the 1-halo term (i.e. the first term in the first equation above) won’t be computed.

  • get_2h (bool) – if False, the 2-halo term (i.e. the second term in the first equation above) won’t be computed.

  • a_arr (array) – an array holding values of the scale factor at which the halo model power spectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of Mpc^-1) at which the halo model power spectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Pk2D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Pk2D.

  • smooth_transition (callable or None) – Modify the halo model 1-halo/2-halo transition region via a time-dependent function \(\alpha(a)\), defined as in HMCODE-2020: \(P(k,a)= (P_{1h}^{\alpha(a)}(k)+P_{2h}^{\alpha(a)}(k))^{1/\alpha}\). If None the extra factor is just 1.

  • suppress_1h (callable or None) –

    Suppress the 1-halo large scale contribution by a time- and scale-dependent function \(k_*(a)\), defined as in HMCODE-2020: \(1/[1+(k_*(a)/k)^4]\). If None the standard 1-halo term is returned with no damping.

  • extrap_pk (bool) – Whether to extrapolate p_of_k_a in case a is out of its support. If `False`, and the queried values are out of bounds, an error is raised.

Returns:

halo model power spectrum.

Return type:

Pk2D

halomod_Tk3D_1h(hmc, prof, *, prof2=None, prof3=None, prof4=None, prof12_2pt=None, prof34_2pt=None, lk_arr=None, a_arr=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False)

Returns a Tk3D object containing the 1-halo trispectrum for four quantities defined by their respective halo profiles. See halomod_trispectrum_1h() for more details about the actual calculation.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile. If None, prof2 will be used as prof4.

  • prof12_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof2. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof34_2pt (Profile2pt) – same as prof12_2pt for prof3 and prof4.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of \({\rm Mpc}^{-1}\)) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

Returns:

1-halo trispectrum.

Return type:

Tk3D

halomod_Tk3D_2h(hmc, prof, prof2=None, prof3=None, prof4=None, prof12_2pt=None, prof13_2pt=None, prof14_2pt=None, prof24_2pt=None, prof32_2pt=None, prof34_2pt=None, p_of_k_a=None, lk_arr=None, a_arr=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False)

Returns a Tk3D object containing the 2-halo trispectrum for four quantities defined by their respective halo profiles. See halomod_trispectrum_1h() for more details about the actual calculation.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • prof12_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof2. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof13_2pt (Profile2pt) – same as prof12_2pt for prof and prof3.

  • prof14_2pt (Profile2pt) – same as prof14_2pt for prof and prof4.

  • prof24_2pt (Profile2pt) – same as prof14_2pt for prof2 and prof4.

  • prof32_2pt (Profile2pt) – same as prof14_2pt for prof3 and prof2.

  • prof34_2pt (Profile2pt) – same as prof34_2pt for prof3 and prof4.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of Mpc^-1) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

Returns:

2-halo trispectrum.

Return type:

Tk3D

halomod_Tk3D_3h(hmc, prof, prof2=None, prof3=None, prof4=None, prof13_2pt=None, prof14_2pt=None, prof24_2pt=None, prof32_2pt=None, lk_arr=None, a_arr=None, p_of_k_a=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False)

Returns a Tk3D object containing the 3-halo trispectrum for four quantities defined by their respective halo profiles. See halomod_trispectrum_3h() for more details about the actual calculation.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • prof13_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof3. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof14_2pt (Profile2pt) – same as prof14_2pt for prof and prof4.

  • prof24_2pt (Profile2pt) – same as prof14_2pt for prof2 and prof4.

  • prof32_2pt (Profile2pt) – same as prof14_2pt for prof3 and prof2.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of Mpc^-1) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

Returns:

3-halo trispectrum.

Return type:

Tk3D

halomod_Tk3D_4h(hmc, prof, prof2=None, prof3=None, prof4=None, lk_arr=None, a_arr=None, p_of_k_a=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False)

Returns a Tk3D object containing the 3-halo trispectrum for four quantities defined by their respective halo profiles. See halomod_trispectrum_4h() for more details about the actual calculation.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of Mpc^-1) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

Returns:

4-halo trispectrum.

Return type:

Tk3D

halomod_Tk3D_SSC(hmc, prof, *, prof2=None, prof3=None, prof4=None, prof12_2pt=None, prof34_2pt=None, p_of_k_a=None, lk_arr=None, a_arr=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False, extrap_pk=False)

Returns a Tk3D object containing the super-sample covariance trispectrum, given by the tensor product of the power spectrum responses associated with the two pairs of quantities being correlated. Each response is calculated as:

\[\frac{\partial P_{u,v}(k)}{\partial\delta_L} = \left(\frac{68}{21}-\frac{d\log k^3P_L(k)}{d\log k}\right) P_L(k)I^1_1(k,|u)I^1_1(k,|v)+I^1_2(k|u,v) - (b_{u} + b_{v}) P_{u,v}(k)\]

where the \(I^a_b\) are defined in the documentation of I_1_1() and I_1_2() and \(b_{u}\) and \(b_{v}\) are the linear halo biases for quantities \(u\) and \(v\), respectively (zero if the profiles are not number counts).

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above). If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above). If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above). If None, prof2 will be used as prof4.

  • prof12_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof2. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof34_2pt (Profile2pt) – same as prof12_2pt for prof3 and prof4. If None, prof12_2pt will be used.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of \({\rm Mpc}^{-1}\)) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

  • extrap_pk (bool) – Whether to extrapolate p_of_k_a in case a is out of its support. If False, and the queried values are out of bounds, an error is raised. The default is False.

Returns:

SSC effective trispectrum.

Return type:

Tk3D

halomod_Tk3D_SSC_linear_bias(hmc, *, prof, bias1=1, bias2=1, bias3=1, bias4=1, is_number_counts1=False, is_number_counts2=False, is_number_counts3=False, is_number_counts4=False, p_of_k_a=None, lk_arr=None, a_arr=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False, extrap_pk=False)

Returns a Tk3D object containing the super-sample covariance trispectrum, given by the tensor product of the power spectrum responses associated with the two pairs of quantities being correlated. Each response is calculated as:

\[\frac{\partial P_{u,v}(k)}{\partial\delta_L} = b_u b_v \left( \left(\frac{68}{21}-\frac{d\log k^3P_L(k)}{d\log k}\right) P_L(k)+I^1_2(k|u,v)\right) - (b_{u} + b_{v}) P_{u,v}(k)\]

where the \(I^1_2\) is defined in the documentation I_1_2() and \(b_{u}\) and \(b_{v}\) are the linear halo biases for quantities \(u\) and \(v\), respectively. The second term is only included if the corresponding profiles do not represent number counts.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – a halo profile representing the matter overdensity.

  • bias1 (float or array) – linear galaxy bias for quantity 1. If an array, it has to have the shape of a_arr.

  • bias2 (float or array) – linear galaxy bias for quantity 2.

  • bias3 (float or array) – linear galaxy bias for quantity 3.

  • bias4 (float or array) – linear galaxy bias for quantity 4.

  • is_number_counts1 (bool) – If True, quantity 1 will be considered number counts and the clustering counter terms computed.

  • is_number_counts2 (bool) – as is_number_counts1 but for quantity 2.

  • is_number_counts3 (bool) – as is_number_counts1 but for quantity 3.

  • is_number_counts4 (bool) – as is_number_counts1 but for quantity 4.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of \({\rm Mpc}^{-1}\)) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

  • extrap_pk (bool) – Whether to extrapolate p_of_k_a in case a is out of its support. If False, and the queried values are out of bounds, an error is raised. The default is False.

Returns:

SSC effective trispectrum.

Return type:

Tk3D

halomod_Tk3D_cNG(hmc, prof, prof2=None, prof3=None, prof4=None, prof12_2pt=None, prof13_2pt=None, prof14_2pt=None, prof24_2pt=None, prof32_2pt=None, prof34_2pt=None, p_of_k_a=None, lk_arr=None, a_arr=None, extrap_order_lok=1, extrap_order_hik=1, use_log=False)

Returns a Tk3D object containing the non-Gaussian covariance trispectrum for four quantities defined by their respective halo profiles. This is the sum of the trispectrum terms 1h + 2h + 3h + 4h.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • prof12_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof2. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof13_2pt (Profile2pt) – same as prof12_2pt for prof and prof3.

  • prof14_2pt (Profile2pt) – same as prof12_2pt for prof and prof4.

  • prof24_2pt (Profile2pt) – same as prof12_2pt for prof2 and prof4.

  • prof32_2pt (Profile2pt) – same as prof12_2pt for prof3 and prof2.

  • prof34_2pt (Profile2pt) – same as prof12_2pt for prof3 and prof4.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • a_arr (array) – an array holding values of the scale factor at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • lk_arr (array) – an array holding values of the natural logarithm of the wavenumber (in units of Mpc^-1) at which the trispectrum should be calculated for interpolation. If None, the internal values used by cosmo will be used.

  • extrap_order_lok (int) – extrapolation order to be used on k-values below the minimum of the splines. See Tk3D.

  • extrap_order_hik (int) – extrapolation order to be used on k-values above the maximum of the splines. See Tk3D.

  • use_log (bool) – if True, the trispectrum will be interpolated in log-space (unless negative or zero values are found).

Returns:

2-halo trispectrum.

Return type:

Tk3D

halomod_bias_1pt(hmc, k, a, prof)

Returns the mass-and-bias-weighted mean halo profile.

\[I^1_1(k,a|u) = \int dM\,n(M,a)\,b(M,a)\, \langle u(k,a|M)\rangle,\]

where \(n(M,a)\) is the halo mass function, \(b(M,a)\) is the halo bias, and \(\langle u(k,a|M)\rangle\) is the halo profile as a function of scale, scale factor and halo mass.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • k (float or array) – comoving wavenumber in \({\rm Mpc}^{-1}\).

  • a (float or array) – scale factor.

  • prof (HaloProfile) – halo profile.

Returns:

integral values evaluated at each combination of k and a. The shape of the output will be (N_a, N_k) where N_k and N_a are the sizes of k and a respectively. If k or a are scalars, the corresponding dimension will be squeezed out on output.

Return type:

(float or array)

halomod_mean_profile_1pt(hmc, k, a, prof)

Returns the mass-weighted mean halo profile.

\[I^0_1(k,a|u) = \int dM\,n(M,a)\,\langle u(k,a|M)\rangle,\]

where \(n(M,a)\) is the halo mass function, and \(\langle u(k,a|M)\rangle\) is the halo profile as a function of scale, scale factor and halo mass.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • k (float or array) – comoving wavenumber in \({\rm Mpc}^{-1}\).

  • a (float or array) – scale factor.

  • prof (HaloProfile) – halo profile.

Returns:

integral values evaluated at each combination of k and a. The shape of the output will be (N_a, N_k) where N_k and N_a are the sizes of k and a respectively. If k or a are scalars, the corresponding dimension will be squeezed out on output.

Return type:

(float or array)

halomod_power_spectrum(hmc, k, a, prof, *, prof2=None, prof_2pt=None, p_of_k_a=None, get_1h=True, get_2h=True, smooth_transition=None, suppress_1h=None, extrap_pk=False)

Computes the halo model power spectrum for two quantities defined by their respective halo profiles. The halo model power spectrum for two profiles \(u\) and \(v\) is:

\[P_{u,v}(k,a) = I^0_2(k,a|u,v) + I^1_1(k,a|u)\,I^1_1(k,a|v)\,P_{\rm lin}(k,a)\]

where \(P_{\rm lin}(k,a)\) is the linear matter power spectrum, \(I^1_1\) is defined in the documentation of I_1_1(), and \(I^0_2\) is defined in the documentation of I_0_2().

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • k (float or array) – comoving wavenumber in Mpc^-1.

  • a (float or array) – scale factor.

  • prof (HaloProfile) – halo profile.

  • prof2 (HaloProfile) – a second halo profile. If None, prof will be used as prof2.

  • prof_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of the two profiles being correlated. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

  • get_1h (bool) – if False, the 1-halo term (i.e. the first term in the first equation above) won’t be computed.

  • get_2h (bool) – if False, the 2-halo term (i.e. the second term in the first equation above) won’t be computed.

  • smooth_transition (callable or None) –

    Modify the halo model 1-halo/2-halo transition region via a time-dependent function \(\alpha(a)\), defined as in HMCODE-2020: \(P(k,a)= (P_{1h}^{\alpha(a)}(k)+P_{2h}^{\alpha(a)}(k))^{1/\alpha}\). If None the extra factor is just 1.

  • suppress_1h (callable or None) –

    Suppress the 1-halo large scale contribution by a time- and scale-dependent function \(k_*(a)\), defined as in HMCODE-2020: \(1/[1+(k_*(a)/k)^4]\). If None the standard 1-halo term is returned with no damping.

  • extrap_pk (bool) – Whether to extrapolate p_of_k_a in case a is out of its support. If `False`, and the queried values are out of bounds, an error is raised.

Returns:

integral values evaluated at each combination of k and a. The shape of the output will be (N_a, N_k) where N_k and N_a are the sizes of k and a respectively. If k or a are scalars, the corresponding dimension will be squeezed out on output.

Return type:

(float or array)

halomod_trispectrum_1h(hmc, k, a, prof, *, prof2=None, prof3=None, prof4=None, prof12_2pt=None, prof34_2pt=None)

Computes the halo model 1-halo trispectrum for four different quantities defined by their respective halo profiles. The 1-halo trispectrum for four profiles \(u_{1,2}\), \(v_{1,2}\) is calculated as:

\[T_{u_1,u_2;v_1,v_2}(k_u,k_v,a) = I^0_{2,2}(k_u,k_v,a|u_{1,2},v_{1,2})\]

where \(I^0_{2,2}\) is defined in the documentation of I_0_22().

Note

This approximation assumes that the 4-point profile cumulant is the same as the product of two 2-point cumulants. We may relax this assumption in future versions of CCL.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • k (float or array) – comoving wavenumber in \({\rm Mpc}^{-1}\).

  • a (float or array) – scale factor.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above).

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above). If None, prof will be used as prof2.

  • prof12_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof2. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • prof34_2pt (Profile2pt) – same as prof12_2pt for prof3 and prof4.

Returns:

1-halo trispectrum evaluated at each combination of k and a. The shape of the output will be (N_a, N_k, N_k) where N_k and N_a are the sizes of k and a respectively. The ordering is such that output[ia, ik2, ik1] = T(k[ik1], k[ik2], a[ia]) If k or a are scalars, the corresponding dimension will be squeezed out on output.

Return type:

(float or array)

halomod_trispectrum_3h(hmc, k, a, prof, *, prof2=None, prof3=None, prof4=None, prof13_2pt=None, prof14_2pt=None, prof24_2pt=None, prof32_2pt=None, p_of_k_a=None)

Computes the isotropized halo model 3-halo trispectrum for four profiles \(u_{1,2}\), \(v_{1,2}\) as

\[\bar{T}^{3h}(k_1, k_2, a) = \int \frac{d\varphi_1}{2\pi} \int \frac{d\varphi_2}{2\pi} T^{2h}_{22}({\bf k_1},-{\bf k_1},{\bf k_2},-{\bf k_2}),\]

with

\[T^{3h}{u_1,u_2;v_1,v_2}(k_u,k_v,a) = B^{PT}({f k_{u_1}}, {f k_{u_2}}, {f k_{v_1}} + {f k_{v_2}}) \, I^1_1(k_{u_1} | u) I^1_1(k_{u_2} | u) I^1_2(k_{v_1}, k_{v_2}|v}) \, + 5 perm\]

where \(I^1_1\) and \(I^1_2\) are defined in the documentation of \(~HMCalculator.I_1_1\) and \(~HMCalculator.I_1_2\), respectively; and \(B^{PT}\) can be found in Eq. 30 of arXiv:1302.6994.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • k (float or array_like) – comoving wavenumber in Mpc^-1.

  • a (float or array_like) – scale factor.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • prof13_2pt (Profile2pt) – a profile covariance object returning the the two-point moment of prof and prof3. If None, the default second moment will be used, corresponding to the products of the means of both profiles.

  • prof14_2pt (Profile2pt) – same as prof14_2pt for prof and prof4.

  • prof24_2pt (Profile2pt) – same as prof14_2pt for prof2 and prof4.

  • prof32_2pt (Profile2pt) – same as prof14_2pt for prof3 and prof2.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

Returns:

integral values evaluated at each combination of k and a. The shape of the output will be (N_a, N_k, N_k) where N_k and N_a are the sizes of k and a respectively. The ordering is such that output[ia, ik2, ik1] = T(k[ik1], k[ik2], a[ia]) If k or a are scalars, the corresponding dimension will be squeezed out on output.

Return type:

float or array_like

halomod_trispectrum_4h(hmc, k, a, prof, prof2=None, prof3=None, prof4=None, p_of_k_a=None)

Computes the isotropized halo model 4-halo trispectrum for four profiles \(u_{1,2}\), \(v_{1,2}\) as

\[\bar{T}^{4h}(k_1, k_2, a) = \int \frac{d\varphi_1}{2\pi} \int \frac{d\varphi_2}{2\pi} T^{4h}({\bf k_1},-{\bf k_1},{\bf k_2},-{\bf k_2}),\]

with

\[T^{4h}{u_1,u_2;v_1,v_2}(k_u,k_v,a) = T^{PT}({f k_{u_1}}, {f k_{u_2}}, {f k_{v_1}}, {f k_{v_2}}) \, I^1_1(k_{u_1} | u) I^1_1(k_{u_2} | u) I^1_1(k_{v_1} | v) \, I^1_1(k_{v_2} | v) \,\]

where \(I^1_1\) is defined in the documentation of \(~HMCalculator.I_1_1\) and \(P^{PT}\) can be found in Eq. 30 of arXiv:1302.6994.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • hmc (HMCalculator) – a halo model calculator.

  • k (float or array_like) – comoving wavenumber in Mpc^-1.

  • a (float or array_like) – scale factor.

  • prof (HaloProfile) – halo profile (corresponding to \(u_1\) above.

  • prof2 (HaloProfile) – halo profile (corresponding to \(u_2\) above. If None, prof will be used as prof2.

  • prof3 (HaloProfile) – halo profile (corresponding to \(v_1\) above. If None, prof will be used as prof3.

  • prof4 (HaloProfile) – halo profile (corresponding to \(v_2\) above. If None, prof2 will be used as prof4.

  • p_of_k_a (Pk2D) – a Pk2D object to be used as the linear matter power spectrum. If None, the power spectrum stored within cosmo will be used.

Returns:

integral values evaluated at each combination of k and a. The shape of the output will be (N_a, N_k, N_k) where N_k and N_a are the sizes of k and a respectively. The ordering is such that output[ia, ik2, ik1] = T(k[ik1], k[ik2], a[ia]) If k or a are scalars, the corresponding dimension will be squeezed out on output.

Return type:

float or array_like

hubble_distance(a)

Hubble distance in \(\rm Mpc\).

\[D_{\rm H} = \frac{cz}{H_0}\]
Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or (na,) array_like) – Scale factor(s) normalized to 1 today.

Returns:

D_H – Hubble distance.

Return type:

float or (na,) numpy.ndarray

kNL(a, *, p_of_k_a='delta_matter:delta_matter')

Non-linear scale \(k_{\rm NL}\). Calculated based on Lagrangian perturbation theory as the inverse of the rms of the displacement field, i.e.:

\[k_{\rm NL}(z) = \left[\frac{1}{6\pi^2} \int dk\,P_L(k,z)\right]^{-1/2}.\]
Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

  • p_of_k_a (Pk2D, str or None) – power spectrum to integrate. If a string, it must correspond to one of the linear power spectra stored in cosmo (e.g. 'delta_matter:delta_matter').

Returns:

\(k_{\rm NL}\).

Return type:

float or array

linear_matter_power(k, a)

The linear matter power spectrum

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • k (float or array) – Wavenumber; \({\rm Mpc}^{-1}\).

  • a (float or array) – Scale factor.

Returns:

Linear matter power spectrum; \({\rm Mpc}^3\).

Return type:

(float or array)

linear_power(k, a, *, p_of_k_a='delta_matter:delta_matter')

The linear power spectrum.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • k (float or array) – Wavenumber; \({\rm Mpc}^{-1}\).

  • a (float or array) – Scale factor.

  • p_of_k_a (str) – string specifying the power spectrum to compute (which should be stored in cosmo). Defaults to the linear matter power spectrum.

Returns:

Linear power spectrum.

Return type:

(float or array)

lookback_time(a)

Difference of the age of the Universe between some scale factor and today, in \(\rm Gyr\).

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or (na,) array_like) – Scale factor(s) normalized to 1 today.

Returns:

t_L – Lookback time at a. nan if a is out of bounds of the spline parametets stored in cosmo.

Return type:

float or (na,) ndarray

luminosity_distance(a)

Luminosity distance.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

Returns:

Luminosity distance; Mpc.

Return type:

(float or array)

mass2radius_lagrangian(M)

Returns Lagrangian radius for a halo of mass \(M\). The Lagrangian radius is defined as that enclosing the mass of the halo assuming a homogeneous Universe.

\[R = \left(\frac{3\,M}{4\pi\,\rho_{M,0}}\right)^{1/3}\]
Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • M (float or array) – halo mass in units of \(M_\odot\).

Returns:

Lagrangian radius in comoving Mpc.

Return type:

(float or array)

nonlin_matter_power(k, a)

The nonlinear matter power spectrum

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • k (float or array) – Wavenumber; \({\rm Mpc}^{-1}\).

  • a (float or array) – Scale factor.

Returns:

Nonlinear matter power spectrum; \({\rm Mpc}^3\).

Return type:

(float or array)

nonlin_power(k, a, *, p_of_k_a='delta_matter:delta_matter')

The non-linear power spectrum.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • k (float or array) – Wavenumber; \({\rm Mpc}^{-1}\).

  • a (float or array) – Scale factor.

  • p_of_k_a (str) – string specifying the power spectrum to compute (which should be stored in cosmo). Defaults to the non-linear matter power spectrum.

Returns:

Non-linear power spectrum.

Return type:

(float or array)

omega_x(a, species)

Density fraction of a given species at a redshift different than z=0.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

  • species (str) –

    species type. Should be one of

    • ’matter’: cold dark matter, massive neutrinos, and baryons

    • ’dark_energy’: cosmological constant or otherwise

    • ’radiation’: relativistic species besides massless neutrinos

    • ’curvature’: curvature density

    • ’neutrinos_rel’: relativistic neutrinos

    • ’neutrinos_massive’: massive neutrinos

    • ’critical’

Returns:

Density fraction of a given species at a scale factor.

Return type:

(float or array)

parse_pk(p_of_k_a=None)

Helper to retrieve the right Pk2D object.

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • p_of_k_a (Pk2D, str or None) – 3D Power spectrum to return. If a Pk2D, it is just returned. If None or “linear”, the linear power spectrum stored in cosmo is returned. If “nonlinear”, the nonlinear matter power spectrum stored in cosmo is returned.

Returns:

Pk2D object corresponding to p_of_k_a.

parse_pk2d(p_of_k_a='delta_matter:delta_matter', *, is_linear=False)

Return the C-level f2d spline associated with a Pk2D object.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • p_of_k_a (Pk2D or str) – if a Pk2D object, its f2d spline will be used. If a string, the linear or non-linear power spectrum stored by cosmo under this name will be used. Defaults to the matter power spectrum stored in cosmo.

  • is_linear (bool) – if True, and if p_of_k_a is a string or None, the linear version of the corresponding power spectrum will be used (otherwise it’ll be the non-linear version).

rho_x(a, species, *, is_comoving=False)

Physical or comoving density as a function of scale factor.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • a (float or array) – Scale factor(s), normalized to 1 today.

  • species (str) –

    species type. Should be one of

    • ’matter’: cold dark matter, massive neutrinos, and baryons

    • ’dark_energy’: cosmological constant or otherwise

    • ’radiation’: relativistic species besides massless neutrinos

    • ’curvature’: curvature density

    • ’neutrinos_rel’: relativistic neutrinos

    • ’neutrinos_massive’: massive neutrinos

    • ’critical’

  • is_comoving (bool) – either physical (False, default) or comoving (True)

Returns:

Physical density of a given species at a scale factor, in units of \(M_\odot / {\rm Mpc}^3\).

Return type:

(float or array)

scale_factor_of_chi(chi)

Scale factor, a, at a comoving radial distance chi.

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • chi (float or array) – Comoving radial distance(s); Mpc.

Returns:

Scale factor(s), normalized to 1 today.

Return type:

(float or array)

sigma2_B_disc(a_arr=None, *, fsky=1.0, p_of_k_a='delta_matter:delta_matter')

Returns the variance of the projected linear density field over a circular disc covering a sky fraction fsky as a function of scale factor. This is given by

\[\sigma^2_B(z) = \int_0^\infty \frac{k\,dk}{2\pi} P_L(k,z)\,\left[\frac{2J_1(k R(z))}{k R(z)}\right]^2,\]

where \(R(z)\) is the corresponding radial aperture as a function of redshift. This quantity can be used to compute the super-sample covariance (see angular_cl_cov_SSC()).

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • a_arr (float, array or None) – an array of scale factor values at which to evaluate the projected variance. If None, a default sampling will be used.

  • fsky (float) – sky fraction.

  • p_of_k_a (Pk2D or str) – Linear power spectrum to use. Defaults to the internal linear power spectrum from cosmo.

Returns:

Tuple containing

  • a_arr (array): an array of scale factor values at which the projected variance has been evaluated. Only returned if a_arr is None on input.

  • sigma2_B (float or array): projected variance.

sigma2_B_from_mask(a_arr=None, *, mask_wl=None, p_of_k_a='delta_matter:delta_matter')

Returns the variance of the projected linear density field, given the angular power spectrum of the footprint mask and scale factor. This is given by

\[\sigma^2_B(z) = \frac{1}{\chi^2(z)}\sum_\ell P_L\left(\frac{\ell+\frac{1}{2}}{\chi(z)},z\right)\, \sum_{m=-\ell}^\ell W^A_{\ell m} {W^B}^*_{\ell m},\]

where \(W^A_{\ell m}\) and \(W^B_{\ell m}\) are the spherical harmonics decomposition of the footprint masks of fields A and B, normalized by their area. This quantity can be used to compute the super-sample covariance (see angular_cl_cov_SSC()).

Parameters:
  • cosmo (Cosmology) – a Cosmology object.

  • a_arr (float, array or None) – an array of scale factor values at which to evaluate the projected variance.

  • mask_wl (array) – Array with the angular power spectrum of the masks. The power spectrum should be given at integer multipoles, starting at \(\ell=0\). The power spectrum is normalized as \({\tt mask\_wl}=\sum_m W^A_{\ell m} {W^B}^*_{\ell m}\). It is the responsibility of the user to the provide the mask power out to sufficiently high ell for their required precision.

  • p_of_k_a (Pk2D or str) – Linear power spectrum to use. Defaults to the internal linear power spectrum from cosmo.

Returns:

Tuple containing

  • a_arr (array): an array of scale factor values at which the projected variance has been evaluated. Only returned if a_arr is None on input.

  • sigma2_B (float or array): projected variance.

sigma8(*, p_of_k_a='delta_matter:delta_matter')

RMS variance in a top-hat sphere of radius \(8\,{\rm Mpc}/h\), (with the value of \(h\) extracted from cosmo) at \(z=0\).

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • p_of_k_a (Pk2D, str or None) – power spectrum to integrate. If a string, it must correspond to one of the linear power spectra stored in cosmo (e.g. 'delta_matter:delta_matter').

Returns:

\(\sigma_8\).

Return type:

float

sigmaM(M, a)

RMS on the scale of a halo of mass \(M\). Calculated as \(\sigma_R\) (see sigmaR()) with \(R\) being the Lagrangian radius of a halo of mass \(M\) (see mass2radius_lagrangian()).

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • M (float or array) – Halo masses.

  • a (float) – scale factor.

Returns:

RMS variance of halo mass.

Return type:

(float or array)

sigmaR(R, a=1, *, p_of_k_a='delta_matter:delta_matter')

RMS of the matter overdensity a top-hat sphere of radius \(R\).

\[\sigma_R^2(z)=\frac{1}{2\pi^2}\int dk\,k^2\,P(k,z)\, |W(kR)|^2,\]

with \(W(x)=(3\sin(x)-x\cos(x))/x^3\).

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • R (float or array) – Radius; Mpc.

  • a (float) – optional scale factor.

  • p_of_k_a (Pk2D, str or None) – power spectrum to integrate. If a string, it must correspond to one of the linear power spectra stored in cosmo (e.g. 'delta_matter:delta_matter').

Returns:

\(\sigma_R\).

Return type:

(float or array)

sigmaV(R, a=1, *, p_of_k_a='delta_matter:delta_matter')

RMS of the linear displacement field in a top-hat sphere of radius R.

\[\sigma_V^2(z)=\frac{1}{6\pi^2}\int dk\,P(k,z)\,|W(kR)|^2,\]

with \(W(x)=(3\sin(x)-x\cos(x))/x^3\).

Parameters:
  • cosmo (Cosmology) – Cosmological parameters.

  • R (float or array) – Radius; Mpc.

  • a (float) – optional scale factor.

  • p_of_k_a (Pk2D, str or None) – power spectrum to integrate. If a string, it must correspond to one of the linear power spectra stored in cosmo (e.g. 'delta_matter:delta_matter').

Returns:

\(\sigma_V\) (\({\rm Mpc}\)).

Return type:

(float or array)

sigma_critical(*, a_lens, a_source)

Returns the critical surface mass density.

\[\Sigma_{\mathrm{crit}} = \frac{c^2}{4\pi G} \frac{D_{\rm{s}}}{D_{\rm{l}}D_{\rm{ls}}},\]

where \(c\) is the speed of light, \(G\) is the gravitational constant, and \(D_i\) is the angular diameter distance. The labels \(i = \{s,\,l,\,ls\}\) denote the distances to the source, lens, and between source and lens, respectively.

Parameters:
  • cosmo (Cosmology) – A Cosmology object.

  • a_lens (float) – lens’ scale factor.

  • a_source (float or array) – source’s scale factor.

Returns:

\(\Sigma_{\mathrm{crit}}\) in units of \(M_{\odot}/{\rm Mpc}^2\)

Return type:

(float or array)

tSZTracer(*, z_max=6.0, n_chi=1024)

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.

translate_IA_norm(*, z, a1=1.0, a1delta=None, a2=None, Om_m2_for_c2=False, Om_m_fid=0.3)

Function to convert from \(A_{ia}\) values to \(c_{ia}\) values, for the intrinsic alignment bias parameters using the standard convention of Blazek et al. 2019 or the variant used by the Dark Energy Survey analysis.

Parameters:
  • cosmo (Cosmology) – cosmology object.

  • z (float or array) – z value(s) where amplitude is evaluated.

  • a1 (float or array) – IA \(A_1\) at input z values.

  • a1delta (float or array) – IA \(A_{1\delta}\) at input z values.

  • a2 (float or array) – IA \(A_2\) at input z values.

  • Om_m2_for_c2 (bool) – True to use the Blazek et al. 2019 convention of \(\Omega_m^2\) scaling.

  • Om_m_fid (float) – Value for Blazek et al. 2019 scaling.

Returns:

Tuple of IA bias parameters

  • c1 (float or array): IA \(C_1\) at input z values.

  • c1delta (float or array): IA \(C_{1\delta}\) at input z values.

  • c2 (float or array): IA \(C_2\) at input z values.

pyccl.cosmology.CosmologyVanillaLCDM(**kwargs)[source]

A cosmology with typical flat Lambda-CDM parameters (Omega_c=0.25, Omega_b = 0.05, Omega_k = 0, sigma8 = 0.81, n_s = 0.96, h = 0.67, no massive neutrinos) for quick instantiation.

Parameters:

**kwargs (dict) – a dictionary of parameters passed as arguments to the Cosmology constructor. It should not contain any of the \(\Lambda\)-CDM parameters (“Omega_c”, “Omega_b”, “n_s”, “sigma8”, “A_s”, “h”), since these are fixed.

class pyccl.cosmology.CosmologyCalculator(self, *, Omega_c=None, Omega_b=None, h=None, n_s=None, sigma8=None, A_s=None, Omega_k=0.0, Omega_g=None, Neff=None, m_nu=0.0, mass_split='normal', w0=-1.0, wa=0.0, T_CMB=2.7255, T_ncdm=0.71611, mg_parametrization=None, background=None, growth=None, pk_linear=None, pk_nonlin=None, nonlinear_model=None)[source]

Bases: Cosmology

A “calculator-mode” CCL Cosmology object. This allows users to build a cosmology from a set of arrays describing the background expansion, linear growth factor and linear and non-linear power spectra, which can then be used to compute more complex observables (e.g. angular power spectra or halo-model quantities). These are stored in background, growth, pk_linear and pk_nonlin.

Note

Although in principle these arrays should suffice to compute most observable quantities some calculations implemented in CCL (e.g. the halo mass function) requires knowledge of basic cosmological parameters such as \(\Omega_M\). For this reason, users must pass a minimal set of \(\Lambda\)-CDM cosmological parameters.

Parameters:
  • Omega_c (float) – Cold dark matter density fraction.

  • Omega_b (float) – Baryonic matter density fraction.

  • h (float) – Hubble constant divided by 100 km/s/Mpc; unitless.

  • A_s (float) – Power spectrum normalization. Exactly one of A_s and sigma_8 is required.

  • sigma8 (float) – Variance of matter density perturbations at an 8 Mpc/h scale. Exactly one of A_s and sigma_8 is required.

  • n_s (float) – Primordial scalar perturbation spectral index.

  • Omega_k (float) – Curvature density fraction. Defaults to 0.

  • Omega_g (float) – Density in relativistic species except massless neutrinos. The default of None corresponds to setting this from the CMB temperature. Note that if a non-None value is given, this may result in a physically inconsistent model because the CMB temperature will still be non-zero in the parameters.

  • Neff (float) – Effective number of massless neutrinos present. Defaults to 3.044.

  • m_nu (float or array) – Mass in eV of the massive neutrinos present. Defaults to 0. If a sequence is passed, it is assumed that the elements of the sequence represent the individual neutrino masses.

  • mass_split (str) – Type of massive neutrinos. Should be one of ‘single’, ‘equal’, ‘normal’, ‘inverted’. ‘single’ treats the mass as being held by one massive neutrino. The other options split the mass into 3 massive neutrinos. Ignored if a sequence is passed in m_nu. Default is ‘normal’.

  • w0 (float) – First order term of dark energy equation of state. Defaults to -1.

  • wa (float) – Second order term of dark energy equation of state. Defaults to 0.

  • T_CMB (float) – The CMB temperature today. The default value is 2.7255.

  • T_ncdm (float) – Non-CDM temperature in units of photon temperature. The default is the same as in the base class

  • mg_parametrization (ModifiedGravity or None) – The modified gravity parametrization to use. Options are None (no MG), or a ModifiedGravity object. Currently, only MuSigmaMG is supported.

  • background (dict) – a dictionary describing the background expansion. It must contain three mandatory entries: 'a': an array of monotonically ascending scale-factor values. 'chi': an array containing the values of the comoving radial distance (in units of Mpc) at the scale factor values stored in a. ‘h_over_h0’: an array containing the Hubble expansion rate at the scale factor values stored in a, divided by its value today (at a=1).

  • growth (dict) – a dictionary describing the linear growth of matter fluctuations. It must contain three mandatory entries: 'a': an array of monotonically ascending scale-factor values. 'growth_factor': an array containing the values of the linear growth factor \(D(a)\) at the scale factor values stored in a. ‘growth_rate’: an array containing the growth rate \(f(a)\equiv d\log D/d\log a\) at the scale factor values stored in a.

  • pk_linear (dict) – a dictionary containing linear power spectra. It must contain the following mandatory entries: 'a': an array of scale factor values. 'k': an array of comoving wavenumbers in units of inverse Mpc. 'delta_matter:delta_matter': a 2D array of shape (n_a, n_k), where n_a and n_k are the lengths of 'a' and 'k' respectively, containing the linear matter power spectrum \(P(k,a)\). This dictionary may also contain other entries with keys of the form 'q1:q2', containing other cross-power spectra between quantities 'q1' and 'q2'.

  • pk_nonlin (dict) – a dictionary containing non-linear power spectra. It must contain the following mandatory entries: 'a': an array of scale factor values. 'k': an array of comoving wavenumbers in units of inverse Mpc. If nonlinear_model is None, it should also contain 'delta_matter:delta_matter': a 2D array of shape (n_a, n_k), where n_a and n_k are the lengths of 'a' and 'k' respectively, containing the non-linear matter power spectrum \(P(k,a)\). This dictionary may also contain other entries with keys of the form 'q1:q2', containing other cross-power spectra between quantities 'q1' and 'q2'.

  • nonlinear_model (str, dict or None) – model to compute non-linear power spectra. If a string, the associated non-linear model will be applied to all entries in pk_linear which do not appear in pk_nonlin. If a dictionary, it should contain entries of the form 'q1:q2': model, where model is a string designating the non-linear model to apply to the 'q1:q2' power spectrum, which must also be present in pk_linear. If model is None, this non-linear power spectrum will not be calculated. If nonlinear_model is None, no additional non-linear power spectra will be computed. The only non-linear model supported is 'halofit', corresponding to the “HALOFIT” transformation of Takahashi et al. 2012.