CAMB

Synopsis:

Managing the CAMB cosmological code

Author:

Jesus Torrado and Antony Lewis

This module imports and manages the CAMB cosmological code. It requires CAMB 1.5 or higher.

Note

If you use this cosmological code, please cite it as:
A. Lewis, A. Challinor, A. Lasenby, Efficient computation of CMB anisotropies in closed FRW (arXiv:astro-ph/9911177)
C. Howlett, A. Lewis, A. Hall, A. Challinor, CMB power spectrum parameter degeneracies in the era of precision cosmology (arXiv:1201.3654)

Usage

If you are using a likelihood that requires some observable from CAMB, simply add CAMB to the theory block.

You can specify any parameter that CAMB understands in the params block:

theory:
  camb:
    extra_args:
      [any param that CAMB understands, for FIXED and PRECISION]

params:
    [any param that CAMB understands, fixed, sampled or derived]

If you want to use your own version of CAMB, you need to specify its location with a path option inside the camb block. If you do not specify a path, CAMB will be loaded from the automatic-install packages_path folder, if specified, or otherwise imported as a globally-installed Python package. If you want to force that the global camb installation is used, pass path='global'. Cobaya will print at initialisation where CAMB was actually loaded from.

Access to CAMB computation products

You can retrieve CAMB computation products within likelihoods (or other pipeline components in general) or manually from a Model as long as you have added them as requisites; see Creating your own cosmological likelihood or Creating your own cosmological likelihood class for the likelihood case, and Using the model wrapper for the manual case.

The products that you can request and later retrieve are listed in must_provide().

If you would like to access a CAMB result that is not accessible that way, you can access the full CAMB results object CAMBdata directly by adding {"CAMBdata": None} to your requisites, and then retrieving it with provider.get_CAMBdata().

In general, the use of CAMBdata should be avoided in public code, since it breaks compatibility with other Boltzmann codes at the likelihood interface level. If you need a quantity for a public code that is not generally interfaced in must_provide(), let us know if you think it makes sense to add it.

Modifying CAMB

If you modify CAMB and add new variables, make sure that the variables you create are exposed in the Python interface (instructions here). If you follow those instructions you do not need to make any additional modification in Cobaya.

If your modification involves new computed quantities, add a retrieving method to CAMBdata, and see Access to CAMB computation products.

You can use the model wrapper to test your modification by evaluating observables or getting derived quantities at known points in the parameter space (set debug: True to get more detailed information of what exactly is passed to CAMB).

In your CAMB modification, remember that you can raise a CAMBParamRangeError or a CAMBError whenever the computation of any observable would fail, but you do not expect that observable to be compatible with the data (e.g. at the fringes of the parameter space). Whenever such an error is raised during sampling, the likelihood is assumed to be zero, and the run is not interrupted.

Note

If your modified CAMB has a lower version number than the minimum required by Cobaya, you will get an error at initialisation. You may still be able to use it by setting the option ignore_obsolete: True in the camb block (though you would be doing that at your own risk; ideally you should translate your modification to a newer CAMB version, in case there have been important fixes since the release of your baseline version).

Installation

Pre-requisites

cobaya calls CAMB using its Python interface, which requires that you compile CAMB using intel’s ifort compiler or the GNU gfortran compiler version 6.4 or later. To check if you have the latter, type gfortran --version in the shell, and the first line should look like

GNU Fortran ([your OS version]) [gfortran version] [release date]

Check that [gfortran's version] is at least 6.4. If you get an error instead, you need to install gfortran (contact your local IT service).

CAMB comes with binaries pre-built for Windows, so if you don’t need to modify the CAMB source code, no Fortran compiler is needed.

If you are using Anaconda you can also install a pre-compiled CAMB package from conda forge using

conda install -c conda-forge camb

Automatic installation

If you do not plan to modify CAMB, the easiest way to install it is using the automatic installation script. Just make sure that theory: camb: appears in one of the files passed as arguments to the installation script.

Manual installation (or using your own version)

If you are planning to modify CAMB or use an already modified version, you should not use the automatic installation script. Use the installation method that best adapts to your needs:

  • [Recommended for staying up-to-date] To install CAMB locally and keep it up-to-date, clone the CAMB repository in GitHub in some folder of your choice, say /path/to/theories/CAMB:

    $ cd /path/to/theories
    $ git clone --recursive https://github.com/cmbant/CAMB.git
    $ cd CAMB
    $ python setup.py build
    

    To update to the last changes in CAMB (master), run git pull from CAMB’s folder and re-build using the last command. If you do not want to use multiple versions of CAMB, you can also make your local installation available to python generally by installing it using

$ python -m pip install -e /path/to/CAMB
  • [Recommended for modifying CAMB] First, fork the CAMB repository in GitHub (follow these instructions) and then follow the same steps as above, substituting the second one with:

    $ git clone --recursive https://[YourGithubUser]@github.com/[YourGithubUser]/CAMB.git
    
  • To use your own version, assuming it’s placed under /path/to/theories/CAMB, just make sure it is compiled.

In the cases above, you must specify the path to your CAMB installation in the input block for CAMB (otherwise a system-wide CAMB may be used instead):

theory:
  camb:
    path: /path/to/theories/CAMB

Note

In any of these methods, if you intend to switch between different versions or modifications of CAMB you should not install CAMB as python package using pip install, as the official instructions suggest. It is not necessary if you indicate the path to your preferred installation as explained above.

Documentation of the camb wrapper class

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

CAMB cosmological Boltzmann code cite{Lewis:1999bs,Howlett:2012mh}.

initialize()

Importing CAMB from the correct path, if given.

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).

initialize_with_provider(provider)

Final initialization after parameters, provider and assigned requirements set. The provider is used to get the requirements of this theory using provider.get_X() and provider.get_param(‘Y’).

Parameters:

provider – the theory.Provider instance that should be used by this component to get computed requirements

get_can_support_params()

Get a list of parameters supported by this component, can be used to support parameters that don’t explicitly appear in the .yaml or class params attribute or are otherwise explicitly supported (e.g. via requirements)

Returns:

iterable of names of parameters

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

set_cl_reqs(reqs)

Sets some common settings for both lensed and unlensed Cl’s.

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.

add_to_redshifts(z)

Adds redshifts to the list of them for which CAMB computes perturbations.

set_collector_with_z_pool(k, zs, method, args=(), kwargs=mappingproxy({}), d=1)

Creates a collector for a z-dependent quantity, keeping track of the pool of z’s.

calculate(state, want_derived=True, **params_values_dict)

Do the actual calculation and store results in state dict

Parameters:
  • state – dictionary to store results

  • want_derived – whether to set state[‘derived’] derived parameters

  • params_values_dict – parameter values

Returns:

None or True if success, False for fail

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

Returns a dictionary of lensed total 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 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_lensed_scal_Cl(ell_factor=False, units='FIRASmuK2')

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

For units options, see get_Cl().

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

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_CAMBdata()

Get the CAMB result object (must have been requested as a requirement).

Returns:

CAMB’s CAMBdata result instance for the current parameters

get_can_provide_params()

Get a list of derived parameters that this component can calculate. The default implementation returns the result based on the params attribute set via the .yaml file or class params (with derived:True for derived parameters).

Returns:

iterable of parameter names

get_version()

Get version information for this component.

Returns:

string or dict of values or None

get_helper_theories()

Transfer functions are computed separately by camb.transfers, then this class uses the transfer functions to calculate power spectra (using A_s, n_s etc).

static get_import_path(path)

Returns the camb module import path if there is a compiled version of CAMB in the given folder. Otherwise, raises FileNotFoundError.