PyCosmo package

Submodules

PyCosmo.Background module

class PyCosmo.Background.Background(cosmo, params, rec, wrapper)[source]

Bases: object

The background calculations are based on a Friedmann-Robertson-Walker model for which the evolution is governed by the Friedmann equation.

H(a)[source]

Calculates the Hubble parameter for a given scale factor by calling the expression from CosmologyCore_model.py.

Parameters:

a – scale factor [1]

Returns:

Hubble parameter \(H(a) [km/s/Mpc]\).

Example:

cosmo.background.H(a)
cs(a)[source]

Returns the photon-baryon fluid sound speed [1] from the Recombination module

Parameters:

a – scale factor [1]

Returns:

cs: photon-baryon sound speed [1]

Example:

cosmo.background.cs(a)
dist_ang_a(a=1.0)[source]

Calculates the angular-diameter distance to a given scale factor.

Parameters:

a – scale factor [1]

Returns:

Angular diameter distance, \(D_A(a)\), in \([Mpc]\).

Example:

cosmo.background.dist_ang_a(a)
dist_lum_a(a=1.0)[source]

Calculates the luminosity distance to a given scale factor.

Parameters:

a – scale factor [1]

Returns:

Luminosity distance, \(D_L(a)\), in \([Mpc]\).

Example:

cosmo.background.dist_lum_a(a)
dist_rad_a(a=1.0)[source]

Calculates the radial comoving distance, also known as the comoving radius.

Parameters:

a – scale factor [1]

Returns:

Radial comoving distance, \(\chi(a)\), in \([Mpc]\).

Example:

cosmo.background.dist_rad_a(a)
dist_trans_a(a=1.0)[source]

Calculates the transverse comoving distance, also known as comoving angular-diameter distance.

Parameters:

a – scale factor [1]

Returns:

Transverse comoving distance, \(r(\chi)\), in \([Mpc]\).

Example:

cosmo.background.dist_trans_a(a)
eta(a)[source]
eta_to_a(eta)[source]
g_a(a=1.0)[source]

Calculates the visibility function for a given scale factor.

Parameters:

a – scale factor [1]

Returns:

visibility function [1/Mpc]

Example:

cosmo.background.g_a(a)
r_s()[source]

Calculates the sound horizon at drag epoch.

Returns:

sound horizon [\(\mathrm{Mpc} / h\)]

tau(a=1.0)[source]

Calculates the optical depth of the photon-baryon fluid, by integrating taudot from the Recombination module.

Parameters:

a – scale factor [1]

Returns:

tau: optical depth [1]

Example:

cosmo.background.tau(a)
tau_b(a=1.0)[source]

Calculates the optical depth of the baryon fluid, by integrating taudot from the Recombination module, weighted with R, the baryon-photon fraction.

Parameters:

a – scale factor [1]

Returns:

tau: optical depth [1]

Example:

cosmo.background.tau(a)
taudot(a=1.0)[source]

Calculates \(\dot{\tau} = \frac{d\tau}{d\eta}\), where \(\tau\) is the optical depth and \(\eta\) is the conformal time.

Parameters:

a – scale factor [1]

Returns:

\(\dot{\tau} = \frac{d\tau}{d\eta}\) in units of [\(h \mathrm{Mpc}^{-1}\)]

Example:

cosmo.background.taudot(a)

PyCosmo.Cosmo module

class PyCosmo.Cosmo.Cosmo(*, model_config=None)[source]

Bases: object

All of main functionalities of PyCosmo are managed by the main (Cosmo) class, which links to the other classes where the bulk of the calculations is performed.

load_params(paramfile=None)[source]

Loads parameter file and additional user-defined parameters to initialize the Cosmo object

Parameters:

paramfile – path to parameter file

print_params(inc_derived=True, inc_internal=True, file=None)[source]

Prints the parameters of PyCosmo instance.

Parameters:
  • inc_derived – prints derived parameters [True or False]

  • inc_internal – prints internal parameterss (e.g. for lin_pert) [True or False]

static recompile_from_cli()[source]
reset_wrapper_globals()[source]
set(**kwargs)[source]

Changes the configuration parameters of the PyCosmo instance. This will also recalculate the derived quantities.

Parameters:

kwargs – keyword arguments for the input parameters to change.

Warning

Always change parameters of an instance of PyCosmo using the set function. Do not change variables directly!

Returns:

If no keywords are passed then a list of options is printed

Example:

cosmo.set(omega_m=0.23)
class PyCosmo.Cosmo.FakeProjection[source]

Bases: object

PyCosmo.Cosmo.make_pickable(checks)[source]
PyCosmo.Cosmo.trace_call(frame, event, arg)[source]

PyCosmo.LinearPerturbationApprox module

class PyCosmo.LinearPerturbationApprox.LinearPerturbationApprox(*a, **kw)[source]

Bases: LinearPerturbationBase

Class created to manage fitting functions for the computation of the linear matter power spectrum \(P_{lin}(k)\). The fitting functions provide transfer functions that are then used to compute the power spectrum.

The different fitting functions can be selected using the set function:

cosmo.set(pk_type = option)

where option can be one the following keywords:

For developers: in order to compare the codes \(\textsf{PyCosmo}\) and \(\texttt{CCL}\) in terms of linear matter power spectrum computed with the BBKS fitting function, the user should choose a routine which is optimized for this purpose. This further option can be selected as:

cosmo.set(pk_type = 'BBKS_CCL')

where BBKS_CCL follows the implementation in the CCL code .

growth_a(a, k=None, norm=0, diag_only=True)[source]

Returns the linear growth factor by calling a function that integrates the growth differential equation. It returns a vector of len(a) if k is None or a number, or if diag_only is set to True. Otherwise it returns a matrix with dimension (len(k), len(a)).

Parameters:
  • a – scale factor [1]

  • k – wavenumber k, necessary for massive neutrino models [\(h Mpc^{-1}\)]

  • norm – normalisation scheme: norm = 0: D(a=1)=1 (default), norm = 1: D(a)=a in matter era

  • diag_only – if set to False, the growth factor is repeated and reshaped to match the shape of k

Returns:

growth factor [1]

max_redshift(k)[source]

computes max redshift for which this model is applicable.

Parameters:

k – wavenumber k, introduced in the abstract base class since this might be needed for some models. [\(h Mpc^{-1}\)]

Returns:

redshift

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Computes the linear total matter power spectrum, \(P_{lin}(k)\), using a choice of fitting functions.

Parameters:
  • a – scale factor [1]

  • k – wavenumber \([Mpc]^{-1}\)

  • diag_only – if set to True: compute powerspectrum for pairs \(a_i, k_i\), else consider all combinations \(a_i, k_j\)

Returns:

Linear matter power spectrum, \(P_{lin}(k)\), in \([Mpc]^3\).

Example:

cosmo.set(pk_type = option)
cosmo.lin_pert.powerspec_a_k(a,k)

where option can be set to one of the fitting functions (EH or BBKS).

print_params()[source]

Print the parameters related to the chosen linear fitting function (EH or BBKS).

Example:

cosmo.lin_pert.print_params()
transfer_k(k)[source]

Computes the linear matter transfer function using a choice of the currently available fitting functions.

Parameters:

k – wavenumber \([Mpc^{-1}]\)

Returns:

Matter transfer function \(T(k)\) in \(Mpc^{3/2}\).

PyCosmo.LinearPerturbationBase module

class PyCosmo.LinearPerturbationBase.LinearPerturbationBase(*a, **kw)[source]

Bases: object

This is the parent class of LinearPerturbationApprox and LinearPerturbationBoltzmann. It provides some abstract methods for growth factors, transfer functions and the linear power spectrum. It also calculates \(\sigma_8\) and \(\sigma_r\).

abstract growth_a(a=1.0, k=None, norm=0, diag_only=True)[source]
abstract max_redshift(k)[source]
min_a(k)[source]
abstract powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]
sigma8(a=1.0, k_grid=None)[source]

Computes sigma8, the rms density contrast fluctuation smoothed with a top hat of radius 8 \(h^{-1} Mpc\). This routine is also used for the normalization of the power spectrum, when pk_type=``sigma8`` is chosen.

Parameters:

a – scale factor [1]

Returns:

sigma8 [1]

sigma_r(r=8.0, a=1.0)[source]

Calculates the rms of the density field on a given scale \(r\). It is a generalization of the function sigma8().

Parameters:

r – scale radius (default set to \(8 h^{-1} Mpc\))

Returns:

\(\sigma_r\) [1]

abstract transfer_k(k)[source]

PyCosmo.LinearPerturbationTable module

class PyCosmo.LinearPerturbationTable.LinearPerturbationTable(*a, **kw)[source]

Bases: LinearPerturbationBase

This class represents the interface to the tables class, which builds tabulated data for the costly calculations. Using this tabulated data this class constructs interpolation functions which can be called instead of the full function in order to increase speed.

growth_a(a=1.0, k=None, norm=0, verbose=False)[source]

Returns an interpolated growth factor computed from a predefined grid with interpolation error smaller than tol which has been set with enrich_params.

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Returns an interpolated linear matter power spectrum computed from a predefined grid with interpolation error smaller than tol which has been set with enrich_params. :param: a: scale factor [1] :param: k: wavenumber [Mpc^-1] :return: P(k, a): interpolated matter power spectrum P(k) [Mpc^3]

transfer_k(k)[source]

Returns an interpolated transfer function computed from a predefined grid with interpolation error smaller than tol which has been set with enrich_params.

PyCosmo.NonLinearPerturbationHaloFit module

class PyCosmo.NonLinearPerturbationHaloFit.NonLinearPerturbationHaloFit(*a, **kw)[source]

Bases: NonLinearPerturbationBase

The class computes the non-linear matter power spectrum, \(P_{nl}(k)\), with a choice of analytic fitting funtions.

Warning

NonLinearPerturbationHalofit not supported by the Boltzmann solver yet. Calling any function from this module while using the Boltzmann solver (pk_type = ‘boltz’) will return AttributeError: ‘NoneType’ object has no attribute ‘function’

k_nl(a)[source]

Computes the non-linear wave number \(k_{nl}\) as a function of the scale factor a.

Parameters:

a – scale factor [1]

Returns:

k_nl: non-linear wave number [\(\mathrm{Mpc}^{-1}\)]

max_redshift(k)[source]

computes max redshift for which this model is applicable.

Parameters:

k – wavenumber k, might be needed for some underlying linear perturbation models. [\(h Mpc^{-1}\)]

Returns:

redshift

min_a(k)[source]
powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Returns the nonlinear matter power spectrum using a choice of fitting functions. Currently both the Halofit fitting function and its revision by Takahashi et al., 2012 are supported. Those can be selected using the set function:

cosmo.set(pk_nonlin_type = option)

where option can be one the following keywords:

After selecting the fitting function, the non-linear matter power spectrum, \(P_{nl}(k)\), can be computed as follows:

cosmo.nonlin_pert.powerspec_a_k(a,k)
Parameters:
  • a – scale factor [1]

  • k – wavenumber \([Mpc^{-1}]\)

Returns:

Nonlinear matter power spectrum, \(P_{nl}(k)\), in \([Mpc^3]\).

PyCosmo.NonLinearPerturbationMead module

class PyCosmo.NonLinearPerturbationMead.NonLinearPerturbationMead(*a, **kw)[source]

Bases: NonLinearPerturbationBase

The class incorporates the implementation of the \(\texttt{HMCode}\) as described in Mead et al. 2015 and Mead et al. 2016 .

Warning

NonLinearPerturbationMead not supported by the Boltzmann solver yet. Calling any function from this module while using the Boltzmann solver (pk_type = ‘boltz’) will return AttributeError: ‘NoneType’ object has no attribute ‘function’

T(x)[source]

Computes the Fourier Transform of a 3D spherical top-hat function.

Parameters:

x – dimensionless parameter (usually equal to \(kr\)) [1]

Returns:

Fourier Transform of a 3D spherical top-hat function.

T_deriv(x)[source]

Analytical computation of the derivative of the Fourier Transform of a 3D spherical top-hat function.

Parameters:

x – dimensionless parameter (usually equal to \(kr\)) [1]

Returns:

Derivative of the Fourier Transform of a 3D spherical top-hat function.

cm(m_msun, a)[source]

Concentration-mass function of Dark Matter Haloes. The implemented fitting function is from Bullock et al, 2001, as described in Eq. (14) of Mead et al, 2015.

Parameters:
  • a – scale factor a [1]

  • m_msun – halo mass in solar masses [Msun]

Returns:

Concentration \(c(a, M)\) of a halo of mass \(m_{msun}\) at scale factor \(a\).

delta_c(a)[source]

Computation of the linear collapse threshold, \(\delta_c\).

Following Mead et al., 2015,, \(\delta_c\) is treated as a fitting parameter (read its fitted value in Table 2 of the paper). The current implementation in \(\textsf{PyCosmo}\) follows the updated formula reported in Eq.(8) of Mead et al., 2016,, where the original prescription from 2015 is augmented by the fitting formula proposed in Eq. C28 of Nakamura & Suto, 1997.

Parameters:

a – scale factor [no dimensions]

Returns:

Linear collapse threshold, \(\delta_c\).

f(nu)[source]

Multiplicity function \(f(\nu)\) as it appears in the calculation of the Halo Mass Function. The available fitting functions are Press & Schechter, 1974, Sheth & Tormen, 1999 , Tinker et al., 2010, and Watson et al., 2013 .

Parameters:

nu – array of \(\nu\) values as in \(\nu = \delta_c / \sigma(M)\) (see eq. 17 in Mead et al., 2015)

Returns:

Multiplicity function \(f(\nu)\) .

Example of setting the multiplicity function:

cosmo.set(multiplicity_fnct = option)

where option can be ‘PS’ for Press & Schechter (1974), ‘ST’ for Sheth & Tormen (1999), ‘Ti’ for Tinker et al., 2010 or ‘Wa’ for Watson et al., 2013.

mass2radius(m_msun)[source]

Converts mass of a sphere in solar masses [Msun] to corresponding radius in [Mpc], assuming homogeneous density corresponding to the critical matter density of the universe (redshift z=0).

Parameters:

m_msun – mass of sphere in solar masses [Msun]

Returns:

Radius of the sphere in [Mpc].

max_redshift(k)[source]

computes max redshift for which this model is applicable.

Parameters:

k – wavenumber k, might be needed for some underlying linear perturbation models. [\(h Mpc^{-1}\)]

Returns:

redshift

min_a(k)[source]
neff(a)[source]

Calculates the effective spectral index of the linear power spectrum at the collapse scale, according to the implementation in the HMCode.

Parameters:

a – scale factor [1]

Returns:

Effective spectral index of the linear power spectrum, \(n_{eff}\), at the collapse scale.

nu2mass(nu, a)[source]

Extracts the halo mass from Eq.(17) in Mead et al., 2015, \(\nu = \delta_c / \sigma\), by converting the \(\nu\)-array into a mass-array by backwards interpolation.

Parameters:

nu\(\nu\)-array as in \(\nu = \delta_c / \sigma\) [1]

Returns:

Mass in solar masses [Msun].

pk_1h(k, a, nu_range)[source]

One-Halo term of the non-linear matter power spectrum as defined in Eq.(8) of Mead et al, 2015.

Parameters:
  • k – wavelength \([Mpc^{-1}]\)

  • a – scale factor [1]

  • nu_range – nu_range-array used for the integration (default set to full lookup table)

  • mode – integration scheme. The values 0, 1 and 2 refer to numpy integration, cython integration and cython with adaptive integration, respectively.

Returns:

One-Halo power spectrum \(P_{1h}(k) \ [Mpc^{3}]\).

pk_2h(k, a)[source]

Two-Halo term of the non-linear matter power spectrum as defined in Eq.(10) of Mead et al, 2015.

Parameters:
  • k – wavelength \([Mpc^{-1}]\)

  • a – scale factor [1]

Returns:

Two-Halo power spectrum \(P_{2h}(k) \ [Mpc^{3}]\).

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Calculates the non-linear matter power spectrum, using the mead model, as described in Mead et al., 2015 and Mead et al., 2016, consisting of the superposition of the pk_1h() and pk_2h() terms.

Parameters:
  • k – wavelength \([Mpc]^{-1}\)

  • a – scale factor [1]

  • nu_range – nu-array used for the integration [1]

  • m_min_msun – minimum halo mass for the integration range [Msun]

  • m_max_msun – maximum halo mass for the integration range [Msun]

Returns:

Halo Model power spectrum, \(P_{nl}(k)\), in \([Mpc^{3}]\).

Example:

cosmo.set(pk_nonlin_type = 'mead')
cosmo.nonlin_pert.powerspec_a_k(a,k)
print_params()[source]

Prints the cosmological setup and the parameters used for the computation of the non-linear matter power spectrum with the Mead et al., 2015 model.

Example:

cosmo.set(pk_nonlin_type = 'mead')
cosmo.nonlin_pert.print_params()
radius2mass(r_mpc)[source]

Converts radius of a sphere in [Mpc] to corresponding mass in solar masses [Msun], assuming homogeneous density corresponding to the critical matter density of the universe (redshift z=0).

Parameters:

r_mpc – radius of sphere [Mpc]

Returns:

Mass of sphere in solar masses [Msun].

rvir(mvir_msun, a)[source]

Calculates the virial radius corresponding to a dark matter halo of mass m_msun in a flat cosmology. See Bryan & Norman, 1998 for more details.

Parameters:

mvir_msun – Halo mass in solar masses [Msun]

Returns:

Virial radius in \([Mpc]\).

sigma(m_msun, a)[source]

Calculates \(\sigma(M, a)\), the RMS of the density field at a given mass.

Parameters:
  • m_msun – mass in solar masses at which sigma is evaluated [Msun]

  • a – scale factor [1]

Returns:

\(\sigma(M, a)\) as the RMS of the density field.

sigma8_a(a)[source]

Compute \(\sigma_8\), the rms density contrast fluctuation, smoothed with a top hat of radius 8 \(h^{-1}\mathrm{Mpc}\) at scale factor a.

Parameters:

a – scale factor [1]

Returns:

\(\sigma_8\) at the desired scale factor.

sigma_d(a, R=100.0)[source]

Computes \(\sigma^2_D(a)\) as defined in Eq. B5 of Mead et al, 2015

Parameters:

a – scale factor [1]

Returns:

\(\sigma^2_D(a)\) at the desired scale factor.

sigma_v_a(a)[source]

Computes \(\sigma^2_V(a)\) as defined in Eq. (22) of Mead et al, 2015 .

Parameters:

a – scale factor [1]

Returns:

\(\sigma^2_v(a)\) at the desired scale factor.

PyCosmo.NonLinearPerturbationMead.integral_halo(k, m_msun, nu_range, rv_mpc, c, a, f, eta)[source]

PyCosmo.NonLinearPerturbationTable module

class PyCosmo.NonLinearPerturbationTable.NonLinearPerturbationTable(*a, **kw)[source]

Bases: NonLinearPerturbationBase

This class represents the interface to the tables class, which builds tabulated data for the costly calculations. Using this tabulated data this class constructs interpolation functions which can be called instead of the full function in order to increase speed.

k_nl(a=1.0)[source]

Returns an interpolated non-linear wave number computed from a predefined grid with interpolation error smaller than tol which has been set with enrich_params.

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Returns an interpolated nonlinear matter power spectrum computed from a predefined grid with interpolation error smaller than tol which has been set with enrich_params. :param: a: scale factor [1] :param: k: wavenumber [Mpc^-1] :return: P(k, a): interpolated nonlinear matter power spectrum P(k) [Mpc^3]

PyCosmo.NonLinearPerturbation_HI module

class PyCosmo.NonLinearPerturbation_HI.NonLinearPerturbation_HI(*a, **kw)[source]

Bases: NonLinearPerturbation_HaloModel

The class incorporates the implementation of the Halo Model for neutral hydrogen (HI) as described in Padmanabhan et al., 2017.

This can be set as:

cosmo.set(pk_nonlin_type='HI')

Warning

NonLinearPerturbation_HI not supported by the Boltzmann solver yet. Calling any function from this module while using the Boltzmann solver (pk_type = ‘boltz’) will return AttributeError: ‘NoneType’ object has no attribute ‘function’

c_hi(m_msun, a)[source]

Concentration-mass function of HI halos.

Parameters:
  • m_msun – Halo mass in solar masses \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor a, scalar or 1d array

Returns:

Concentration \(c_\text{HI}(M, a)\) of a halo mass \([M_{\odot}]\) at scale factor \(a\)

mean_hi_temp(a)[source]

Calculates the mean HI brightness temperature \(T_\text{b}(a)\) in \([K]\) as a function of the scale factor.

Parameters:

a – Scale factor, scalar or 1d array

Returns:

Mean HI brightness temperature, in \([K]\)

mhi_of_m(m_msun, a, diag_only=False)[source]

HI-halo mass relation: Calculates the amount of HI mass in a dark matter halo of mass \(M\).

Parameters:
  • m_msun – Halo mass in solar masses \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor a, scalar or 1d array

Returns:

HI mass in solar masses \([M_{\odot}]\)

pk_HI_1h(k, a)[source]

One-halo term \(P_\text{HI,1h}(k, a)\) of the non-linear HI power spectrum as described in the Halo Model for HI.

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

One-Halo power spectrum \(P_\text{HI,1h}(k, a) \ [Mpc^{3}]\)

Example of setting the multiplicity function, linear halo bias, halo profile, and the assumed minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation:

cosmo.set(multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,
          halo_profile=option_p,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

where the possible option_mf are given in the documentation to \(f(\nu, a)\) and option_b in the documentation to the linear halo bias in the NonLinearPerturbation_HaloModel module. option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

pk_HI_2h(k, a)[source]

Two-halo term \(P_\text{HI,2h}(k, a)\) of the non-linear power spectrum as described in the Halo Model for HI.

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

Two-Halo power spectrum \(P_\text{HI,2h}(k, a) \ [Mpc^{3}]\)

Example of setting the multiplicity function, linear halo bias, halo profile, and the assumed minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation:

cosmo.set(multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,
          halo_profile=option_p,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

where the possible option_mf are given in the documentation to \(f(\nu, a)\) and option_b in the documentation to the linear halo bias in the NonLinearPerturbation_HaloModel module. option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Calculates the non-linear matter power spectrum using the Halo Model, described e.g. in Mo et al. “Galaxy Formation and Evolution” (2010) chapter 7.6, as the sum of the pk_1h() and pk_2h() terms.

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

Halo Model power spectrum, \(P_\text{nl}(k, a)\), in \([Mpc^{3}]\)

Example on how to use the Halo model, set the multiplicity function, linear halo bias, halo profile, and the minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation, and then calculate the power spectrum:

cosmo.set(pk_nonlin_type='HaloModel',
          multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,
          halo_profile=option_p,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

pk = cosmo.nonlin_pert.powerspec_a_k(a,k)

where the possible option_mf are given in the documentation to \(f(\nu, a)\), option_b in the documentation to the linear halo bias and option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

powerspec_a_k_HI(a=1.0, k=0.1, diag_only=False)[source]

Calculates the non-linear HI power spectrum using the Halo Model for HI, as the sum of the pk_HI_1h() and pk_HI_2h() terms.

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

HI Halo Model power spectrum, \(P_\text{HI,nl}(k, a)\), in \([Mpc^{3}]\)

Example on how to use the HI Halo model, set the multiplicity function, linear halo bias, halo profile, and the minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation, and then calculate the power spectrum:

cosmo.set(pk_nonlin_type='HI',
          multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,
          halo_profile=option_p,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

pk = cosmo.nonlin_pert.powerspec_a_k_HI(a,k)

where the possible option_mf are given in the documentation to \(f(\nu, a)\) and option_b in the documentation to the linear halo bias in the NonLinearPerturbation_HaloModel module. option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

rho_HI(a)[source]

Calculates the mean HI density as a function of the scale factor.

Parameters:

a – Scale factor, scalar or 1d array

Returns:

Mean HI density, in \([M_{\odot} \ Mpc^{-3}]\)

Example of setting the multiplicity function and the assumed minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation:

cosmo.set(multiplicity_fnct=option,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

where option can be ‘PS’ for Press & Schechter (1974), ‘ST’ for Sheth & Tormen (1999), ‘Ti’ for Tinker et al., 2010, or ‘Wa’ for Watson et al., 2013. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

PyCosmo.NonLinearPerturbation_HaloModel module

class PyCosmo.NonLinearPerturbation_HaloModel.NonLinearPerturbation_HaloModel(*a, **kw)[source]

Bases: NonLinearPerturbationBase

The class incorporates the implementation of the Halo Model.

This can be set as:

cosmo.set(pk_nonlin_type='HaloModel')

Warning

NonLinearPerturbation_HaloModel not supported by the Boltzmann solver yet. Calling any function from this module while using the Boltzmann solver (pk_type=’boltz’) will return AttributeError: ‘NoneType’ object has no attribute ‘function’

T(x)[source]

Computes the Fourier Transform of a 3D spherical top-hat function.

Parameters:

x – Dimensionless parameter (usually equal to \(kr\)) [1]

Returns:

Fourier Transform of a 3D spherical top-hat function

cm(m_msun, a)[source]

Concentration-mass function of dark matter haloes. The implemented fitting function is from Bullock et al., 2001.

Parameters:
  • m_msun – Halo mass in solar masses \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor a, scalar or 1d array

Returns:

Concentration \(c(M, a)\) of a halo mass \([M_{\odot}]\) at scale factor \(a\)

delta_c(a)[source]

Computation of the cosmology dependent linear collapse threshold, \(\delta_c\). For the mass functions of Press & Schechter, 1974 , Sheth & Tormen, 1999 , and Watson et al., 2013 , we use the cosmology dependent expression from Nakamura & Suto, 1997 , or Mo & White, 2010 . For the mass function of Tinker et al., 2010 , we assume \(\delta_c=1.686\).

Parameters:

a – Scale factor, scalar or 1d array

Returns:

Linear collapse threshold \(\delta_c\)

delta_vir(a)[source]

Calculates the mean overdensity of a virialised dark matter halo for different cosmologies. We use equation (6) of Bryan & Norman, 1998 modified by a \(1/\Omega_m\) term, because we work with respect to the matter density, rather than the critical density.

Parameters:

a – Scale factor, scalar or 1d array

Returns:

Mean overdensity of a dark matter halo \(\Delta_{vir}\)

dn_dlnm_of_m(m_msun, a, dm=100000.0)[source]

Calculates the logarithmic mass function \(\frac{dn}{dln M}\) as a function of the halo mass. The functional form depends on the selected multiplicity function \(f(\nu, a)\), see the corresponding documentation.

Parameters:
  • m_msun – Halo mass \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

  • dm – Accuracy for triangle approximation of derivative \([M_{\odot}]\), scalar

Return dn_dlnm:

Logarithmic mass function \(\frac{dn}{d\ln M}\) \([Mpc^{-3}]\)

dn_dlnm_of_nu(nu, a, dm=100000.0)[source]

Calculates logarithmic mass function \(\frac{dn}{d\ln M}\) as a function of the peak height \(\nu=\delta_c / \sigma(M, a)\). The functional form depends on the selected multiplicity function \(f(\nu, a)\), see the corresponding documentation.

Parameters:
  • nu – Peak height \(\nu=\delta_c / \sigma(M, a)\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

  • dm – Accuracy for triangle approximation of derivative \([M_{\odot}]\), scalar

Return dn_dlnm:

Logarithmic mass function \(\frac{dn}{d\ln M}\) \([Mpc^{-3}]\)

dn_dlogm_of_m(m_msun, a, dm=100000.0)[source]

Calculates logarithmic mass function \(\frac{dn}{d\log_{10} M}\) as a function of the halo mass. The functional form depends on the selected multiplicity function \(f(\nu, a)\), see the corresponding documentation.

Parameters:
  • m_msun – Halo mass \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

  • dm – Accuracy for triangle approximation of derivative \([M_{\odot}]\), scalar

Return dn_dlogm:

Logarithmic mass function \(\frac{dn}{d\log_{10} M}\) \([Mpc^{-3}]\)

dn_dlogm_of_nu(nu, a, dm=100000.0)[source]

Calculates logarithmic mass function \(\frac{dn}{d\log_{10} M}\) as a function of the peak height \(\nu=\delta_c / \sigma(M, a)\). The functional form depends on the selected multiplicity function \(f(\nu, a)\), see the corresponding documentation.

Parameters:
  • nu – Peak height \(\nu=\delta_c / \sigma(M, a)\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

  • dm – Accuracy for triangle approximation of derivative \([M_{\odot}]\), scalar

Return dn_dlogm:

Logarithmic mass function \(\frac{dn}{d\log_{10} M}\) \([Mpc^{-3}]\)

dn_dm_of_m(m_msun, a, dm=100000.0)[source]

Calculates mass function \(\frac{dn}{dM}\) as a function of the halo mass. The functional form depends on the selected multiplicity function \(f(\nu, a)\), see the corresponding documentation.

Parameters:
  • m_msun – Halo mass \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

  • dm – Accuracy for triangle approximation of derivative \([M_{\odot}]\), scalar

Return dn_dm:

Mass function \(\frac{dn}{dM}\) \([M_{\odot}^{-1} Mpc^{-3}]\)

dn_dm_of_nu(nu, a, dm=100000.0)[source]

Calculates mass function \(\frac{dn}{dM}\) as a function of the peak height \(\nu=\delta_c / \sigma(M, a)\). The functional form depends on the selected multiplicity function \(f(\nu, a)\), see the corresponding documentation.

Parameters:
  • nu – Peak height \(\nu=\delta_c / \sigma(M, a)\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

  • dm – Accuracy for triangle approximation of derivative \([M_{\odot}]\), scalar

Return dn_dm:

Mass function \(\frac{dn}{dM}\) \([M_{\odot}^{-1} Mpc^{-3}]\)

f(nu, a=1.0)[source]

Multiplicity function \(f(\nu, a)\) as it appears in the calculation of the Halo Mass Function. The available fitting functions are Press & Schechter (‘PS’), 1974, Sheth & Tormen (‘ST’), 1999 , Tinker et al. (‘Ti’), 2010 , and Watson et al. (‘Wa’), 2013 .

Parameters:
  • nu – Peak height \(\nu=\delta_c / \sigma(M, a)\), scalar or 1d array

  • a – Scale factor, a=1.0 for multiplicity function=’PS’, ‘ST’, and ‘Wa’, scalar or 1d array for ‘Ti’

Returns:

Multiplicity function \(f(\nu, a)\)

Example of setting the multiplicity function:

cosmo.set(multiplicity_fnct=option)

where option can be ‘PS’ for Press & Schechter (1974), ‘ST’ for Sheth & Tormen (1999), ‘Ti’ for Tinker et al., 2010, or ‘Wa’ for Watson et al., 2013.

lin_halo_bias_of_nu(nu, a)[source]

Calculates the linear halo bias for a given peak height \(\nu=\delta_c / \sigma(M, a)\). The available biases are Mo & White (‘MW’), 1996 , Sheth & Tormen (‘ST’), 1999 , Sheth, Mo & Tormen (‘SMT’), 2001 , and Tinker et al. (‘Ti’), 2010 .

Parameters:
  • nu – Peak height \(\nu=\delta_c / \sigma(M, a)\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Return lin_halo_bias:

Linear halo bias

Example of setting the linear halo bias:

cosmo.set(lin_halo_bias_type=option)

where option can be ‘MW’ for Mo & White (1996), ‘ST’ for Sheth & Tormen (1999), ‘SMT’ for Sheth, Mo & Tormen (2001), or ‘Ti’ for Tinker et al., 2010.

mass2nu(m_msun, a)[source]

Calculates the peak height \(\nu(M, a)=\delta_c / \sigma(M, a)\) as a function of the halo mass and scale factor.

Parameters:
  • m_msun – Halo mass \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Return nu:

Peak height

mass2radius(m_msun)[source]

Converts mass of a sphere in solar masses \([M_{\odot}]\) to the corresponding radius in \([Mpc]\), assuming homogeneous density corresponding to the matter density of the universe at redshift \(z=0\).

Parameters:

m_msun – Mass of sphere in solar masses \([M_{\odot}]\), scalar or 1d array

Returns:

Radius of the sphere in \([Mpc]\)

max_redshift(k)[source]

Computes max redshift for which this model is applicable.

Parameters:

k – Wavenumber k, might be needed for some underlying linear perturbation models. \([h Mpc^{-1}]\)

Returns:

Redshift

nu2mass(nu, a)[source]

Extracts the halo mass from Eq.(17) in Mead et al., 2015, \(\nu(M, a)=\delta_c / \sigma(M, a)\), by converting the \(\nu\)-array into a mass-array by backwards interpolation.

Parameters:
  • nu – Peak height \(\nu=\delta_c / \sigma(M, a)\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

Mass in solar masses \([M_{\odot}]\)

pk_1h(k, a)[source]

One-Halo term \(P_\text{1h}(k, a)\) of the non-linear matter power spectrum as described in Eq. (7.191) of Mo et al. “Galaxy Formation and Evolution” (2010).

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

One-Halo power spectrum \(P_\text{1h}(k, a) \ [Mpc^{3}]\)

Example of setting the multiplicity function, linear halo bias, halo profile, and the assumed minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation:

cosmo.set(multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,s
          halo_profile=option_p,
          min_halo_mass=1e8, max_halo_mass=1e14)

where the possible option_mf are given in the documentation to \(f(\nu, a)\), option_b in the documentation to the linear halo bias and option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

pk_2h(k, a)[source]

Two-Halo term \(P_\text{2h}(k, a)\) of the non-linear matter power spectrum as described in Eq. (7.194) of Mo et al. “Galaxy Formation and Evolution” (2010).

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

Two-Halo power spectrum \(P_\text{2h}(k, a) \ [Mpc^{3}]\)

Example of setting the multiplicity function, linear halo bias, halo profile, and the assumed minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation:

cosmo.set(multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,
          halo_profile=option_p,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

where the possible option_mf are given in the documentation to \(f(\nu, a)\), option_b in the documentation to the linear halo bias and option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]

Calculates the non-linear matter power spectrum using the Halo Model, described e.g. in Mo et al. “Galaxy Formation and Evolution” (2010) chapter 7.6, as the sum of the pk_1h() and pk_2h() terms.

Parameters:
  • k – Wavelength \([Mpc^{-1}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

Halo Model power spectrum, \(P_\text{nl}(k, a)\), in \([Mpc^{3}]\)

Example on how to use the Halo model, set the multiplicity function, linear halo bias, halo profile, and the minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation, and then calculate the power spectrum:

cosmo.set(pk_nonlin_type='HaloModel',
          multiplicity_fnct=option_mf,
          lin_halo_bias_type=option_b,
          halo_profile=option_p,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

pk = cosmo.nonlin_pert.powerspec_a_k(a,k)

where the possible option_mf are given in the documentation to \(f(\nu, a)\), option_b in the documentation to the linear halo bias and option_p is either True for a NFW-profile or False for assumed point sources. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

print_params()[source]

Prints the cosmological setup and the parameters used for the computation of the non-linear matter power spectrum with the halo model.

Example:

cosmo.set(pk_nonlin_type='HaloModel')
cosmo.nonlin_pert.print_params()
radius2mass(r_mpc)[source]

Converts radius of a sphere in \([Mpc]\) to corresponding mass in solar masses \([M_{\odot}]\), assuming homogeneous density corresponding to the matter density of the universe at redshift \(z=0\).

Parameters:

r_mpc – Radius of sphere \([Mpc]\), scalar or 1d array

Returns:

Mass of sphere in solar masses \([M_{\odot}]\)

rho_m(a)[source]

Calculates the mean dark matter density as a function of the scale factor.

Parameters:

a – Scale factor, scalar or 1d array

Returns:

Mean matter density, in \([M_{\odot} \ Mpc^{-3}]\)

Example of setting the multiplicity function and the assumed minimal and maximal halo masses in solar masses \([M_{\odot}]\) considered in the calculation:

cosmo.set(multiplicity_fnct=option,
          min_halo_mass=1e8,
          max_halo_mass=1e14)

where option can be ‘PS’ for Press & Schechter (1974), ‘ST’ for Sheth & Tormen (1999), ‘Ti’ for Tinker et al., 2010, or ‘Wa’ for Watson et al., 2013. Default values for the minimal and maximal halo masses are min_halo_mass=1 and max_halo_mass=1e+20.

rvir(mvir_msun, a)[source]

Calculates the virial radius corresponding to a dark matter halo of a given mass for different cosmologies. See Bryan & Norman, 1998 for more details.

Parameters:
  • mvir_msun – Halo mass in solar masses \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, 1d array

Returns:

Virial radius in \([Mpc]\)

sigma(m_msun, a)[source]

Calculates \(\sigma(M, a)\), the RMS of the density field at a given mass.

Parameters:
  • m_msun – Mass in solar masses at which \(\sigma\) is evaluated \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, scalar or 1d array

Returns:

\(\sigma(M, a)\) as the RMS of the density field

sigma8_a(a)[source]

Computes \(\sigma_8\), the RMS density contrast fluctuation, smoothed with a top hat of radius 8 \(h^{-1}Mpc\) at scale factor \(a\).

Parameters:

a – Scale factor, scalar or 1d array

Returns:

\(\sigma_8\) at the desired scale factor

subhalos_dn_dm_of_m(subhalo_m_msun, mean_subhalo_mass_fraction=0.5, host_halo_mass_msun=1000000000000.0, beta=0.3, gamma=0.9)[source]

Compute subhalo mass function as a function of subhalo masses and host halo mass.

References:

Chapter 7.5 of Mo, van den Bosch and White (Equations 7.149 and 7.150) Left panel of Figure 2 from Jiang et al., 2017 informs the value for mean_subhalo_mass_fraction ~ 0.1

Recommended values:
mean_subhalo_mass_fraction: By definition, between 0 and 1.

Default=0.1, informed by arxiv:1610.02399

beta ~ 0.1 - 0.5, always < 1 (according to Mo, van den Bosch and White)

gamma ~ 0.8 - 1.0 (according to Mo, van den Bosch and White)

Parameters:
  • subhalo_m_msun – Subhalo mass \([M_{\odot}]\), scalar or 1d array

  • mean_subhalo_mass_fraction – Mean mass fraction contained in subhalos, scalar

  • host_halo_mass_msun – Mass of the host halo for which this subhalo mass function is computed \([M_{\odot}]\), scalar

  • beta – Exponent, scalar

Return dn_dm:

Subhalo mass function \(\frac{dn}{dM}\) \([M_{\odot}^{-1} Mpc^{-3}]\)

vvir(mvir_msun, a, diag_only=False)[source]

Calculates the virial velocity corresponding to a dark matter halo of a given mass as described in Barnes & Haehnelt, 2014 .

Parameters:
  • mvir_msun – Halo mass in solar masses \([M_{\odot}]\), scalar or 1d array

  • a – Scale factor, 1d array

Returns:

Virial velocity in \([km/s]\)

PyCosmo.NonLinearPerturbation_HaloModel.fft_top_hat(x)[source]

Fourier transform of 3D spherical top hat.

Parameters:

x – Dimensionless parameter (usually k r)

Return t:

Fourier transform of a 3D spherical top hat.

PyCosmo.Obs module

class PyCosmo.Obs.Obs[source]

Bases: object

The class is used to create an instance that is needed to calculate the observables, containing information about experiments and surveys. In this way the PyCosmo.Obs class acts at the same level as the PyCosmo.Cosmo class and is initialised using a parameter file or a set of parameters. As shown also in the Section Usage of the documentation, the following example illustrates how to instantiate the PyCosmo.Obs class to compute the lensing power spectrum \(C_{\ell}^{\gamma \gamma}\):

# PyCosmo instance
cosmo = PyCosmo.Cosmo()

# Set the fitting functions for the power spectrum
cosmo.set(pk_type = 'EH')
cosmo.set(pk_nonlin_type = 'rev_halofit')

# Information about the observables
# Note: 'perturb': 'nonlinear' can be replaced with 'linear' by the user.
clparams = {'nz': ['smail', 'smail'],
            'z0': 1.13,
            'beta': 2.,
            'alpha': 2.,
            'probes': ['gamma', 'gamma'],
            'perturb': 'nonlinear',
            'normalised': False,
            'bias': 1.,
            'm': 0,
            'z_grid': np.array([0.0001, 5., 1000])
            }

# Obs instance
obs = PyCosmo.Obs()

# Cls computation
cls_haloEH = obs.cl(ells, cosmo, clparams)

where the function cl() is documented below.

Warning

Obs not supported by the Boltzmann solver yet, since it uses the Projection class. Calling any function from this module while using the Boltzmann solver (pk_type = ‘boltz’) will return ValueError: you must set pk_nonlin_type to access projections

C1(cosmo)[source]

Returns th:e linear intrinsic alignment amplitude C1.

Parameters:

cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

Returns:

Linear intrinsic alignment amplitude for current value of h

F(a, cosmo, obsparams)[source]

Returns the \(F\) function of the linear alignment model. This follows Eq. (8) of Hildebrandt et al., 2016. Here we set \(\eta = \beta = 0\), i.e., no redshift and luminosity dependence is considered.

Parameters:
  • a – array of scale factor values. \(a\)

  • cosmo – instance of PyCosmo.Cosmo (specifies the cosmological model)

  • obsparams – dictionary containing specifications for the observables

Returns:

Value of \(F\) as a function of the scale factor \(a\).

bin_setup(obsparams)[source]

Set up the redshift selection functions, i.e. each redshift selection function is either calculated from the default hard coded ones or the tabulated selection functions are read in and saved as attributes in self.n.

Parameters:

obsparams – dictionary containing specifications for observable

Returns:

bin_setup_old(obsparams)[source]

Set up the redshift selection functions, i.e. each redshift selection function is either calculated from the default hard coded ones or the tabulated selection functions are read in and saved as attributes in self.n.

Parameters:

obsparams – dictionary containing specifications for observable

Returns:

cl(ells, cosmo, obsparams)[source]

Wrapper around all the other cl routines which determines from the obsparams dictionary for which cosmological probes the angular power spectrum has to be computed. The provided cosmological probes are galaxy/halo overdensity ‘deltag’, neutral hydrogen HI overdensity ‘HI’, cosmic shear ‘gamma’, and CMB lensing convergence ‘cmbkappa’.

This is especially useful when interfacing with CosmoHammer.

Parameters:
  • ells – array of angular multipoles \(\ell\)

  • cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

  • obsparams – dictionary containing specifications for observable

Returns:

Array of values of the angular power spectrum for the desired probes for \(\ell\)’s

cl_IG(ells, cosmo, obsparams)[source]

Calculates the angular power spectrum for GI intrinsic alignments (IAs) using the NLA or LA model for the multipole array specified by the angular multipoles \(\ell\). The NLA implementation follows Hildebrandt et al., 2016.

Parameters:
  • ells – array of angular multipoles \(\ell\)

  • cosmo – instance of PyCosmo.Cosmo (it specifies the cosmological model)

  • obsparams – dictionary containing specifications for the observables

Returns:

Array of values of the angular power spectrum for GI IAs.

cl_II(ells, cosmo, obsparams)[source]

Calculates the angular power spectrum for II intrinsic alignments (IAs) using the NLA or LA model for the multipole array specified by \(\ell\). The NLA implementation follows Hildebrandt et al., 2016.

Parameters:
  • ells – array of angular multipoles \(\ell\)

  • cosmo – instance of PyCosmo.Cosmo (it specifies the cosmological model)

  • obsparams – dictionary containing specifications for the observables

Returns:

Array of values of the angular power spectrum for II IAs.

cl_multi(ells, cosmo, obsparams_list)[source]

Calculates the \(C_{\ell}\) for the desired probes, as specified by the “probes” parameter in obsparams.

Parameters:
  • ells – array of angular multipoles \(\ell\)

  • cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

  • obsparams_list – list of dictionaries containing specifications for observables

growth_ISW(a, projection, mode='num')[source]

Returns the generalised growth factor used in the ISW angular power spectrum integrand both for the analytic approximation and the numerical derivative.

Parameters:
  • a – scale factor

  • cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

  • mode – numerical or analytic

Returns:

The value of the generalised growth factor at redshift x

growth_suba(cosmo, a=1.0)[source]

Returns the growth factor divided by the scale factor a.

Parameters:
  • cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

  • a – scale factor a

Returns:

Growth factor divided by scale factor

growth_suba_deriv(cosmo, a=1.0)[source]

Returns the derivative of growth_suba by scale factor a.

Parameters:
  • cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

  • a – scale factor a

Returns:

Derivative of growth_suba w.r.t. z

nz(obsparams, nzmode=None, path2zdist=None)[source]

Computes and normalises the redshift selection-function of a survey.

Parameters:
  • obsparams – dictionary containing specifications for the observables

  • nzmode – string tag in [smail, cfhtlens, custom, custom_array, uniform]

  • path2zdist – if \(n(z)\) comes from a file, then this is the path to the file (nzmode = custom) if \(n(z)\) comes from an array, this is the array (nzmode = custom_array)

  • z – redshift grid

Returns:

Corresponding normalised redshift selection-function, \(P(z)_{|norm}\).

print_params()[source]
rhocrit(cosmo)[source]

Return the critical density of the Universe at a_0 = 1 (today).

Parameters:

cosmo – instance of PyCosmo.Cosmo (specifies cosmological model)

Returns:

critical density today, [rhoc] = kg/Mpc^3

setup(obsparams)[source]

Sets the multiplicative calibration parameters as defined in the dictionary obsparams.

Parameters:

obsparams – dictionary containing specifications for the observables

windows(obsparams)[source]

Sets up the radial window functions for the surveys and the probes, i.e., computes the window functions appropriate for the galaxy overdensity, HI, \(\gamma\), CMB temperature and CMB \(- \kappa\) and saves them in the attribute self.window.

Parameters:

obsparams – dictionary containing specifications for observable

Return window:

radial window function

PyCosmo.PerturbationBase module

class PyCosmo.PerturbationBase.ClassContractMeta(name, bases, d, **kwargs)[source]

Bases: object

class PyCosmo.PerturbationBase.NonLinearPerturbationBase(*a, **kw)[source]

Bases: object

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]
PyCosmo.PerturbationBase.check_protypes(clz)[source]
PyCosmo.PerturbationBase.prototype(func)[source]

PyCosmo.PerturbationTable module

class PyCosmo.PerturbationTable.PerturbationTable(cosmo, perturbation)[source]

Bases: object

This class has been written to handle operations involving building tabulated data. For some of the calculations this is often useful. These tables can then be accessed through interpolation routines rather than doing the full calculations repeatedly.

powerspec_a_k(a=1.0, k=0.1, diag_only=False)[source]
PyCosmo.PerturbationTable.check_diff(diff, thresh)[source]
PyCosmo.PerturbationTable.optimize_grid(perturbation, initial_k_grid, tolerance, width, pk_tobe=None)[source]

PyCosmo.Projection module

class PyCosmo.Projection.Projection(params, background, lin_pert, nonlin_pert)[source]

Bases: object

Core projection functions to be used by the Obs class to predict observables.

Warning

Projection not supported by the Boltzmann solver yet, since it requires nonlinear perturbations. Calling any function from this module while using the Boltzmann solver (pk_type = ‘boltz’) will return ValueError: you must set pk_nonlin_type to access projections

cl_limber(ells, weight_function1, weight_function2, a_grid, probes, perturb='nonlinear')[source]

Computes the angular power spectrum of the auto- or cross-correlation between two LSS probes i.e. galaxy/halo overdensity, HI overdensity, cosmic shear or CMB lensing in the Limber Approximation.

Parameters:
  • ells – array of angular multipoles \(\ell\)

  • weight_function1 – radial weight function for first probe

  • weight_function2 – radial weight function for second probe

  • a_grid – integration grid in scale factor a

  • probes – list with name of probes

  • perturb – string tag denoting if we use linear or nonlinear power spectrum. It can be either linear or nonlinear.

Returns:

Angular power spectrum \(C_{\ell}\) at the multipoles \(\ell\).

cl_limber_IG(ell, weight_function1, weight_function2, F, a_grid, IAmodel='NLA')[source]

Computes the angular power spectrum of the cross correlation between intrinsic galaxy ellipticities and tracers of the LSS.

Parameters:
  • ell – array of angular multipole \(\ell\)

  • weight_function1 – radial weight function for first probe -> this needs to be the weight function

  • weight_function2 – radial weight function for second probe -> this needs to be n(z)

  • F – IA bias function

  • a_grid – integration grid in scale factor a

  • IAmodel – string tag denoting if we use NLA or LA IA model. It can be either NLA or IA.

Returns:

Angular power spectrum \(C_{\ell}\) at the multipoles \(\ell\).

cl_limber_IG_int(a, ells, weight_function1, weight_function2, F, IAmodel='NLA')[source]

Returns the integrand for the angular power spectrum of intrinsic alignments (IAs) computed using the NLA or LA model.

Parameters:
  • a – array of scale factor values a

  • ells – array of angular multipoles \(\ell\)

  • weight_function1 – radial weight function for first probe -> this needs to be the weight function

  • weight_function2 – radial weight function for second probe -> this needs to be n(z)

  • F – IA bias function

  • IAmodel – string tag denoting if we use NLA or LA IA model. It can be either NLA or IA.

Returns:

Integrand at \(\\ell\) for the values of a.

cl_limber_II(ell, weight_function1, weight_function2, F, a_grid, IAmodel='NLA')[source]

Computes the angular power spectrum of the auto power spectrum of intrinsic galaxy ellipticities.

Parameters:
  • ell – array of angular multipole \(\ell\)

  • weight_function1 – redshift selection function for first probe

  • weight_function2 – redshift selection function for second probe

  • F – IA bias function

  • a_grid – integration grid in scale factor a

  • IAmodel – string tag denoting if we use NLA or LA IA model. It can be either NLA or IA.

Returns:

Angular power spectrum \(C_{\ell}\) at the multipoles \(\ell\).

cl_limber_II_int(a, ells, weight_function1, weight_function2, F, IAmodel='NLA')[source]

Returns the integrand for the angular power spectrum of the auto correlation of intrinsic alignments (IAs) computed using the NLA or LA model.

Parameters:
  • a – array of scale factor values a

  • ells – array of angular multipoles \(\ell\)

  • weight_function1 – redshift selection function for first probe

  • weight_function2 – redshift selection function for second probe

  • F – IA bias function

  • IAmodel – string tag denoting if we use NLA or LA IA model. It can be either NLA or IA.

Returns:

Integrand at \(\ell\) for the values of a.

cl_limber_ISW(ell, weight_function1, weight_function2, growth_ISW, a_grid, perturb='linear')[source]

Computes the angular power spectrum of the cross correlation between the CMB temperature anisotropies and the galaxy overdensity/cosmic shear.

Parameters:
  • ell – array of angular multipole \(\ell\)

  • weight_function1 – radial weight function for first probe

  • weight_function2 – radial weight function for second probe

  • growth_ISW – growth function for ISW

  • a_grid – integration grid in scale factor a

  • linear – string tag denoting if we use linear or nonlinear power spectrum. It can be either linear or nonlinear

Returns:

Angular power spectrum \(.C_{\ell}\) at the multipoles \(\ell\).

cl_limber_ISW_int(a, ells, weight_function1, weight_function2, growth_ISW, perturb='linear')[source]

Returns the integrand needed to compute the angular power spectrum of the cross correlation between the CMB temperature anisotropies and the galaxy overdensity/cosmic shear/CMB lensing.

Parameters:
  • a – array of scale factor values a

  • ells – array of angular multipoles \(\ell\)

  • weight_function1 – radial weight function for first probe

  • weight_function2 – radial weight function for second probe

  • growth_ISW – growth function for ISW

  • perturb – string tag denoting if we use linear or nonlinear power spectrum. It can be either linear or nonlinear.

Returns:

Integrand at \(\ell\) for the values of a.

cl_limber_int(a, ells, weight_function1, weight_function2, probes, perturb='nonlinear')[source]

Returns the integrand needed to compute the angular power spectrum of the auto- or cross-correlation between two LSS probes in the Limber approximation.

Parameters:
  • a – array of scale factor values a

  • ells – array of angular multipoles \(\ell\)

  • weight_function1 – radial weight function for first probe

  • weight_function2 – radial weight function for second probe

  • probes – list with name of probes

  • perturb – string tag denoting if we use linear or nonlinear power spectrum. It can be either linear or nonlinear.

Returns:

Integrand at \(\ell\) for the values of a.

PyCosmo.TheoryPredTables module

class PyCosmo.TheoryPredTables.Spline(x, y)[source]

Bases: object

class PyCosmo.TheoryPredTables.TheoryPredTables(lin_pert)[source]

Bases: object

Class for calculation of tabulated cosmological quantities from PyCosmo.

growth_deriv(a=1.0, mode='fivepoint')[source]

Returns the derivative of the growth function w.r.t. a. :param a: scale factor :param mode: 5-point derivative or 2-point derivative :return growth_deriv: derivative of the growth factor w.r.t. a

growth_suba(a=1.0)[source]

Returns the growth factor divided by the scale factor a :param a: scale factor a :return growth_suba_temp: growth factor divided by scale factor

growth_suba_deriv(a=1.0)[source]

Returns the derivative of growth_suba by redshift z :param a: scale factor a :return growth_deriv: derivative of growth_suba w.r.t. z

growth_tab_a(a=1.0)[source]

Returns the interpolated growth factor from PyCosmo. :param a: array of scale factor values a :return growth_interp: interpolated growth factor at a

inv_growth_tab_a(g=1.0)[source]

Returns the interpolated inverse of the growth factor from PyCosmo. :param g: array of growth factor values g :return inv_growth_interp: interpolated inverse growth factor at g

powerspec_tab_a_k(a: object = 1.0, k: object = 0.1, fields: object = None) object[source]

Returns the interpolated linear matter power spectrum from PyCosmo. :param a: array of scale factor values a :param k: array of wave vector k values :return pk_lin_interp: interpolated linear power spectrum at a and k

PyCosmo.TheoryPredTables.Timer(*a, **kw)[source]

PyCosmo.build module

PyCosmo.build.build(model_name='lcdm', **extra_settings)[source]
PyCosmo.build.update_nested(original, updates)[source]

PyCosmo.core_file_handling module

PyCosmo.core_file_handling.extract_code_from_nb(path, cosmology_config)[source]
PyCosmo.core_file_handling.extract_code_from_py(paths, cosmology_config)[source]
PyCosmo.core_file_handling.import_equations_file(paths, cosmology)[source]
PyCosmo.core_file_handling.load_core_files(paths, cosmology_config)[source]
PyCosmo.core_file_handling.remove_globals(f)[source]
PyCosmo.core_file_handling.setup_wrappers(module)[source]

PyCosmo.disable_multithreading module

PyCosmo.disable_multithreading.disable_multithreading()[source]

try to disable multithreading in numpy et al.

source: https://stackoverflow.com/questions/30791550/

PyCosmo.ini_handling module

class PyCosmo.ini_handling.Bunch(dict_=None)[source]

Bases: object

as_dict()[source]
PyCosmo.ini_handling.best_convert(txt)[source]
PyCosmo.ini_handling.load_ini(path)[source]

PyCosmo.model_config module

PyCosmo.model_config.model_config(model_name, extra_settings)[source]

PyCosmo.optimize_nq module

PyCosmo.optimize_nq.gl_quadrature(mnu_relerr)[source]

PyCosmo.patches module

PyCosmo.predefined_k_grids module

PyCosmo.predefined_k_grids.k_precise = array([1.00000000e-04, 1.02804473e-04, 1.08651577e-04, 1.11698682e-04,        1.14831241e-04, 1.16430313e-04, 1.18051653e-04, 1.21362380e-04,        1.24765955e-04, 1.28264983e-04, 1.31862140e-04, 1.35560179e-04,        1.39361927e-04, 1.41302599e-04, 1.43270295e-04, 1.47288272e-04,        1.51418933e-04, 1.60031031e-04, 1.69132952e-04, 1.78752553e-04,        1.88919278e-04, 1.94217468e-04, 2.05263775e-04, 2.29276931e-04,        2.42317279e-04, 2.70665207e-04, 3.37698031e-04, 5.25679112e-04,        8.18300682e-04, 1.27381132e-03, 1.58928287e-03, 1.77520801e-03,        1.98288395e-03, 2.21485523e-03, 2.34082728e-03, 2.61467321e-03,        4.07014245e-03, 5.07815211e-03, 6.33580499e-03, 7.90492762e-03,        9.86265846e-03, 1.23052400e-02, 1.53527503e-02, 1.71488197e-02,        1.91550056e-02, 2.13958887e-02, 2.38989257e-02, 2.66947849e-02,        2.98177229e-02, 3.33060034e-02, 3.72023668e-02, 4.15545533e-02,        4.64158883e-02, 4.90558371e-02, 5.18459354e-02, 5.47947234e-02,        5.79112265e-02, 6.12049837e-02, 6.46860766e-02, 7.22534949e-02,        7.63629826e-02, 8.07062014e-02, 8.52964450e-02, 9.01477631e-02,        9.52750047e-02, 1.00693863e-01, 1.06420924e-01, 1.12473718e-01,        1.15628013e-01, 1.18870770e-01, 1.22204469e-01, 1.29154967e-01,        1.32777083e-01, 1.36500781e-01, 1.40328908e-01, 1.44264395e-01,        1.48310251e-01, 1.56745541e-01, 1.61141428e-01, 1.65660596e-01,        1.70306503e-01, 1.79992851e-01, 1.85040702e-01, 1.90230119e-01,        1.95565072e-01, 2.01049642e-01, 2.06688025e-01, 2.12484535e-01,        2.18443607e-01, 2.24569800e-01, 2.30867799e-01, 2.37342425e-01,        2.50841506e-01, 2.57876289e-01, 2.72543253e-01, 2.80186656e-01,        2.88044415e-01, 2.96122544e-01, 3.04427221e-01, 3.21741815e-01,        3.30764978e-01, 3.40041193e-01, 3.49577557e-01, 3.69460121e-01,        3.79821531e-01, 3.90473524e-01, 4.01424249e-01, 4.24255643e-01,        4.36153779e-01, 4.48385595e-01, 4.60960449e-01, 5.00840799e-01,        5.14886745e-01, 5.44171429e-01, 5.75121707e-01, 6.07832313e-01,        7.58367791e-01, 1.18051653e+00, 1.47288272e+00, 1.83765620e+00,        2.56099310e+00, 3.56904935e+00, 9.66017480e+00, 2.34082728e+01,        3.26222201e+01, 4.07014245e+01, 5.07815211e+01, 5.67222897e+01,        6.69616005e+01, 7.90492762e+01, 9.33189772e+01, 1.00000000e+02])
>>> len(k2)
106
>>> e1
7.999793047824832e-05
>>> e2
0.0004712228618850567
>>>

Module contents

This is the PyCosmo package.

class PyCosmo.Cosmo(*, model_config=None)[source]

Bases: object

All of main functionalities of PyCosmo are managed by the main (Cosmo) class, which links to the other classes where the bulk of the calculations is performed.

load_params(paramfile=None)[source]

Loads parameter file and additional user-defined parameters to initialize the Cosmo object

Parameters:

paramfile – path to parameter file

print_params(inc_derived=True, inc_internal=True, file=None)[source]

Prints the parameters of PyCosmo instance.

Parameters:
  • inc_derived – prints derived parameters [True or False]

  • inc_internal – prints internal parameterss (e.g. for lin_pert) [True or False]

static recompile_from_cli()[source]
reset_wrapper_globals()[source]
set(**kwargs)[source]

Changes the configuration parameters of the PyCosmo instance. This will also recalculate the derived quantities.

Parameters:

kwargs – keyword arguments for the input parameters to change.

Warning

Always change parameters of an instance of PyCosmo using the set function. Do not change variables directly!

Returns:

If no keywords are passed then a list of options is printed

Example:

cosmo.set(omega_m=0.23)