PyCosmo package

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]

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

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]

Calculates conformal time \(\eta\) to a given scale factor.

Parameters

a – scale factor [1]

Returns

\(\eta\) conformal time [\(\mathrm{Mpc}/h\)]

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)

LinearPerturbationApprox module

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

Bases: PyCosmo.LinearPerturbationBase.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}\).

BoltzmannSolver.LinearPerturbationBoltzmann module

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

Bases: PyCosmo.LinearPerturbationBase.LinearPerturbationBase

Class for computing linear perturbations by solving the Einstein-Boltzmann ODE system.

The Boltzmann solver is selected using the set function:

cosmo.set(pk_type = "boltz")
fields(k, grid, sec_factor=5, keep_lna0=False, initial_conditions=None, enable_fast_solver=True, enable_sparse_lu_solver=True)[source]

Solves the Einstein-Boltzmann ODE system for the evolution of the linear order perturbation of the fields.

Parameters
  • k – wavenumber \([h/Mpc]\)

  • grid – ln(a) values at which to output fields [1]

  • sec_factor – relaxes row permutation criterium in optimized LUP solver

  • keep_lna0 – if True includes the fields at initial time a_0

  • initial_conditions – can pass a_0, y_0 (vector of initial conditions)

  • enable_fast_solver – if set to False always use standard LUP solver for full matrix

  • enable_sparse_lu_solver – if set to True: avoid iterating over known zero-entries in fallback LUP solver

Returns

Linear order perturbations, accessed with fields.a, fields.Phi etc.

growth_a(a=1.0, k=None)[source]

Returns the linear growth factor computed from the Boltzmann solver at a given k value.

Parameters
  • a – scale factor [1]

  • k – wavenumber used to compute the growth factor (default:k=1 \(Mpc^{-1}\)) [\(Mpc^{-1}\)]

Returns

D(a): growth factor normalised to 1 at a=1 [1]

max_redshift(k)[source]

computes max redshift for which this model is applicable. uses the implemented initial conditions to determine this.

Parameters

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

Returns

redshift

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

Returns the linear total matter power spectrum computed from the Boltzmann solver, this includes cold dark matter, baryons and massive neutrinos for an array of a and k values.

Parameters
  • a – scale factor [1] (default:a=1)

  • k – wavenumber used to compute the power spectrum (default:k=1 Mpc^1) [Mpc^-1]

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

  • pool – instance of multiprocessing.pool.Pool for parallel computation

Warning

If pk_norm_type = ‘A_s’ this will compute the total matter power spectrum following the general relativistic treatment (using delta_m_tot). If pk_norm_type = ‘deltah’ the Poisson equation will be assumed and the evolution will be split in growth factor and transfer function. We recommend the use of A_s or sigma8 normalizations.

Returns

P(a,k): total matter power spectrum [\(Mpc^3\)]

powerspec_cb_a_k(a=1.0, k=0.1, diag_only=False, pool=None)[source]

Returns the linear matter power spectrum of cold dark matter and baryons computed from the Boltzmann solver for an array of a and k values.

Parameters
  • a – scale factor [1] (default:a=1)

  • k – wavenumber used to compute the power spectrum (default:k=1 \(Mpc^{-1}\)) [\(Mpc^{-1}\)]

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

  • pool – instance of multiprocessing.pool.Pool for parallel computation

Returns

P_cb(a,k): power spectrum of CDM+baryons [\(Mpc^3\)]

transfer_k(k, pool=None)[source]

Computes the linear matter transfer function using the Boltzmann solver assuming the Poisson equation.

Parameters

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

Returns

Transfer function \([ Mpc^{3/2} ]\)

PyCosmo.BoltzmannSolver.LinearPerturbationBoltzmann.print_new_traces(cache_file)[source]

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 max_redshift(k)[source]
sigma8(a=1.0)[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]

NonLinearPerturbationHaloFit module

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

Bases: PyCosmo.PerturbationBase.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

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

NonLinearPerturbationMead module

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

Bases: PyCosmo.PerturbationBase.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 , and Tinker et al., 2010 .

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) or ‘Ti’ for Tinker et al., 2010.

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

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

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, cosmo, perturb='nonlinear')[source]

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

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

  • 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, cosmo, 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

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

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

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.

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 ells

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]

  • path2zdist – if \(n(z)\) comes from a file, then this is the path to the file

  • 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, \(\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