API Reference

cctoolkit

cosmology

cosmology.py

This module defines the CosmologyCalculator class, which manages cosmological parameters and calculations related to the power spectrum and other cosmological quantities using CAMB.

class cctoolkit.cosmology.CosmologyCalculator(params=None, power_spectrum=None, zmax=2.0, nz=100, var='cdm+b')[source]

Bases: object

A class to handle cosmological calculations using CAMB.

This class facilitates the calculation of various cosmological quantities, including power spectra, density parameters, and functions related to the Halo Mass Function (HMF). It includes methods for calculating the Peak Background Split (PBS) bias and corrected bias, utilizing the CAMB library.

Attributes:

paramsdict

Cosmological parameters for the model, including ‘H0’, ‘Ob0’, ‘Om0’, ‘sigma8’ or ‘As’, ‘ns’, ‘TCMB’, ‘mnu’, ‘number of massive neutrinos species’, ‘w0’, and ‘wa’.

knp.ndarray

Array of wavenumbers used in the power spectrum calculations.

z_valsnp.ndarray

Array of redshift values used for power spectrum calculations.

varstr

Variable used to compute the matter power spectrum, default is ‘cdm+b’.

Methods:

get_power_spectrum(z):

Returns the wavenumbers and dimensionless power spectrum at a given redshift.

set_cosmology(new_params):

Updates the cosmological parameters and recalculates the power spectrum normalization.

Omega_m(z):

Returns the matter density parameter Omega_m as a function of redshift.

Omega_nu(z):

Returns the neutrino density parameter Omega_nu as a function of redshift.

Omega_neutrino(z):

Returns the massless neutrino density parameter Omega_neutrino as a function of redshift.

Omega_photon(z):

Returns the photon density parameter Omega_photon as a function of redshift.

Omega_k(z):

Returns the curvature density parameter Omega_k as a function of redshift.

Omega_DE(z):

Returns the dark energy density parameter Omega_DE as a function of redshift.

wz(z):

Returns the dark energy equation of state at a given redshift.

critical_density(z):

Returns the critical density at a given redshift.

H(z):

Returns the Hubble parameter H at a given redshift.

peak_height(M, z):

Calculates the peak-height delta_c / sigma(R(M), z) for a given mass and redshift.

lagrangian_radius(M):

Calculates the radius of the Lagrangian patch corresponding to a given mass M.

sigma(R, z):

Calculates the RMS fluctuation sigma(R, z) for a given radius and redshift.

dlnsigma_dlnR(R, z):

Calculates the derivative of the logarithm of sigma with respect to the logarithm of R.

vfv(M, z, return_variables=False, halo_finder=’ROCKSTAR’, model=’castro23’):

Calculates the multiplicity function νf(ν) and related variables for a given mass and redshift.

dndlnM(M, z, halo_finder=’ROCKSTAR’, model=’castro23’):

Calculates the halo mass function dn/dlnM, representing the number density of halos per logarithmic mass interval.

pbs_bias(M, z, halo_finder=’ROCKSTAR’, hmf_model=’castro23’, return_variables=False):

Calculates the PBS bias for halos of mass M at redshift z.

bias(M, z, halo_finder=’ROCKSTAR’, hmf_model=’castro23’):

Calculates the corrected linear halo bias for given mass and redshift, using the PBS prediction and correction factors.

Notes:

  • This class requires CAMB to be installed and available for the power spectrum calculations.

  • The methods assume that the input masses (M) are in units of solar masses over h (M_sun/h), and distances (R) are in Mpc/h.

  • The halo finder parameters must be provided correctly to ensure accurate mass function and bias calculations.

H(z)[source]

Get the Hubble parameter H at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The Hubble parameter H(z) in km/s/Mpc at redshift z.

Omega_DE(z)[source]

Get the dark energy density parameter Omega_DE at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The dark energy density parameter Omega_DE at redshift z.

Omega_k(z)[source]

Get the curvature density parameter Omega_k at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The curvature density parameter Omega_k at redshift z.

Omega_m(z)[source]

Get the matter density parameter Omega_m at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The matter density parameter Omega_m at redshift z.

Omega_neutrino(z)[source]

Get the massless neutrino density parameter Omega_neutrino at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The massless neutrino density parameter Omega_neutrino at redshift z.

Omega_nu(z)[source]

Get the massive neutrino density parameter Omega_nu at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The neutrino density parameter Omega_nu at redshift z.

Omega_photon(z)[source]

Get the photon density parameter Omega_photon at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The photon density parameter Omega_neutrino at redshift z.

bias(M, z, halo_finder='ROCKSTAR', hmf_model='castro23')[source]

Calculate the corrected linear halo bias.

Parameters:

Mnp.ndarray

Mass of the halos in solar masses over h. Can be a scalar or array.

zfloat

Redshift.

halo_finderstr

The halo finder used to determine the best-fit parameters. Default is ‘ROCKSTAR’.

hmf_modelstr

The hmf model used for the multiplicity function. Default is ‘castro23’.

Returns:

np.ndarray

The corrected linear halo bias.

critical_density(z)[source]

Get the critical density at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The critical density at redshift z in units of 10^10 M_sun / (h^2 Mpc^3).

dlnsigma_dlnR(R, z)[source]

Calculate the derivative of the logarithm of sigma with respect to the logarithm of R.

Parameters:

Rfloat

Radius in Mpc/h.

zfloat

Redshift.

Returns:

float

The derivative dln(sigma)/dln(R).

dndlnM(M, z, halo_finder='ROCKSTAR', model='castro23', PkDE=None)[source]

Calculate the halo mass function, dn/dlnM, which gives the number density of halos per logarithmic mass interval.

Parameters:

Mnp.ndarray

Mass of the halos in solar masses over h.

zfloat

Redshift.

halo_finderstr, optional

Only used for model castro23. Descriptor for the best-fit values for the mass function parameters. Default is ROCKSTAR. Options are: AHF, ROCKSTAR, SUBFIND, and VELOCIraptor.

modelstr, optional

Descriptor for the multiplicity function model. Default is castro23. Options are castro23 and castro25.

PkDEnp.ndarray, optional

Tabulated clustering dark energy power-spectrum at turn-around z_ta=(1+z)/0.5^(2/3)-1. Default is None which means that the modification presented in Castro et al 2025b is not used.

Returns:

np.ndarray

The halo mass function dn/dlnM, with units of number density per logarithmic mass interval.

Raises:

ValueError

If an unsupported model is specified.

Notes:

For the ‘castro25’ model, the halo_finder parameter is ignored.

get_power_spectrum(z=0)[source]

Get the wavenumber and dimensionless power spectrum arrays at a given redshift.

Parameters:

zfloat, optional

Redshift at which to get the power spectrum. Default is 0.

Returns:

tuple

(k, Dk) where k is the array of wavenumbers and Dk is the dimensionless power spectrum.

lagrangian_radius(M)[source]

Calculate the radius of the Lagrangian patch corresponding to a given mass M.

Parameters:

Mfloat

Mass in solar masses over h.

Returns:

float

The Lagrangian radius in Mpc/h.

pbs_bias(M, z, halo_finder='ROCKSTAR', return_variables=False, hmf_model='castro23')[source]

Calculate the PBS bias for halos of mass M at redshift z.

Parameters:

Mnp.ndarray

Mass of the halos in solar masses over h.

zfloat

Redshift.

halo_finderstr

The halo finder used to determine the best-fit parameters. Default is ‘ROCKSTAR’.

return_variablesbool, optional

If True, returns additional intermediate variables used in the calculation. If False, returns only the bias. Default is False.

hmf_modelstr

The hmf model used for the multiplicity function. Default is ‘castro23’.

Returns:

np.ndarray or tuple

If return_variables is False, returns an array of the PBS bias. If return_variables is True, returns a tuple containing: - M : Mass array - R : Lagrangian radius array corresponding to M - v : Peak-height array - dlnsdlnR : Derivative of log sigma with respect to log R - vfv : Multiplicity function values - bias : PBS bias values

peak_height(M, z)[source]

Calculate the peak-height delta_c / sigma(R(M), z).

Parameters:

Mfloat

Mass in solar masses over h.

zfloat

Redshift.

Returns:

float

The peak-height.

set_cosmology(new_params)[source]

Update cosmological parameters and recalculate the power spectrum normalization.

Parameters:

new_paramsdict

Dictionary of new cosmological parameters.

sigma(R, z)[source]

Calculate the RMS fluctuation sigma(R, z).

Parameters:

Rfloat

Radius in Mpc/h.

zfloat

Redshift.

Returns:

float

The RMS fluctuation sigma.

vfv(M, z, return_variables=False, halo_finder='ROCKSTAR', model='castro23')[source]

Calculate the multiplicity function, defined as νf(ν), where ν is the peak-height, as well as related variables for a given mass and redshift.

Parameters:

Mnp.ndarray

Mass of the halos in solar masses over h.

zfloat

Redshift.

return_variablesbool, optional

If True, returns additional intermediate variables used in the calculation. If False, returns only the multiplicity function. Default is False.

halo_finderstr, optional

Only used for model castro23. Descriptor for the best-fit values for the mass function parameters. Default is ROCKSTAR. Options are: AHF, ROCKSTAR, SUBFIND, and VELOCIraptor.

modelstr, optional

Descriptor for the multiplicity function model. Default is castro23. Options are castro23 and castro25.

Returns:

np.ndarray or tuple

If return_variables is False, returns an array of the multiplicity function νf(ν). If return_variables is True, returns a tuple containing: - M : Mass array - R : Lagrangian radius array corresponding to M - ν : Peak-height array - dlnσ/dlnR : Derivative of log sigma with respect to log R - νf(ν) : Multiplicity function values

Raises:

ValueError

If an unsupported model is specified.

Notes:

For the ‘castro25’ model, the halo_finder parameter is ignored, and appropriate best-fit values are used.

wz(z)[source]

Get the dark energy equation of state at a given redshift.

Parameters:

zfloat

Redshift.

Returns:

float

The dark energy equation of state at redshift z.

zta(z)[source]

Get the turnaround redshift for a given collapse redshift z.

Parameters:

zfloat

Redshift.

Returns:

float

EdS solution for the turnaroudn redshift.

bias

bias.py

This module contains functions for calculating the halo bias, including a correction model for the Peak Background Split (PBS) prescription, as described in the provided cosmological study. The correction accounts for dependencies on Omega_m(z), the slope of the power spectrum, and the clustering amplitude S_8.

cctoolkit.bias.bias_correction_PBS(Omega_m_z, dlnsigma_dlnR, S8, A0=1.15, a1=0.0929, b1=0.256, b2=0.173, c1=-0.0372)[source]

Compute the correction factor for the linear halo bias based on the PBS prescription.

Parameters:

Omega_m_zfloat

The matter density parameter at redshift z.

dlnsigma_dlnRfloat

The logarithmic derivative of the variance with respect to the radius.

S8float

The clustering amplitude parameter.

Returns:

float

The correction factor for the halo bias.

References

  • [arxiv:2208.02174](https://arxiv.org/pdf/2311.01465): “Euclid preparation. XXIV. Calibration of the halo mass function in Λ(ν)CDM cosmologies”, Castro et al., 2023.

  • [arxiv:2409.01877](https://arxiv.org/pdf/2409.01877): “Euclid preparation. Calibration of the linear halo bias in Λ(ν)CDM cosmologies”, Castro et al., 2024.

cctoolkit.bias.corrected_bias(b_PBS, Omega_m_z, dlnsigma_dlnR, S8)[source]

Calculate the corrected linear halo bias using the PBS prediction and the correction factor.

Parameters:

b_PBSfloat

The linear halo bias predicted by the PBS model.

Omega_m_zfloat

The matter density parameter at redshift z.

dlnsigma_dlnRfloat

The logarithmic derivative of the variance with respect to the radius.

S8float

The clustering amplitude parameter.

Returns:

float

The corrected linear halo bias.

References

  • [arxiv:2208.02174](https://arxiv.org/pdf/2311.01465): “Euclid preparation. XXIV. Calibration of the halo mass function in Λ(ν)CDM cosmologies”, Castro et al., 2023.

  • [arxiv:2409.01877](https://arxiv.org/pdf/2409.01877): “Euclid preparation. Calibration of the linear halo bias in Λ(ν)CDM cosmologies”, Castro et al., in prep.

hmf

hmf.py

This module contains functions and data related to the Halo Mass Function (HMF), including the multiplicity function and best-fit parameters for different halo finders.

cctoolkit.hmf.multiplicity_function(peak_height, dlns_dlnR, Omega_m_z, best_fit_values, model='castro23', Omega_de_zta=None, w_de_zta=None)[source]

Handler for the multiplicity function.

Parameters:

peak_heightfloat

The peak height (ν).

dlns_dlnRfloat

The logarithmic slope of the variance with respect to the radius.

Omega_m_zfloat

The matter density parameter at redshift z.

best_fit_valuesdict

Dictionary containing the best-fit parameters for the chosen halo finder.

modelstr, optional

The model to use for the multiplicity function. Options are ‘castro23’ (default) and ‘castro25’.

Omega_de_ztafloat, optional

The dark energy density parameter at turn-around. Required if model is ‘castro25’.

w_de_ztafloat, optional

The dark energy equation of state parameter at turn-around. Required if model is ‘castro25’.

Returns:

float

The value of the multiplicity function.

Raises:

ValueError

If required parameters are not provided for the selected model.

References

  • [arXiv:2208.02174](https://arxiv.org/pdf/2208.02174): “Euclid preparation. XXIV. Calibration of the halo mass function in Λ(ν)CDM cosmologies”, Castro et al., 2023.

  • [arXiv:2504.07608](https://arxiv.org/pdf/2504.07608): “DUCA: Dynamic Universe Cosmological Analysis. The halo mass function in dynamic dark energy cosmologies”, Castro et al., 2025.

cctoolkit.hmf.multiplicity_function_castro23(peak_height, dlns_dlnR, Omega_m_z, best_fit_values)[source]

Calculate the multiplicity function using the Castro et al. 2023 model.

Parameters:

peak_heightfloat

The peak height (ν).

dlns_dlnRfloat

The logarithmic slope of the variance with respect to the radius.

Omega_m_zfloat

The matter density parameter at redshift z.

best_fit_valuesdict

Dictionary containing the best-fit parameters for the chosen halo finder.

Returns:

float

The value of the multiplicity function.

References

  • [arXiv:2208.02174](https://arxiv.org/pdf/2208.02174): “Euclid preparation. XXIV. Calibration of the halo mass function in Λ(ν)CDM cosmologies”, Castro et al., 2023.

cctoolkit.hmf.multiplicity_function_castro25(peak_height, dlns_dlnR, Omega_m_z, Omega_de_zta, w_de_zta, best_fit_values)[source]

Calculate the multiplicity function using the Castro et al. 2025 model, which accounts for dynamic dark energy cosmologies.

Parameters:

peak_heightfloat

The peak height (ν).

dlns_dlnRfloat

The logarithmic slope of the variance with respect to the radius.

Omega_m_zfloat

The matter density parameter at redshift z.

Omega_de_ztafloat

The dark energy density parameter at turn-around.

w_de_ztafloat

The dark energy equation of state parameter at turn-around.

best_fit_valuesdict

Dictionary containing the best-fit parameters for the chosen halo finder.

Returns:

float

The value of the multiplicity function.

References

  • [arXiv:2504.07608](https://arxiv.org/pdf/2504.07608): “DUCA: Dynamic Universe Cosmological Analysis. I. The halo mass function in dynamic dark energy cosmologies”, Castro et al., 2025.

baryons

baryons.py

This module implements the model for the impact of baryonic effects on the masses of clusters and groups presented in Castro et al. 2024 (https://inspirehep.net/literature/2718844).

cctoolkit.baryons.baryon_fraction_magneticum(M_vir, z)[source]

Mean baryon fraction contained inside the virial radius for Magneticum simulations.

Parameters

M_virfloat

Virial mass.

zfloat

Redshift.

Returns

fb_fidfloat

The mean baryon fraction inside the virial radius for Magneticum clusters.

cctoolkit.baryons.baryon_fraction_tng300(M_vir, z)[source]

Mean baryon fraction contained inside the virial R for TNG300 simulations.

Parameters

M_virfloat

Virial mass.

zfloat

Redshift.

Returns

fb_fidfloat

The mean baryon fraction inside the virial radius for TNG300 clusters.

cctoolkit.baryons.compute_dmo_mass(M_vir_hydro, z, fb_cosmic, relation='magneticum')[source]

Compute the SO equivalent DMO mass and overdensity Delta (virial units) to the virial object definition from hydrodynamic simulations.

Parameters

M_vir_hydrofloat

Virial mass from hydrodynamic simulations.

zfloat

Redshift.

fb_cosmicfloat

Cosmic baryon fraction.

relationstr or function, optional

The baryon fraction relation to use, either ‘magneticum’ or ‘tng300’ or a custom function with signature (M, z). Default is ‘magneticum’.

Returns

M_delta_dmofloat

The equivalent DMO mass at the Delta overdensity.

Deltafloat

The overdensity of equivalence in virial units.

Raises

ValueError

If an invalid relation is specified.

cctoolkit.baryons.compute_rec_mass(cosmo, M_Delta_dmo, Delta, z, c=<function concentration_duffy_et_al_2008>)[source]

Reconstruct the dark matter-only (DMO) mass of a hydro object using the model presented in [arxiv:2311.01465](https://arxiv.org/pdf/2311.01465).

This function calculates the reconstructed mass of a hydro object, denoted as ( M_{Delta, ext{dmo}} ), using the provided concentration-mass relation and cosmological parameters.

Parameters

cosmoobject

Cosmology object containing cosmological parameters and methods to compute relevant quantities, such as the critical density and Omega_m(z).

M_Delta_dmofloat

Mass of the dark matter-only object within radius defined by (Delta) times the critical density, in solar masses.

Deltafloat

Overdensity parameter, defining the radius within which the mass ( M_{Delta, ext{dmo}} ) is calculated.

zfloat

Redshift at which the calculation is performed.

cfunction, optional

Function to compute the concentration parameter given mass and redshift. Default is concentration_duffy_et_al_2008.

Returns

float

Reconstructed virial mass of the object at redshift z.

Notes

The function first computes the virial overdensity parameter, then iteratively solves for the virial radius ( R_{ ext{vir}} ) using the auxiliary NFW function. The virial mass is finally obtained using the computed radius and density.

References

  • [arxiv:2311.01465](https://arxiv.org/pdf/2311.01465): “Euclid preparation XXXIX. The effect of baryons on the halo mass function”, Castro et al., 2024.

cctoolkit.baryons.concentration_duffy_et_al_2008(M, z)[source]

Calculate the concentration-mass relation according to Duffy et al. (2008) for virial masses.

This function implements the formula provided by Duffy et al. (2008) for determining the concentration parameter as a function of mass and redshift. The formula is given by:

\[c(M, z) = A_{vir} \left(\frac{M}{M_{piv}}\right)^{B_{vir}} (1+z)^{C_{vir}}\]

where the parameters are defined as follows: - \(A_{vir} = 7.85\) - \(B_{vir} = -0.081\) - \(C_{vir} = -0.71\) - \(M_{piv} = 2 \times 10^{12} M_\odot\)

Parameters

Mfloat

Virial mass of the halo in solar masses.

zfloat

Redshift at which the concentration is calculated.

Returns

float

Concentration parameter for the given mass and redshift.

References

Duffy, A. R., et al. (2008). “Dark matter halo concentrations in the Wilkinson Microwave Anisotropy Probe year 5 cosmology.” Monthly Notices of the Royal Astronomical Society, 390(2), L64-L68.

utils

utils.py

This module can contain utility functions or constants that may be useful across the package.

cctoolkit.utils.compute_sigma8_norm(k, Dk, sigma8)[source]

Normalize the dimensionless power spectrum to the given sigma8.

Parameters:

knp.ndarray

Wavenumber array.

Dknp.ndarray

Dimensionless power spectrum array.

sigma8float

The value of sigma8.

Returns:

float

The normalization factor for the power spectrum.

cctoolkit.utils.d_s_given_Dk(r, k, Dk)[source]

Compute the derivative of sigma_r with respect to r.

Parameters:

rfloat

Radius of the spherical top-hat window.

knp.ndarray

Wavenumber array.

Dknp.ndarray

Dimensionless power spectrum array.

Returns:

float

The derivative of sigma_r with respect to r.

cctoolkit.utils.d_s_rsq_integrand_given_Dk(r, k, Dk)[source]

Integrand for the derivative of sigma_r^2 with respect to r.

Parameters:

rfloat

Radius of the spherical top-hat window.

knp.ndarray

Wavenumber array.

Dknp.ndarray

Dimensionless power spectrum array.

Returns:

np.ndarray

The integrand for the derivative of sigma_r^2 with respect to r.

cctoolkit.utils.delta_c(Om, z=None)[source]

Calculate the critical density contrast (delta_c) using eq. A.6 of Kitayama & Suto.

Parameters:

Omfloat or callable

Matter density parameter. If callable, it should return Om as a function of z.

zfloat or array-like, optional

Redshift at which to evaluate Om if Om is a function.

Returns:

float

The critical density contrast delta_c.

cctoolkit.utils.dw_tophat_sq_dr(r, k)[source]

Derivative of the square of the Fourier Transform of the top-hat function with respect to r.

Parameters:

rfloat

Radius of the spherical top-hat window.

kfloat or np.ndarray

Wavenumber(s).

Returns:

float or np.ndarray

The derivative with respect to r.

cctoolkit.utils.s_rsq_integrand_given_Dk(r, k, Dk)[source]

Integrand for the calculation of sigma_r^2.

Parameters:

rfloat

Radius of the spherical top-hat window.

knp.ndarray

Wavenumber array.

Dknp.ndarray

Dimensionless power spectrum array.

Returns:

np.ndarray

The integrand for sigma_r^2.

cctoolkit.utils.ssq_given_Dk(r, k, Dk)[source]

Compute sigma_r^2 for a given radius and power spectrum.

Parameters:

rfloat

Radius of the spherical top-hat window.

knp.ndarray

Wavenumber array.

Dknp.ndarray

Dimensionless power spectrum array.

Returns:

float

The value of sigma_r^2.

cctoolkit.utils.virial_Delta(Om, z=None)[source]

Calculate the virial overdensity (Delta_vir) in critical density unit using Bryan and Norman fit.

Parameters:

Omfloat or callable

Matter density parameter. If callable, it should return Om as a function of z.

zfloat or array-like, optional

Redshift at which to evaluate Om if Om is a function.

Returns:

float

The virial overdensity.

cctoolkit.utils.w_tophat(r, k)[source]

Fourier Transform of the spherical top-hat function.

Parameters:

rfloat

Radius of the spherical top-hat window.

kfloat or np.ndarray

Wavenumber(s).

Returns:

float or np.ndarray

The Fourier Transform of the top-hat function.