pyccl.core module

The core functionality of ccl, including the core data types. This includes the cosmology and parameters objects used to instantiate a model from which one can compute a set of theoretical predictions.

Supported Models for the Power Spectrum, Mass Function, etc.

The classes in this module accept strings indicating which model to use for various physical quantities (e.g., the transfer function). The various options are as follows.

transfer_function options
  • ‘emulator’: the transfer function defined by the Comsic Emu
  • ‘fitting_function’: the Eisenstein and Hu (1998) fitting function
  • ‘eisenstein_hu’: the Eisenstein and Hu (1998) fitting function
  • ‘bbks’: the BBKS approximation
  • ‘boltzmann’: use CLASS to compute the transfer function
  • ‘boltzmann_class’: use CLASS to compute the transfer function
  • ‘class’: use CLASS to compute the transfer function
  • ‘boltzmann_camb’: not implemented
  • ‘camb’: not implemented
matter_power_spectrum options
  • ‘halo_model’: use a halo model
  • ‘halofit’: use HALOFIT
  • ‘linear’: neglect non-linear power spectrum contributions
  • ‘emu’: use the Cosmic Emu
baryons_power_spectrum options
  • ‘nobaryons’: neglect baryonic contributions to the power spectrum
  • ‘bcm’: use the baryonic correction model
mass_function options
  • ‘tinker’: the Tinker et al. (2008) mass function
  • ‘tinker10’: the Tinker et al. (2010) mass function
  • ‘watson’: the Watson et al. mass function
  • ‘angulo’: the Angulo et al. mass function
  • ‘shethtormen’: the Sheth and Tormen mass function
halo_concentration options
  • ‘bhattacharya2011’: Bhattacharya et al. (2011) relation
  • ‘duffy2008’: Duffy et al. (2008) relation
  • ‘constant_concentration’: use a constant concentration
mnu_type options

This parameter specifies the model for massive neutrinos.

  • ‘list’: specify each mass yourself in eV
  • ‘sum’: use the normal hierarchy to convert total mass to individual masses
  • ‘sum_inverted’: use the inverted hierarchy to convert total mass to individual masses
  • ‘sum_equal’: assume equal masses when converting the total mass to individual masses
emulator_neutrinos options

This parameter specifies how to handle inconsistencies in the treatment of neutrinos between the Cosmic Emu (equal masses) and other models.

  • ‘strict’: fail unless things are absolutely consistent
  • ‘equalize’: redistribute the total mass equaly before using the Cosmic Emu. This option may result in slight internal inconsistencies in the physical model assumed for neutrinos.

Controlling Splines and Numerical Accuracy

The internal splines and integration accuracy are controlled by the attributes of Cosmology.cosmo.spline_params and Cosmology.cosmo.gsl_params. These should be set after instantiation, but before the object is used. For example, you can set the generic relative accuracy for integration by executing c = Cosmology(...); c.cosmo.gsl_params.INTEGRATION_EPSREL = 1e-5. The default values for these parameters are located in src/ccl_core.c.

The intrnal splines are controlled by the following parameters.

  • A_SPLINE_NLOG: the number of logarithmically spaced bins between A_SPLINE_MINLOG and A_SPLINE_MIN.
  • A_SPLINE_NA: the number of linearly spaced bins between A_SPLINE_MIN and A_SPLINE_MAX.
  • A_SPLINE_MINLOG: the minimum value of the scale factor splines used for distances, etc.
  • A_SPLINE_MIN: the transition scale factor between logarithmically spaced spline points and linearly spaced spline points.
  • A_SPLINE_MAX: the the maximum value of the scale factor splines used for distances, etc.
  • LOGM_SPLINE_NM: the number of logarithmically spaced values in mass for splines used in the computation of the halo mass function.
  • LOGM_SPLINE_MIN: the base-10 logarithm of the minimum halo mass for splines used in the computation of the halo mass function.
  • LOGM_SPLINE_MAX: the base-10 logarithm of the maximum halo mass for splines used in the computation of the halo mass function.
  • LOGM_SPLINE_DELTA: the step in base-10 logarithmic units for computing finite difference derivatives in the computation of the mass function.
  • A_SPLINE_NLOG_PK: the number of logarithmically spaced bins between A_SPLINE_MINLOG_PK and A_SPLINE_MIN_PK.
  • A_SPLINE_NA_PK: the number of linearly spaced bins between A_SPLINE_MIN_PK and A_SPLINE_MAX.
  • A_SPLINE_MINLOG_PK: the minimum value of the scale factor used for the power spectrum splines.
  • A_SPLINE_MIN_PK: the transition scale factor between logarithmically spaced spline points and linearly spaced spline points for the power spectrum.
  • K_MIN: the minimum wavenumber for the power spectrum splines for analytic models (e.g., BBKS, Eisenstein & Hu, etc.).
  • K_MAX: the maximum wavenumber for the power spectrum splines for analytic models (e.g., BBKS, Eisenstein & Hu, etc.).
  • K_MAX_SPLINE: the maximum wavenumber for the power spectrum splines for numerical models (e.g., ComsicEmu, CLASS, etc.).
  • N_K: the number of spline nodes per decade for the power spectrum splines.
  • N_K_3DCOR: the number of spline points in wavenumber per decade used for computing the 3D correlation function.
  • ELL_MIN_CORR: the minimum value of the spline in angular wavenumber for correlation function computations with FFTlog.
  • ELL_MAX_CORR: the maximum value of the spline in angular wavenumber for correlation function computations with FFTlog.
  • N_ELL_CORR: the number of logarithmically spaced bins in angular wavenumber between ELL_MIN_CORR and ELL_MAX_CORR.

The numrical accuracy of GSL computations are controlled by the following parameters.

  • N_ITERATION: the size of the GSL workspace for numerical integration.
  • INTEGRATION_GAUSS_KRONROD_POINTS: the Gauss-Kronrod quadrature rule used for adaptive integrations.
  • INTEGRATION_EPSREL: the relative error tolerance for numerical integration; used if not specified by a more specific parameter.
  • INTEGRATION_LIMBER_GAUSS_KRONROD_POINTS: the Gauss-Kronrod quadrature rule used for adaptive integrations on subintervals for Limber integrals.
  • INTEGRATION_LIMBER_EPSREL: the relative error tolerance for numerical integration of Limber integrals.
  • INTEGRATION_DISTANCE_EPSREL: the relative error tolerance for numerical integration of distance integrals.
  • INTEGRATION_SIGMAR_EPSREL: the relative error tolerance for numerical integration of power spectrum variance intrgals for the mass function.
  • ROOT_EPSREL: the relative error tolerance for root finding used to invert the relationship between comoving distance and scale factor.
  • ROOT_N_ITERATION: the maximum number of iterations used to for root finding to invert the relationship between comoving distance and scale factor.
  • ODE_GROWTH_EPSREL: the relative error tolerance for integrating the linear growth ODEs.
  • EPS_SCALEFAC_GROWTH: 10x the starting step size for integrating the linear growth ODEs and the scale factor of the initial condition for the linear growth ODEs.
  • HM_MMIN: the minimum mass for halo model integrations.
  • HM_MMAX: the maximum mass for halo model integrations.
  • HM_EPSABS: the absolute error tolerance for halo model integrations.
  • HM_EPSREL: the relative error tolerance for halo model integrations.
  • HM_LIMIT: the size of the GSL workspace for halo moodel integrations.
  • HM_INT_METHOD: the Gauss-Kronrod quadrature rule used for adaptive integrations for the halo model comptutations.

Specifying Physical Constants

The values of physical constants are set globally. These can be changed by assigning a new value to the attributes of pyccl.physical_constants. The following constants are defined and their default values are located in src/ccl_core.c. Note that the neutrino mass splittings are taken from Lesgourgues & Pastor (2012; 1212.6154).

basic physical constants
  • CLIGHT_HMPC: speed of light / H0 in units of Mpc/h
  • GNEWT: Newton’s gravitational constant in units of m^3/Kg/s^2
  • SOLAR_MASS: solar mass in units of kg
  • MPC_TO_METER: conversion factor for Mpc to meters.
  • PC_TO_METER: conversion factor for parsecs to meters.
  • RHO_CRITICAL: critical density in units of M_sun/h / (Mpc/h)^3
  • KBOLTZ: Boltzmann constant in units of J/K
  • STBOLTZ: Stefan-Boltzmann constant in units of kg/s^3 / K^4
  • HPLANCK: Planck’s constant in units kg m^2 / s
  • CLIGHT: speed of light in m/s
  • EV_IN_J: conversion factor between electron volts and Joules
  • TCMB: temperature of the CMB in K
  • TNCDM: temperature of the cosmological neutrino background in K
neutrino mass splittings
  • DELTAM12_sq: squared mass difference between eigenstates 2 and 1.
  • DELTAM13_sq_pos: squared mass difference between eigenstates 3 and 1 for the normal hierarchy.
  • DELTAM13_sq_neg: squared mass difference between eigenstates 3 and 1 for the inverted hierarchy.
class pyccl.core.Cosmology(Omega_c=None, Omega_b=None, h=None, n_s=None, sigma8=None, A_s=None, Omega_k=0.0, Omega_g=None, Neff=3.046, m_nu=0.0, mnu_type=None, w0=-1.0, wa=0.0, bcm_log10Mc=14.079181246047625, bcm_etab=0.5, bcm_ks=55.0, z_mg=None, df_mg=None, transfer_function='boltzmann_class', matter_power_spectrum='halofit', baryons_power_spectrum='nobaryons', mass_function='tinker10', halo_concentration='duffy2008', emulator_neutrinos='strict')[source]

Bases: object

A cosmology including parameters and associated data.

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. Note however that this approximation is common for late-time LSS computations.

Note

BCM stands for the “baryonic correction model” of Schneider & Teyssier (2015; https://arxiv.org/abs/1510.06034). See the DESC Note for details.

Note

After instantiation, you can set parameters related to the internal splines and numerical integration accuracy by setting the values of the attributes of Cosmology.cosmo.spline_params and Cosmology.cosmo.gsl_params. For example, you can set the generic relative accuracy for integration by executing c = Cosmology(...); c.cosmo.gsl_params.INTEGRATION_EPSREL = 1e-5. See the module level documetaion of pyccl.core for details.

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, optional) – Curvature density fraction. Defaults to 0.
  • Omega_g (float, optional) – 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, optional) – Effective number of massless neutrinos present. Defaults to 3.046.
  • m_nu (float, optional) – Total mass in eV of the massive neutrinos present. Defaults to 0.
  • mnu_type (str, optional) – The type of massive neutrinos.
  • w0 (float, optional) – First order term of dark energy equation of state. Defaults to -1.
  • wa (float, optional) – Second order term of dark energy equation of state. Defaults to 0.
  • bcm_log10Mc (float, optional) – One of the parameters of the BCM model. Defaults to np.log10(1.2e14).
  • bcm_etab (float, optional) – One of the parameters of the BCM model. Defaults to 0.5.
  • bcm_ks (float, optional) – One of the parameters of the BCM model. Defaults to 55.0.
  • df_mg (array_like, optional) – Perturbations to the GR growth rate as a function of redshift \(\Delta f\). Used to implement simple modified growth scenarios.
  • z_mg (array_like, optional) – Array of redshifts corresponding to df_mg.
  • transfer_function (str, optional) – The transfer function to use. Defaults to ‘boltzmann_class’.
  • matter_power_spectrum (str, optional) – The matter power spectrum to use. Defaults to ‘halofit’.
  • baryons_power_spectrum (str, optional) – The correction from baryonic effects to be implemented. Defaults to ‘nobaryons’.
  • mass_function (str, optional) – The mass function to use. Defaults to ‘tinker10’ (2010).
  • halo_concentration (str, optional) – The halo concentration relation to use. Defaults to Duffy et al. (2008) ‘duffy2008’.
  • emulator_neutrinosstr, optional): If using the emulator for the power spectrum, specified treatment of unequal neutrinos. Options are ‘strict’, which will raise an error and quit if the user fails to pass either a set of three equal masses or a sum with mnu_type = ‘sum_equal’, and ‘equalize’, which will redistribute masses to be equal right before calling the emualtor but results in internal inconsistencies. Defaults to ‘strict’.
compute_distances()[source]

Compute the distance splines.

compute_growth()[source]

Compute the growth function.

compute_power()[source]

Compute the power spectrum.

has_distances()[source]

Checks if the distances have been precomputed.

Returns:True if precomputed, False otherwise.
Return type:bool
has_growth()[source]

Checks if the growth function has been precomputed.

Returns:True if precomputed, False otherwise.
Return type:bool
has_power()[source]

Checks if the power spectra have been precomputed.

Returns:True if precomputed, False otherwise.
Return type:bool
has_sigma()[source]

Checks if sigma8 has been computed.

Returns:True if precomputed, False otherwise.
Return type:bool
classmethod read_yaml(filename)[source]

Read the parameters from a YAML file.

Parameters:filename (str) –
status()[source]

Get error status of the ccl_cosmology object.

Note

error statuses are currently under development.

Returns:str containing the status message.
write_yaml(filename)[source]

Write a YAML representation of the parameters to file.

Parameters:filename (str) –
pyccl.core.check(status, cosmo=None)[source]

Check the status returned by a ccllib function.

Parameters:
  • status (int or core.error_types) – Flag or error describing the success of a function.
  • cosmo (Cosmology, optional) – A Cosmology object.