nbodykit.cosmology.cosmology

Functions

astropy_to_dict(cosmo) Convert an astropy cosmology object to a dictionary of parameters suitable for initializing a Cosmology object.
check_args(args)
check_deprecated_init(args, kwargs) Check if kwargs uses the (now deprecated) signature of Cosmology prior to version 0.2.6.
compile_args(args) Compile the input args of Cosmology object to the input parameters (pars) to a Cosmology object.
find_eqcls(key)
merge_args(args, moreargs) merge moreargs into args.
store_user_kwargs() Decorator that adds the _user_kwargs attribute to the class to track which arguments the user actually supplied.

Classes

Cosmology([h, T0_cmb, Omega0_b, Omega0_cdm, …]) A cosmology calculator based on the CLASS binding in classylss.
nbodykit.cosmology.cosmology.astropy_to_dict(cosmo)[source]

Convert an astropy cosmology object to a dictionary of parameters suitable for initializing a Cosmology object.

nbodykit.cosmology.cosmology.check_deprecated_init(args, kwargs)[source]

Check if kwargs uses the (now deprecated) signature of Cosmology prior to version 0.2.6.

If using the deprecated syntax, this returns the necessary arguments for the new signature, and None otherwise.

nbodykit.cosmology.cosmology.compile_args(args)[source]

Compile the input args of Cosmology object to the input parameters (pars) to a Cosmology object.

A variety of defaults are set to tune CLASS for quantities used in large scale structures.

Difference between pars and args:
  • anything that is valid pars is also valid args.
  • after replacing our customizations in args, we get pars.

Note that CLASS will check for additional conflicts.

see merge_args()

nbodykit.cosmology.cosmology.merge_args(args, moreargs)[source]

merge moreargs into args.

Those defined in moreargs takes priority than those defined in args.

see compile_args()

nbodykit.cosmology.cosmology.store_user_kwargs()[source]

Decorator that adds the _user_kwargs attribute to the class to track which arguments the user actually supplied.

class nbodykit.cosmology.cosmology.Cosmology(h=0.67556, T0_cmb=2.7255, Omega0_b=0.0482754208891869, Omega0_cdm=0.26377065934278865, N_ur=None, m_ncdm=[0.06], P_k_max=10.0, P_z_max=100.0, gauge='synchronous', n_s=0.9667, nonlinear=False, verbose=False, **kwargs)[source]

A cosmology calculator based on the CLASS binding in classylss.

It is a collection of all method provided by the CLASS interfaces. The object is immutable. To obtain an instance with a new set of parameters use clone() or match().

The individual interfaces can be accessed too, such that c.Spectra.get_transfer and c.get_transfer are identical.

Important

A default set of units is assumed. Those units are:

  • temperature: \(\mathrm{K}\)
  • distance: \(h^{-1} \mathrm{Mpc}\)
  • wavenumber: \(h \mathrm{Mpc}^{-1}\)
  • power: \(h^{-3} \mathrm{Mpc}^3\)
  • density: \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\)
  • neutrino mass: \(\mathrm{eV}\)
  • time: \(\mathrm{Gyr}\)
  • \(H_0\): \((\mathrm{km} \ \mathrm{s^{-1}}) / (h^{-1} \ \mathrm{Mpc})\)

Notes

  • The default configuration assumes a flat cosmology, \(\Omega_{0,k}=0\). Pass Omega0_k as a keyword to specify the desired non-flat curvature.
  • For consistency of variable names, the present day values can be passed with or without ‘0’ postfix, e.g., Omega0_cdm is translated to Omega_cdm as CLASS always uses the names without 0 as input parameters.
  • By default, a cosmological constant (Omega0_lambda) is assumed, with its density value inferred by the curvature condition.
  • Non-cosmological constant dark energy can be used by specifying the w0_fld, wa_fld, and/or Omega_fld values.
  • To pass in CLASS parameters that are not valid Python argument names, use the dictionary/keyword arguments trick, e.g. Cosmology(..., **{'temperature contributions': 'y'})
  • Cosmology(**dict(c)) is not supposed to work; use Cosmology.from_dict(dict(c)).
Parameters:
  • h (float) – the dimensionless Hubble parameter
  • T0_cmb (float) – the temperature of the CMB in Kelvins
  • Omega0_b (float) – the current baryon density parameter, \(\Omega_{b,0}\). Currently unrealistic cosmology where Omega_b == 0 is not supported.
  • Omega0_cdm (float) – the current cold dark matter density parameter, \(\Omega_{cdm,0}\)
  • N_ur (float) – the number of ultra-relativistic (massless neutrino) species; the default number is inferred based on the number of massive neutrinos via the following logic: if you have respectively 1,2,3 massive neutrinos and use the default T_ncdm value (0.71611 K), designed to give m/omega of 93.14 eV, and you wish to have N_eff=3.046 in the early universe, then N_ur is set to 2.0328, 1.0196, 0.00641, respectively.
  • m_ncdm (list, None) – the masses (in eV) for all massive neutrino species; an empty list should be passed for no massive neutrinos. The default is a single massive neutrino with mass of 0.06 eV
  • P_k_max (float) – the maximum k value to compute power spectrum results to, in units of \(h/Mpc\)
  • P_z_max (float) – the maximum redshift to compute power spectrum results to
  • gauge (str,) – either synchronous or newtonian
  • n_s (float) – the tilt of the primordial power spectrum
  • nonlinear (bool) – whether to compute nonlinear power spectrum results via HaloFit
  • verbose (bool) – whether to turn on the default CLASS logging for all submodules
  • **kwargs – extra keyword parameters to pass to CLASS. Mainly used to pass-in parameter names that are not valid Python function argument names, e.g. temperature contributions, or number count contributions. Users should be wary of configuration options that may conflict with the base set of parameters. To override parameters, chain the result with clone().

Attributes

Omega0_cb The total density of CDM and Baryon.
sigma8 The amplitude of matter fluctuations at \(z=0\) in a sphere of radius \(r = 8 \ h^{-1}\mathrm{Mpc}\).

Methods

clone(**kwargs) Create a new cosmology based on modification of self, with the input keyword parameters changed.
from_astropy(cosmo, **kwargs) Initialize and return a Cosmology object from a subclass of astropy.cosmology.FLRW.
from_dict(pars) Creates a Cosmology from a pars dictionary.
from_file(filename, **kwargs) Initialize a Cosmology object from the CLASS parameter file
match([sigma8, Omega0_cb, Omega0_m]) Creates a new cosmology that matches a derived parameter.
to_astropy() Initialize and return a subclass of astropy.cosmology.FLRW from the Cosmology class.

Attributes

A_s The scalar amplitude of the primordial power spectrum at \(k_\mathrm{pivot}\).
C The speed of light in units of km/s.
G The gravitational constant in units of \(10^{-10} \ (M_\odot/h)^{-1} (\mathrm{Mpc}/h) \mathrm{km}^2 \mathrm{s}^{-2}\).
H0 Current Hubble parameter in units of \(\mathrm{km/s} (\mathrm{Mpc}/h)^{-1}.\)
N_ncdm The number of distinguishable ncdm (massive neutrino) species.
N_ur The number of ultra-relativistic species.
Neff Effective number of relativistic species, summed over ultra-relativistic and ncdm species.
Ob0 Returns Omega0_b.
Ode0 Returns the sum of Omega0_lambda and Omega0_fld.
Odm0 Returns Omega0_cdm.
Ogamma0 Returns Omega0_g.
Ok0 Returns Omega0_k.
Om0 Returns Omega0_m.
Omega0_b Current density parameter for photons, \(\Omega_{b,0}\).
Omega0_cb The total density of CDM and Baryon.
Omega0_cdm Current density parameter for cold dark matter, \(\Omega_{cdm,0}\).
Omega0_dcdm Current density parammeter for decaying cold dark matter, \(\Omega_{dcdm,0}\).
Omega0_fld Current density parameter for dark energy (fluid) \(\Omega_{fld, 0}\).
Omega0_g Current density parameter for photons, \(\Omega_{g,0}\).
Omega0_k Current density parameter for curvaturve, \(\Omega_{k,0}\).
Omega0_lambda Current density parameter for cosmological constant, \(\Omega_{\Lambda,0}\).
Omega0_m The sum of density parameters for all non-relativistic components, \(\Omega_{0,m}\).
Omega0_ncdm Current density parameter for distinguishable (massive) neutrinos for each species as an array, \(\Omega_{0, ncdm}\).
Omega0_ncdm_tot Current total density parameter of all distinguishable (massive) neutrinos.
Omega0_pncdm The pressure contribution to the current density parameter for the non-relativatistic part of massive neutrinos (an array holding all species).
Omega0_pncdm_tot The sum of \(\Omega_{0,pncdm}\) for all species.
Omega0_r Current density parameter of radiation, \(\Omega_{0,r}\).
Omega0_ur Current density parameter of ultra-relativistic (massless) neutrinos, \(\Omega_{0,\nu_r}\).
Onu0 Returns the sum of Omega0_ncdm_tot and Omega0_ur.
P_k_max The maximum k value measured for power spectra in \(h \mathrm{Mpc}^{-1}\).
P_k_min The minimum k value for which power spectra have been computed in \(h \mathrm{Mpc}^{-1}\).
P_z_max The input parameter specifying the maximum redshift measured for power spectra.
T0_cmb The current CMB temperature in Kelvins.
T0_ncdm An array holding the current ncdm temperature in Kelvins for each species.
Tcmb0 Returns T0_cmb.
Tnu0 Returns T0_ncdm.
a_max The maximum scale factor for which results can be computed; it can be greater than 1.0.
a_today An arbitrary number that sets the reference scaling factor.
age0 The current age of the universe in gigayears.
gauge The gauge name as a string, either ‘newtonian’ or ‘synchronous’.
h The dimensionless Hubble parameter.
has_massive_nu Returns True if N_ncdm is greater than zero.
has_pk_matter Boolean flag specifying whether matter power spectra have been requested as output.
k_max_for_pk The input parameter specifying the maximum k value to compute spectra for; units of \(h \mathrm{Mpc}^{-1}\).
k_pivot The primordial power spectrum pivot scale, where the primordial power is equal to \(A_s\).
ln_1e10_A_s Return \(\log(10^{10}A_s)\).
m_ncdm The masses of the distinguishable ncdm (massive neutrino) species, in units of eV.
n_s The tilt of the primordial power spectrum.
nonlinear Boolean flag specifying whether the power spectrum is nonlinear.
parameter_file A string holding the parameter names and values as loaded by CLASS.
rs_drag The comoving sound horizon at baryon drag, in \(\mathrm{Mpc}/h\).
rs_rec The comoving sound horizon at recombination, \(z=z_\mathrm{rec}\).
sigma8 The amplitude of matter fluctuations at \(z=0\) in a sphere of radius \(r = 8 \ h^{-1}\mathrm{Mpc}\).
tau_reio The reionization optical depth.
theta_s The sound horizon angle at recombination, equal to \(r_s(z_\mathrm{rec}) / D_a(z_\mathrm{rec})\).
w0 Returns w0_fld.
w0_fld Current fluid equation of state parameter, \(w_{0,fld}\).
wa Returns wa_fld.
wa_fld Fluid equation of state derivative, \(w_{a,fld}\).
z_drag The baryon drag redshift.
z_rec The redshift at which the visibility reaches its maximum; equals the recombination redshift.
z_reio The reionization redshift.

Methods

Ob(z) Returns Omega_b().
Ode(z) Returns the sum of Omega_lambda() and Omega_fld().
Odm(z) Returns Omega_cdm().
Ogamma(z) Returns Omega_g().
Ok(z) Returns Omega_k().
Om(z) Returns Omega_m().
Omega_b(self, z) Density parameter of baryons.
Omega_cdm(self, z) Density parameter of cold dark matter.
Omega_fld(self, z) Density parameter of dark energy (fluid).
Omega_g(self, z) Density parameter of photons.
Omega_k(self, z) Density parameter of curvature.
Omega_lambda(self, z) Density of dark energy (cosmological constant).
Omega_m(self, z) Density parameter of non-relativistic (matter-like) component, including non-relativistic part of massive neutrino.
Omega_ncdm(self, z[, species]) Density parameter of massive neutrinos.
Omega_pncdm(self, z[, species]) Return \(\Omega_{pncdm}\) as a function redshift.
Omega_r(self, z) Density parameter of relativistic (radiation-like) component, including relativistic part of massive neutrino and massless neutrino.
Omega_ur(self, z) Density parameter of ultra relativistic neutrinos.
Onu(z) Returns the sum of Omega_ncdm() and Omega_ur().
T_cmb(self, z) The CMB temperature as a function of redshift.
T_ncdm(self, z) The ncdm temperature (massive neutrinos) as a function of redshift.
Tcmb(z) Returns \((1+z)\) T0_cmb.
Tnu(z) Returns \((1+z)\) T0_ncdm.
angular_diameter_distance(self, z) Angular diameter distance in \(\mathrm{Mpc}/h\) at a given redshift.
clone(**kwargs) Create a new cosmology based on modification of self, with the input keyword parameters changed.
comoving_distance(self, z) Comoving line-of-sight distance in \(\mathrm{Mpc}/h\) at a given redshift.
comoving_transverse_distance(self, z) Comoving transverse distance in \(\mathrm{Mpc}/h\) at a given redshift.
compute_for_z(self, z, int column) Internal function to compute the background module at a specific redshift.
efunc(self, z) Function giving \(E(z)\), where the Hubble parameter is defined as \(H(z) = H_0 E(z)\).
efunc_prime(self, z) Function giving \(dE(z) / da\).
from_astropy(cosmo, **kwargs) Initialize and return a Cosmology object from a subclass of astropy.cosmology.FLRW.
from_dict(pars) Creates a Cosmology from a pars dictionary.
from_file(filename, **kwargs) Initialize a Cosmology object from the CLASS parameter file
get_pk(self, k, z) The primary power spectrum result (nonlinear if enabled) on k and z array.
get_pklin(self, k, z) Linear power spectrum result (linear even if nonlinear is enabled) on k and z array.
get_primordial(self) Return the primordial scalar and/or tensor spectrum depending on ‘modes’.
get_transfer(self, z[, output_format]) Return the density and/or velocity transfer functions for all initial conditions today.
hubble_function(self, z) The Hubble function in CLASS units, returning ba.index_bg_H.
hubble_function_prime(self, z) Derivative of Hubble function: \(dH/d\tau\), where \(d\tau/da = 1 / (a^2 H)\) in CLASS units.
luminosity_distance(self, z) Luminosity distance in \(\mathrm{Mpc}/h\) at redshift z.
match([sigma8, Omega0_cb, Omega0_m]) Creates a new cosmology that matches a derived parameter.
nu_relative_density(z) Returns Onu() / Ogamma().
p_ncdm(self, z[, species]) Pressure of non-relative part of massive neutrino.
rho_b(self, z) Density of baryons \(\rho_b\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_cdm(self, z) Density of cold dark matter \(\rho_{cdm}\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_crit(self, z) Critical density excluding curvature \(\rho_c\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_fld(self, z) Density of dark energy fluid \(\rho_{fld}\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_g(self, z) Density of photons \(\rho_g\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_k(self, z) Density of curvature \(\rho_k\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_lambda(self, z) Density of cosmological constant \(\rho_\Lambda\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_m(self, z) Density of matter \(\rho_b\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_ncdm(self, z[, species]) Density of non-relativistic part of massive neutrinos \(\rho_{ncdm}\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_r(self, z) Density of radiation \(\rho_r\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_tot(self, z) Total density \(\rho_\mathrm{tot}\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
rho_ur(self, z) Density of ultra-relativistic radiation (massless neutrinos) \(\rho_{ur}\) as a function of redshift, in units of \(10^{10} (M_\odot/h) (\mathrm{Mpc}/h)^{-3}\).
scale_independent_growth_factor(self, z) Return the scale invariant growth factor \(D(a)\) for CDM perturbations.
scale_independent_growth_rate(self, z) The scale invariant growth rate \(d\mathrm{ln}D/d\mathrm{ln}a\) for CDM perturbations.
sigma8_z(self, z) Return \(\sigma_8(z)\).
tau(self, z) Conformal time, equal to comoving distance when K = 0.0 (flat universe).
time(self, z) Proper time (age of universe) in gigayears.
to_astropy() Initialize and return a subclass of astropy.cosmology.FLRW from the Cosmology class.
Omega0_cb

The total density of CDM and Baryon.

This is not an input CLASS parameter. To scale Omega0_cb, use match().

__dir__()[source]

a list of all members from all delegate classes

__getattr__(name)[source]

Find the proper delegate, initialize it, and run the method

__iter__()[source]

Allows dict() to be used on class. Use from_dict() to reconstruct an instance.

__str__()[source]

Return a dict string when printed

clone(**kwargs)[source]

Create a new cosmology based on modification of self, with the input keyword parameters changed.

Parameters:**kwargs – keyword parameters to adjust
Returns:a copy of self, with the input kwargs adjusted
Return type:Cosmology
classmethod from_astropy(cosmo, **kwargs)[source]

Initialize and return a Cosmology object from a subclass of astropy.cosmology.FLRW.

Parameters:
  • cosmo (subclass of astropy.cosmology.FLRW.) – the astropy cosmology instance
  • **kwargs – extra keyword parameters to pass when initializing; they shall not be in conflict with the parameters inferred from cosmo. To override parameters, chain the result with clone().
Returns:

the initialized cosmology object

Return type:

Cosmology

classmethod from_dict(pars)[source]

Creates a Cosmology from a pars dictionary.

This is a rather ‘raw’ API. The dictionary must be readable by ClassEngine. Unlike Cosmology(**args), pars must not contain any convenient names defined here.

classmethod from_file(filename, **kwargs)[source]

Initialize a Cosmology object from the CLASS parameter file

Parameters:
  • filename (str) – the name of the parameter file to read
  • **kwargs – extra keyword parameters to pass when initializing; they shall not be in conflict with the parameters inferred from cosmo. To override parameters, chain the result with clone().
match(sigma8=None, Omega0_cb=None, Omega0_m=None)[source]

Creates a new cosmology that matches a derived parameter. This is different from clone, where CLASS parameters are used.

Note that we only supoort matching one derived parameter at a time, because the matching is in general non-commutable.

Parameters:
  • sigma8 (float or None) – We scale the scalar amplitude A_s to achieve the desired sigma8.
  • Omega0_cb (float or None) – Desired total energy density of CDM and baryon.
  • Omega0_m (float or None) – Desired total energy density of matter-like components (included ncdm)
Returns:

Return type:

A new cosmology parameter where the derived parameter matches the given constrain.

sigma8

The amplitude of matter fluctuations at \(z=0\) in a sphere of radius \(r = 8 \ h^{-1}\mathrm{Mpc}\).

This is not an input CLASS parameter. To scale sigma8, use match(), which adjusts scalar amplitude A_s to achieve the desired sigma8.

to_astropy()[source]

Initialize and return a subclass of astropy.cosmology.FLRW from the Cosmology class.

Returns:the astropy class holding the cosmology values
Return type:subclass of astropy.cosmology.FLRW