Cosmological theory codes and likelihoods

Models in Cosmology are usually split in two: \(\mathcal{M}=\mathcal{T}+\mathcal{E}\), where

  • \(\mathcal{T}\), the theoretical model, is used to compute observable quantities \(\mathcal{O}\)

  • \(\mathcal{E}\), the experimental model, accounts for instrumental errors, foregrounds… when comparing the theoretical observable with some data \(\mathcal{D}\).

In practice the theoretical model is encapsulated in one or more theory codes (CLASS, CAMB…) and the experimental model in a likelihood, which gives the probability of the data being a realization of the given observable in the context of the experiment:

\[\mathcal{L}\left[\mathcal{D}\,|\,\mathcal{M}\right] = \mathcal{L}\left[\mathcal{D}\,|\,\mathcal{O},\mathcal{E}\right]\]

Each iteration of a sampler reproduces the model using the following steps:

  1. A new set of theory+experimental parameters is proposed by the sampler.

  2. The theory parameters are passed to the theory codes, which compute one or more observables.

  3. The experimental parameters are passed to each of the likelihoods, which in turn ask the theory code for the current value of the observables they need, and use all of that to compute a log-probability of the data.

cobaya wraps the most popular cosmological codes under a single interface, documented below. The codes themselves are documented in the next sections, followed by the internal likelihoods included in cobaya.

Cosmological theory code

class theories.cosmo.BoltzmannBase(info=mappingproxy({}), name=None, timing=None, packages_path=None, initialize=True, standalone=True)
initialize()

Initializes the class (called from __init__, before other initializations).

initialize_with_params()

Additional initialization after requirements called and input_params and output_params have been assigned (but provider and assigned requirements not yet set).

get_allow_agnostic()

Whether it is allowed to pass all unassigned input parameters to this component (True) or whether parameters must be explicitly specified (False).

Returns:

True or False

get_param(p)

Interface function for likelihoods and other theory components to get derived parameters.

must_provide(**requirements)

Specifies the quantities that this Boltzmann code is requested to compute.

Typical requisites in Cosmology (as keywords, case-insensitive):

  • Cl={...}: CMB lensed power spectra, as a dictionary {[spectrum]: l_max}, where the possible spectra are combinations of "t", "e", "b" and "p" (lensing potential). Get with get_Cl().

  • unlensed_Cl={...}: CMB unlensed power spectra, as a dictionary {[spectrum]: l_max}, where the possible spectra are combinations of "t", "e", "b". Get with get_unlensed_Cl().

  • [BETA: CAMB only; notation may change!] source_Cl={...}: \(C_\ell\) of given sources with given windows, e.g.: source_name: {"function": "spline"|"gaussian", [source_args]; for now, [source_args] follows the notation of CAMBSources. It can also take "lmax": [int], "limber": True if Limber approximation desired, and "non_linear": True if non-linear contributions requested. Get with get_source_Cl().

  • Pk_interpolator={...}: Matter power spectrum interpolator in \((z, k)\). Takes "z": [list_of_evaluated_redshifts], "k_max": [k_max], "nonlinear": [True|False], "vars_pairs": [["delta_tot", "delta_tot"], ["Weyl", "Weyl"], [...]]}. Notice that k_min cannot be specified. To reach a lower one, use extra_args in CAMB to increase accuracyboost and TimeStepBoost, or in CLASS to decrease k_min_tau0. The method get_Pk_interpolator() to retrieve the interpolator admits extrapolation limits extrap_[kmax|kmin]. It is recommended to use extrap_kmin to reach the desired k_min, and increase precision parameters only as much as necessary to improve over the interpolation. All \(k\) values should be in units of \(1/\mathrm{Mpc}\).

    Non-linear contributions are included by default. Note that the non-linear setting determines whether non-linear corrections are calculated; the get_Pk_interpolator() method also has a non-linear argument to specify if you want the linear or non-linear spectrum returned (to have both linear and non-linear spectra available request a tuple (False,True) for the non-linear argument).

  • Pk_grid={...}: similar to Pk_interpolator except that rather than returning a bicubic spline object it returns the raw power spectrum grid as a (k, z, P(z,k)) set of arrays. Get with get_Pk_grid().

  • sigma_R={...}: RMS linear fluctuation in spheres of radius \(R\) at redshifts \(z\). Takes "z": [list_of_evaluated_redshifts], "k_max": [k_max], "vars_pairs": [["delta_tot", "delta_tot"],  [...]], "R": [list_of_evaluated_R]. Note that \(R\) is in \(\mathrm{Mpc}\), not \(h^{-1}\,\mathrm{Mpc}\). Get with get_sigma_R().

  • Hubble={'z': [z_1, ...]}: Hubble rates at the requested redshifts. Get it with get_Hubble().

  • Omega_b={'z': [z_1, ...]}: Baryon density parameter at the requested redshifts. Get it with get_Omega_b().

  • Omega_cdm={'z': [z_1, ...]}: Cold dark matter density parameter at the requested redshifts. Get it with get_Omega_cdm().

  • Omega_nu_massive={'z': [z_1, ...]}: Massive neutrinos’ density parameter at the requested redshifts. Get it with get_Omega_nu_massive().

  • angular_diameter_distance={'z': [z_1, ...]}: Physical angular diameter distance to the redshifts requested. Get it with get_angular_diameter_distance().

  • angular_diameter_distance_2={'z_pairs': [(z_1a, z_1b), (z_2a, z_2b)...]}: Physical angular diameter distance between the pairs of redshifts requested. If a 1d array of redshifts is passed as z_pairs, all possible combinations of two are computed and stored (not recommended if only a subset is needed). Get it with get_angular_diameter_distance_2().

  • comoving_radial_distance={'z': [z_1, ...]}: Comoving radial distance from us to the redshifts requested. Get it with get_comoving_radial_distance().

  • sigma8_z={'z': [z_1, ...]}: Amplitude of rms fluctuations \(\sigma_8\) at the redshifts requested. Get it with get_sigma8().

  • fsigma8={'z': [z_1, ...]}: Structure growth rate \(f\sigma_8\) at the redshifts requested. Get it with get_fsigma8().

  • Other derived parameters that are not included in the input but whose value the likelihood may need.

requested()

Returns the full set of requested cosmological products and parameters.

check_no_repeated_input_extra()

Checks that there are no repeated parameters between input and extra.

Should be called at initialisation, and at the end of every call to must_provide()

get_Cl(ell_factor=False, units='FIRASmuK2')

Returns a dictionary of lensed CMB power spectra and the lensing potential pp power spectrum.

Set the units with the keyword units=number|'muK2'|'K2'|'FIRASmuK2'|'FIRASK2'. The default is FIRASmuK2, which returns CMB \(C_\ell\)’s in FIRAS-calibrated \(\mu K^2\), i.e. scaled by a fixed factor of \((2.7255\cdot 10^6)^2\) (except for the lensing potential power spectrum, which is always unitless). The muK2 and K2 options use the model’s CMB temperature.

If ell_factor=True (default: False), multiplies the spectra by \(\ell(\ell+1)/(2\pi)\) (or by \([\ell(\ell+1)]^2/(2\pi)\) in the case of the lensing potential pp spectrum, and \([\ell(\ell+1)]^{3/2}/(2\pi)\) for the the cross spectra tp and ep).

get_unlensed_Cl(ell_factor=False, units='FIRASmuK2')

Returns a dictionary of unlensed CMB power spectra.

For units options, see get_Cl().

If ell_factor=True (default: False), multiplies the spectra by \(\ell(\ell+1)/(2\pi)\).

get_Hubble(z, units='km/s/Mpc')

Returns the Hubble rate at the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called. The available units are "km/s/Mpc" (i.e. \(cH(\mathrm(Mpc)^{-1})\)) and 1/Mpc.

get_Omega_b(z)

Returns the Baryon density parameter at the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called.

get_Omega_cdm(z)

Returns the Cold Dark Matter density parameter at the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called.

get_Omega_nu_massive(z)

Returns the Massive neutrinos’ density parameter at the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called.

get_angular_diameter_distance(z)

Returns the physical angular diameter distance in \(\mathrm{Mpc}\) to the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called.

get_angular_diameter_distance_2(z_pairs)

Returns the physical angular diameter distance between pairs of redshifts z_pairs in \(\mathrm{Mpc}\).

The redshift pairs must be a subset of those requested when must_provide() was called.

Return zero for pairs in which z1 > z2.

get_comoving_radial_distance(z)

Returns the comoving radial distance in \(\mathrm{Mpc}\) to the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called.

get_Pk_interpolator(var_pair=('delta_tot', 'delta_tot'), nonlinear=True, extrap_kmin=None, extrap_kmax=None)

Get a \(P(z,k)\) bicubic interpolation object (PowerSpectrumInterpolator).

In the interpolator returned, the input \(k\) and resulting \(P(z,k)\) are in units of \(1/\mathrm{Mpc}\) and \(\mathrm{Mpc}^3\) respectively (not in \(h^{-1}\) units).

Parameters:
  • var_pair – variable pair for power spectrum

  • nonlinear – non-linear spectrum (default True)

  • extrap_kmin – use log linear extrapolation from extrap_kmin up to min \(k\).

  • extrap_kmax – use log linear extrapolation beyond max \(k\) computed up to extrap_kmax.

Returns:

PowerSpectrumInterpolator instance.

get_Pk_grid(var_pair=('delta_tot', 'delta_tot'), nonlinear=True)

Get matter power spectrum, e.g. suitable for splining. Returned arrays may be bigger or more densely sampled than requested, but will include required values.

In the grid returned, \(k\) and \(P(z,k)\) are in units of \(1/\mathrm{Mpc}\) and \(\mathrm{Mpc}^3\) respectively (not in \(h^{-1}\) units), and \(z\) and \(k\) are in ascending order.

Parameters:
  • nonlinear – whether the linear or non-linear spectrum

  • var_pair – which power spectrum

Returns:

k, z, Pk, where k and z are 1-d arrays, and the 2-d array Pk[i,j] is the value of \(P(z,k)\) at z[i], k[j].

get_sigma_R(var_pair=('delta_tot', 'delta_tot'))

Get \(\sigma_R(z)\), the RMS power in a sphere of radius \(R\) at redshift \(z\).

Note that \(R\) is in \(\mathrm{Mpc}\), not \(h^{-1}\,\mathrm{Mpc}\), and \(z\) and \(R\) are returned in ascending order.

You may get back more values than originally requested, but the requested \(R\) and \(z\) should be in the returned arrays.

Parameters:

var_pair – which two fields to use for the RMS power

Returns:

z, R, sigma_R, where z and R are arrays of computed values, and sigma_R[i,j] is the value \(\sigma_R(z)\) for z[i], R[j].

get_source_Cl()

Returns a dict of power spectra of for the computed sources, with keys a tuple of sources ([source1], [source2]), and an additional key ell containing the multipoles.

get_sigma8_z(z)

Present day linear theory root-mean-square amplitude of the matter fluctuation spectrum averaged in spheres of radius 8 h^{−1} Mpc.

The redshifts must be a subset of those requested when must_provide() was called.

get_fsigma8(z)

Structure growth rate \(f\sigma_8\), as defined in eq. 33 of Planck 2015 results. XIII. Cosmological parameters, at the given redshift(s) z.

The redshifts must be a subset of those requested when must_provide() was called.

get_auto_covmat(params_info, likes_info)

Tries to get match to a database of existing covariance matrix files for the current model and data.

params_info should contain preferably the slow parameters only.

class theories.cosmo.PowerSpectrumInterpolator(z, k, P_or_logP, extrap_kmin=None, extrap_kmax=None, logP=False, logsign=1)

2D spline interpolation object (scipy.interpolate.RectBivariateSpline) to evaluate matter power spectrum as function of z and k.

This class is adapted from CAMB’s own P(k) interpolator, by Antony Lewis; it’s mostly interface-compatible with the original.

Parameters:
  • z – values of z for which the power spectrum was evaluated.

  • k – values of k for which the power spectrum was evaluated.

  • P_or_logP – Values of the power spectrum (or log-values, if logP=True).

  • logP – if True (default: False), log of power spectrum are given and used for the underlying interpolator.

  • logsign – if logP is True, P_or_logP is log(logsign*Pk)

  • extrap_kmax – if set, use power law extrapolation beyond kmax up to extrap_kmax; useful for tails of integrals.

property input_kmin

Minimum k for the interpolation (not incl. extrapolation range).

property input_kmax

Maximum k for the interpolation (not incl. extrapolation range).

property kmin

Minimum k of the interpolator (incl. extrapolation range).

property kmax

Maximum k of the interpolator (incl. extrapolation range).

check_ranges(z, k)

Checks that we are not trying to extrapolate beyond the interpolator limits.

P(z, k, grid=None)

Get the power spectrum at (z,k).

logP(z, k, grid=None)

Get the log power spectrum at (z,k). (or minus log power spectrum if islog and logsign=-1)

Cosmological theory code inheritance

Inheritance diagram of theorys