API

Cosmology

Cosmology class

cosmoprimo.cosmology.Background(cosmology, engine=None, set_engine=True, **extra_params)

Return Background calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

class cosmoprimo.cosmology.BaseBackground(engine)

Bases: BaseSection

Base background engine, including a few definitions.

Omega_Lambda(z)

Density of cosmological constant, unitless.

Omega_b(z)

Density parameter of baryons, unitless.

Omega_cdm(z)

Density parameter of cold dark matter, unitless.

Omega_de(z)

Density of total dark energy (fluid + cosmological constant), unitless.

Omega_fld(z)

Density of cosmological constant, unitless.

Omega_g(z)

Density parameter of photons, unitless.

Omega_k(z)

Density parameter of curvature, unitless.

Omega_m(z)

Density parameter of non-relativistic (matter-like) component, including non-relativistic part of massive neutrino, unitless.

Omega_ncdm(z, species=None)

Density parameter of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_ncdm_tot(z)

Total density parameter of massive neutrinos, unitless.

Omega_pncdm(z, species=None)

Density parameter of pressure of non-relativistic part of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_pncdm_tot(z)

Total density parameter of pressure of non-relativistic part of massive neutrinos, unitless.

Omega_r(z)

Density parameter of radiation, including photons and relativistic part of massive and massless neutrinos, unitless.

Omega_ur(z)

Density parameter of massless neutrinos, unitless.

T_cmb(z)

The CMB temperature, in \(K\).

T_ncdm(z, species=None)

Return the ncdm temperature (massive neutrinos), in \(K\). Returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)).

angular_diameter_distance(z)

Proper angular diameter distance, in \(\mathrm{Mpc}/h\).

See eq. 18 of astro-ph/9905116 for \(D_{A}(z)\).

angular_diameter_distance_2(z1, z2)

Angular diameter distance of object at \(z_{2}\) as seen by observer at \(z_{1}\), that is, \(S_{K}((\chi(z_{2}) - \chi(z_{1})) \sqrt{|K|}) / \sqrt{|K|} / (1 + z_{2})\), where \(S_{K}\) is the identity if \(K = 0\), \(\sin\) if \(K < 0\) and \(\sinh\) if \(K > 0\). camb’s angular_diameter_distance2(z1, z2) is not used as it returns 0 when z2 < z1.

comoving_angular_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

comoving_transverse_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

efunc(z)

Function giving \(E(z)\), where the Hubble parameter is defined as \(H(z) = H_{0} E(z)\), unitless.

hubble_function(z)

Hubble function ba.index_bg_H, in \(\mathrm{km}/\mathrm{s}/\mathrm{Mpc}\).

luminosity_distance(z)

Luminosity distance, in \(\mathrm{Mpc}/h\).

See eq. 21 of astro-ph/9905116 for \(D_{L}(z)\).

p_ncdm(z, species=None)

Pressure of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return pressure for this species.

p_ncdm_tot(z)

Total pressure of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_Lambda(z)

Comoving density of cosmological constant \(\rho_{\Lambda}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_b(z)

Comoving density of baryons \(\rho_{b}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_cdm(z)

Comoving density of cold dark matter \(\rho_{cdm}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_crit(z)

Comoving critical density excluding curvature \(\rho_{c}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

This is defined as:

\[\rho_{\mathrm{crit}}(z) = \frac{3 H(z)^{2}}{8 \pi G}.\]
rho_de(z)

Total comoving density of dark energy \(\rho_{\mathrm{de}}\) (fluid + cosmological constant), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_fld(z)

Comoving density of dark energy fluid \(\rho_{\mathrm{fld}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_g(z)

Comoving density of photons \(\rho_{g}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_k(z)

Comoving density of curvature \(\rho_{k}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_m(z)

Comoving density of matter \(\rho_{m}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ncdm(z, species=None)

Comoving density of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

rho_ncdm_tot(z)

Total comoving density of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_r(z)

Comoving density of radiation \(\rho_{r}\), including photons and relativistic part of massive and massless neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_tot(z)

Comoving total density \(\rho_{\mathrm{tot}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ur(z)

Comoving density of massless neutrinos \(\rho_{ur}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

class cosmoprimo.cosmology.BaseCosmoParams

Bases: BaseClass

get(*args, **kwargs)

Return an input (or easily derived) parameter.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

class cosmoprimo.cosmology.BaseEngine(cosmo, **extra_params)

Bases: BaseCosmoParams

Base engine for cosmological calculation.

Initialize engine.

Parameters:
  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Engine parameters.

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background()

Return Background calculations.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier()

Return Fourier calculations.

get_harmonic()

Return Harmonic calculations.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations()

Return Perturbations calculations.

get_primordial()

Return Primordial calculations.

get_thermodynamics()

Return Thermodynamics calculations.

get_transfer()

Return Transfer calculations.

class cosmoprimo.cosmology.BaseSection(engine)

Bases: object

Base section.

class cosmoprimo.cosmology.Cosmology(engine=None, extra_params=None, **params)

Bases: BaseCosmoParams

Cosmology, defined as a set of parameters (and possibly a current engine attached to it).

Initialize Cosmology.

Note

If Omega_m (or omega_m) is provided, Omega_cdm is infered by subtracting Omega_b and the non-relativistic part of Omega_ncdm from Omega_m. Massive neutrinos can be provided e.g. through m_ncdm or Omega_ncdm/omega_ncdm with their temperatures w.r.t. CMB T_ncdm_over_cmb. In the case of Omega_ncdm, the neutrino energy density (see _compute_ncdm_momenta()) will be inverted to recover m_ncdm. If a single value for m_ncdm or Omega_ncdm is provided, neutrino_hierarchy can be set to None (default, single massive neutrino) or ‘normal’, ‘inverted’, ‘degenerate’ (all neutrinos with same mass), which will determine the masses of the 3 neutrinos. If the number of relativistic species N_ur is not provided (or None), it will be determined from the desired effective number of neutrinos N_eff (typically kept at 3.044 for 3 neutrinos whatever m_ncdm or Omega_ncdm/omega_ncdm) and the number of massless neutrinos (\(m \leq 0.00017 \; \mathrm{eV}\)), which are then removed from the list m_ncdm. Parameter Omega_ncdm/omega_ncdm (accessed as cosmo['Omega_ncdm']/cosmo['omega_ncdm']) will always provide the total energy density of neutrinos (single value). The pivot scale k_pivot is in \(\mathrm{Mpc}^{-1}\). If ‘non_linear’ is required, we recommend “mead” for class and camb matching.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, no engine is set.

  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Cosmological and calculation parameters which take priority over the default ones.

clone(base='input', engine=None, extra_params=None, **params)

Clone current cosmology instance, optionally updating engine and parameters.

Parameters:
  • base (string, default='input') – If ‘internal’ or None, update parameters in the internal \(h, \Omega, m_{cdm}\) basis. If, e.g. input parameters are \(h, \omega_{b}, \omega_{cdm}\), clone(base='internal', h=0.7) returns the same cosmology, but with \(h = 0.7\); since \(\Omega_{b}, \Omega_{cdm}\) are kept fixed, \(\omega_{b}, \omega_{cdm}\) are modified. If ‘input’, update input parameters. With clone(base='input', h=0.7) \(\omega_{b}, \omega_{cdm}\) are left unchanged, but \(\Omega_{b}, \Omega_{cdm}\) are modified.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, use same engine (class) as current instance.

  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Cosmological and calculation parameters which take priority over the current ones.

Returns:

new – Copy of current instance, with updated engine and parameters.

Return type:

Cosmology

classmethod from_state(state)

Instantiate and initalise class with state dictionary.

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background(engine=None, set_engine=True, **extra_params)

Get background.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier(engine=None, set_engine=True, **extra_params)

Get fourier.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

get_harmonic(engine=None, set_engine=True, **extra_params)

Get harmonic.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='input') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations(engine=None, set_engine=True, **extra_params)

Get perturbations.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

get_primordial(engine=None, set_engine=True, **extra_params)

Get primordial.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

get_thermodynamics(engine=None, set_engine=True, **extra_params)

Get thermodynamics.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

get_transfer(engine=None, set_engine=True, **extra_params)

Get transfer.

Parameters:
  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_variants’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

classmethod load(filename)

Load class from disk.

save(filename)

Save class to disk.

set_engine(engine, set_engine=True, **extra_params)

Set engine for cosmological calculation.

Parameters:
  • engine (string) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’].

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology. (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

solve(param, func, target=0.0, limits=None, xtol=1e-06, rtol=1e-06, maxiter=100)

Return cosmology cosmo that verifies func(cosmo) == target, by varying parameter param.

Parameters:
  • param (string) – Input parameter name, e.g. ‘h’.

  • func (callable, string) – Function that takes a Cosmology instance (clone of self) and returns a value, e.g. lambda cosmo: cosmo.get_thermodynamics().theta_star. If ‘theta_MC_100’, match 100 * cosmo['theta_cosmomc'] to target (and engine should be defined).

  • target (float, default=0.) – Target value.

  • limits (tuple, list, default=None) – Variation range for param.

  • xtol (float, default=1e-6) – Absolute tolerance on the value of param. See scipy.optimize.bisect().

  • rtol (float, default=1e-6) – Relative tolerance on the value of param. See scipy.optimize.bisect().

  • maxiter (int, default=100) – If convergence is not achieved in maxiter iterations, an error is raised. Must be >= 0.

Returns:

new

Return type:

Cosmology

exception cosmoprimo.cosmology.CosmologyComputationError

Bases: CosmologyError

Exception raised when error in cosmology computation.

with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

exception cosmoprimo.cosmology.CosmologyError

Bases: Exception

Exception raised by Cosmology.

with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

exception cosmoprimo.cosmology.CosmologyInputError

Bases: CosmologyError

Exception raised when error in the value of input parameters.

with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

class cosmoprimo.cosmology.DefaultBackground(engine)

Bases: BaseBackground

Omega_Lambda(z)

Density of cosmological constant, unitless.

Omega_b(z)

Density parameter of baryons, unitless.

Omega_cdm(z)

Density parameter of cold dark matter, unitless.

Omega_de(z)

Density of total dark energy (fluid + cosmological constant), unitless.

Omega_fld(z)

Density of cosmological constant, unitless.

Omega_g(z)

Density parameter of photons, unitless.

Omega_k(z)

Density parameter of curvature, unitless.

Omega_m(z)

Density parameter of non-relativistic (matter-like) component, including non-relativistic part of massive neutrino, unitless.

Omega_ncdm(z, species=None)

Density parameter of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_ncdm_tot(z)

Total density parameter of massive neutrinos, unitless.

Omega_pncdm(z, species=None)

Density parameter of pressure of non-relativistic part of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_pncdm_tot(z)

Total density parameter of pressure of non-relativistic part of massive neutrinos, unitless.

Omega_r(z)

Density parameter of radiation, including photons and relativistic part of massive and massless neutrinos, unitless.

Omega_ur(z)

Density parameter of massless neutrinos, unitless.

T_cmb(z)

The CMB temperature, in \(K\).

T_ncdm(z, species=None)

Return the ncdm temperature (massive neutrinos), in \(K\). Returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)).

property age

The current age of the Universe, in \(\mathrm{Gy}\).

angular_diameter_distance(z)

Proper angular diameter distance, in \(\mathrm{Mpc}/h\).

See eq. 18 of astro-ph/9905116 for \(D_{A}(z)\).

angular_diameter_distance_2(z1, z2)

Angular diameter distance of object at \(z_{2}\) as seen by observer at \(z_{1}\), that is, \(S_{K}((\chi(z_{2}) - \chi(z_{1})) \sqrt{|K|}) / \sqrt{|K|} / (1 + z_{2})\), where \(S_{K}\) is the identity if \(K = 0\), \(\sin\) if \(K < 0\) and \(\sinh\) if \(K > 0\). camb’s angular_diameter_distance2(z1, z2) is not used as it returns 0 when z2 < z1.

comoving_angular_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

comoving_radial_distance(z)

Comoving radial distance, in \(mathrm{Mpc}/h\).

See eq. 15 of astro-ph/9905116 for \(D_C(z)\).

comoving_transverse_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

efunc(z)

Function giving \(E(z)\), where the Hubble parameter is defined as \(H(z) = H_{0} E(z)\), unitless.

hubble_function(z)

Hubble function ba.index_bg_H, in \(\mathrm{km}/\mathrm{s}/\mathrm{Mpc}\).

luminosity_distance(z)

Luminosity distance, in \(\mathrm{Mpc}/h\).

See eq. 21 of astro-ph/9905116 for \(D_{L}(z)\).

p_ncdm(z, species=None)

Pressure of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return pressure for this species.

p_ncdm_tot(z)

Total pressure of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_Lambda(z)

Comoving density of cosmological constant \(\rho_{\Lambda}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_b(z)

Comoving density of baryons \(\rho_{b}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_cdm(z)

Comoving density of cold dark matter \(\rho_{cdm}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_crit(z)

Comoving critical density excluding curvature \(\rho_{c}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

This is defined as:

\[\rho_{\mathrm{crit}}(z) = \frac{3 H(z)^{2}}{8 \pi G}.\]
rho_de(z)

Total comoving density of dark energy \(\rho_{\mathrm{de}}\) (fluid + cosmological constant), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_fld(z)

Comoving density of dark energy fluid \(\rho_{\mathrm{fld}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_g(z)

Comoving density of photons \(\rho_{g}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_k(z)

Comoving density of curvature \(\rho_{k}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_m(z)

Comoving density of matter \(\rho_{m}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ncdm(z, species=None)

Comoving density of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

rho_ncdm_tot(z)

Total comoving density of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_r(z)

Comoving density of radiation \(\rho_{r}\), including photons and relativistic part of massive and massless neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_tot(z)

Comoving total density \(\rho_{\mathrm{tot}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ur(z)

Comoving density of massless neutrinos \(\rho_{ur}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

time(z)

Proper time (age of universe), in \(\mathrm{Gy}\).

cosmoprimo.cosmology.Fourier(cosmology, engine=None, set_engine=True, **extra_params)

Return Fourier calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

cosmoprimo.cosmology.Harmonic(cosmology, engine=None, set_engine=True, **extra_params)

Return Harmonic calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

class cosmoprimo.cosmology.MetaSection(name, bases, class_dict)

Bases: type

Metaclass registering BaseEngine-derived classes.

mro()

Return a type’s method resolution order.

cosmoprimo.cosmology.Perturbations(cosmology, engine=None, set_engine=True, **extra_params)

Return Perturbations calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

cosmoprimo.cosmology.Primordial(cosmology, engine=None, set_engine=True, **extra_params)

Return Primordial calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

class cosmoprimo.cosmology.RegisteredEngine(name, bases, class_dict)

Bases: type

Metaclass registering BaseEngine-derived classes.

mro()

Return a type’s method resolution order.

cosmoprimo.cosmology.Thermodynamics(cosmology, engine=None, set_engine=True, **extra_params)

Return Thermodynamics calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

cosmoprimo.cosmology.Transfer(cosmology, engine=None, set_engine=True, **extra_params)

Return Transfer calculations.

Parameters:
  • cosmology (Cosmology) – Current cosmology.

  • engine (string, default=None) – Engine name, one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘bbks’]. If None, returns current Cosmology.engine.

  • set_engine (bool, default=True) – Whether to attach returned engine to cosmology (Set False if e.g. you want to use this engine for a single calculation).

  • extra_params (dict) – Extra engine parameters, typically precision parameters.

Returns:

engine

Return type:

BaseEngine

cosmoprimo.cosmology.check_params(args, **kwargs)

Check for conflicting parameters in args parameter dictionary.

class cosmoprimo.cosmology.class_or_instancemethod

Bases: classmethod

cosmoprimo.cosmology.compute_ncdm_momenta(T_eff, m, z, method='laguerre', epsabs=1e-07, epsrel=1e-07, out='rho')

Return momenta of non-CDM components (massive neutrinos) by integrating over the phase-space distribution (frozen since CMB).

Parameters:
  • T_eff (float) – Effective temperature; typically T_cmb * T_ncdm_over_cmb.

  • m (float) – Mass in \(\mathrm{eV}\).

  • z (float, array) – Redshift.

  • epsrel (float, default=1e-7) – Relative precision (for scipy.integrate.quad() integration).

  • out (string, default='rho') – If ‘rho’, return energy density. If ‘drhodm’, return derivative of energy density w.r.t. to mass m. If ‘p’, return pressure.

Returns:

out – For each input redshift, required momentum, in units of \(10^{10} M_{\odot} / \mathrm{Mpc}^{3}\) (/ \(\mathrm{eV}\) if out is ‘drhodm’)

Return type:

float, array

cosmoprimo.cosmology.find_conflicts(name, conflicts=())

Return conflicts corresponding to input parameter name.

Parameters:

name (string) – Parameter name.

Returns:

conflicts – Conflicting parameter names.

Return type:

tuple

cosmoprimo.cosmology.get_engine(engine)

Return engine (class) for cosmological calculation.

Parameters:

engine (type, string) – Engine or one of [‘class’, ‘camb’, ‘eisenstein_hu’, ‘eisenstein_hu_nowiggle’, ‘eisenstein_hu_nowiggle_variants’, ‘bbks’].

Returns:

engine

Return type:

BaseEngine

cosmoprimo.cosmology.merge_params(args, moreargs, **kwargs)

Merge moreargs parameters into args. moreargs parameters take priority over those defined in args.

Note

args is modified in-place.

Parameters:
  • args (dict) – Base parameter dictionary.

  • moreargs (dict) – Parameter dictionary to be merged into args.

Returns:

args – Merged parameter dictionary.

Return type:

dict

Engines

CLASS

Cosmological calculation with the Boltzmann code CLASS.

class cosmoprimo.classy.Background(*args: Any, **kwargs: Any)

Bases: BaseClassBackground, Background

class cosmoprimo.classy.ClassEngine(*args, **kwargs)

Bases: BaseEngine

Engine for the Boltzmann code CLASS.

Initialize engine.

Parameters:
  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Engine parameters.

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background()

Return Background calculations.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier()

Return Fourier calculations.

get_harmonic()

Return Harmonic calculations.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations()

Return Perturbations calculations.

get_primordial()

Return Primordial calculations.

get_thermodynamics()

Return Thermodynamics calculations.

get_transfer()

Return Transfer calculations.

class cosmoprimo.classy.Fourier(*args: Any, **kwargs: Any)

Bases: BaseClassFourier, Fourier

pk_interpolator(non_linear=False, of='delta_m', **kwargs)

Return PowerSpectrumInterpolator2D instance.

Parameters:
  • non_linear (bool, default=False) – Whether to return the non_linear power spectrum (if requested in parameters, with ‘non_linear’: ‘halofit’ or ‘mead’). Computed only for of == 'delta_m' or ‘delta_cb’.

  • of (string, tuple, default='delta_m') – Perturbed quantities. Either ‘delta_m’ for matter perturbations or ‘delta_cb’ for cold dark matter + baryons perturbations will use precomputed spectra. Else, e.g. (‘delta_m’, ‘theta_cb’) for the cross matter density - cold dark matter + baryons velocity power spectra, are computed on-the-fly.

  • kwargs (dict) – Arguments for PowerSpectrumInterpolator2D.

pk_kz(k, z, non_linear=False, of='m')

Return power spectrum, in \((\mathrm{Mpc}/h)^{3}\), using original CLASS routine.

Parameters:
  • k (array_like) – Wavenumbers, in \(h/\mathrm{Mpc}\).

  • z (array_like) – Redshifts.

  • non_linear (bool, default=False) – Whether to return the non_linear power spectrum (if requested in parameters, with ‘non_linear’: ‘halofit’ or ‘mead’).

  • of (string, default='delta_m') – Perturbed quantities. Either ‘delta_m’ for matter perturbations or ‘delta_cb’ for cold dark matter + baryons perturbations.

Returns:

pk – Power spectrum array of shape (len(k), len(z)).

Return type:

array

property sigma8_cb

Current r.m.s. of cold dark matter + baryons perturbations in a sphere of \(8 \mathrm{Mpc}/h\) unitless.

property sigma8_m

Current r.m.s. of matter perturbations in a sphere of \(8 \mathrm{Mpc}/h\), unitless.

sigma8_z(z, of='delta_m')

Return the r.m.s. of of perturbations in sphere of \(8 \mathrm{Mpc}/h\).

sigma_rz(r, z, of='delta_m', **kwargs)

Return the r.m.s. of of perturbations in sphere of \(r \mathrm{Mpc}/h\).

table(non_linear=False, of='delta_m')

Return power spectrum table, in \((\mathrm{Mpc}/h)^{3}\).

Parameters:
  • non_linear (bool, default=False) – Whether to return the non_linear power spectrum (if requested in parameters, with ‘non_linear’: ‘halofit’ or ‘mead’). Computed only for of == 'delta_m' or ‘delta_cb’.

  • of (string, tuple, default='delta_m') – Perturbed quantities. Either ‘delta_m’ for matter perturbations or ‘delta_cb’ for cold dark matter + baryons perturbations will use precomputed spectra. Else, e.g. (‘delta_m’, ‘theta_cb’) for the cross matter density - cold dark matter + baryons velocity power spectra, are computed on-the-fly.

Returns:

  • k (array) – Wavenumbers.

  • z (array) – Redshifts.

  • pk (array) – Power spectrum array of shape (len(k), len(z)).

class cosmoprimo.classy.Harmonic(*args: Any, **kwargs: Any)

Bases: BaseClassHarmonic, Harmonic

lensed_table(ellmax=-1, of=None)

Return table of lensed \(C_{\ell}\), unitless.

Parameters:
  • ellmax (int, default=-1) – Maximum \(\ell\) desired. If negative, is relative to the requested maximum \(\ell\).

  • of (list, default=None) – List of outputs, [‘tt’, ‘ee’, ‘bb’, ‘pp’, ‘te’, ‘tp’]. If None, return all computed outputs.

Returns:

cell – Structured array.

Return type:

array

unlensed_table(ellmax=-1, of=None)

Return table of unlensed \(C_{\ell}\) (i.e. CMB power spectra without lensing and lensing potentials), unitless.

Parameters:
  • ellmax (int, default=-1) – Maximum \(\ell\) desired. If negative, is relative to the requested maximum \(\ell\).

  • of (list, default=None) – List of outputs, [‘tt’, ‘ee’, ‘bb’, ‘te’, ‘pp’, ‘tp’, ‘ep’]. If None, return all computed outputs.

Returns:

cell – Structured array.

Return type:

array

Note

Normalisation is \(C_{\ell}\) rather than \(\ell(\ell+1)C_{\ell}/(2\pi)\) (or \(\ell^{2}(\ell+1)^{2}/(2\pi)\) in the case of the lensing potential pp spectrum). Usually multiplied by CMB temperature in \(\mu K\).

class cosmoprimo.classy.Perturbations(*args: Any, **kwargs: Any)

Bases: BaseClassPerturbations, Perturbations

class cosmoprimo.classy.Primordial(*args: Any, **kwargs: Any)

Bases: BaseClassPrimordial, Primordial

property A_s

Scalar amplitude of the primordial power spectrum at \(k_\mathrm{pivot}\), unitless.

property ln_1e10_A_s

\(\ln(10^{10}A_s)\), unitless.

pk_interpolator(mode='scalar')

Return power spectrum interpolator.

Parameters:

mode (string, default='scalar') – ‘scalar’, ‘vector’ or ‘tensor’ mode.

Returns:

interpPowerSpectrumInterpolator1D instance if only one type of initial conditions (typically adiabatic), else dictionary of class:PowerSpectrumInterpolator1D corresponding to the tuples of initial conditions.

Return type:

PowerSpectrumInterpolator1D, dict

pk_k(k, mode='scalar')

The primordial spectrum of curvature perturbations at k, generated by inflation, in \((\mathrm{Mpc}/h)^{3}\). For scalar perturbations this is e.g. defined as:

\[\mathcal{P_R}(k) = A_s \left (\frac{k}{k_\mathrm{pivot}} \right )^{n_s - 1 + 1/2 \alpha_s \ln(k/k_\mathrm{pivot})}\]

See also: eq. 2 of this reference.

Parameters:
  • k (array_like) – Wavenumbers, in \(h/\mathrm{Mpc}\).

  • mode (string, default='scalar') – ‘scalar’, ‘vector’ or ‘tensor’ mode.

Returns:

pk – The primordial power spectrum if only one type of initial conditions (typically adiabatic), else dictionary of primordial power spectra corresponding to the tuples of initial conditions.

Return type:

array, dict

table()

Return primordial table.

Returns:

data – Structured array containing primordial data.

Return type:

array

class cosmoprimo.classy.Thermodynamics(*args: Any, **kwargs: Any)

Bases: BaseClassThermodynamics, Thermodynamics

class cosmoprimo.classy.Transfer(*args: Any, **kwargs: Any)

Bases: BaseClassTransfer, Transfer

CAMB

Cosmological calculation with the Boltzmann code CAMB.

class cosmoprimo.camb.Background(engine)

Bases: BaseBackground

Omega_Lambda(z)

Density of cosmological constant, unitless.

Omega_b(z)

Density parameter of baryons, unitless.

Omega_cdm(z)

Density parameter of cold dark matter, unitless.

Omega_de(z)

Total density of dark energy (fluid + cosmological constant), unitless.

Omega_fld(z)

Density of cosmological constant, unitless.

Omega_g(z)

Density parameter of photons, unitless.

Omega_k(z)

Density parameter of curvature, unitless.

Omega_m(z)

Density parameter of non-relativistic (matter-like) component, including non-relativistic part of massive neutrino, unitless.

Omega_ncdm(z, species=None)

Density parameter of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_ncdm_tot(z)

Total density parameter of massive neutrinos, unitless.

Omega_pncdm(z, species=None)

Density parameter of pressure of non-relativistic part of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_pncdm_tot(z)

Total density parameter of pressure of non-relativistic part of massive neutrinos, unitless.

Omega_r(z)

Density parameter of radiation, including photons and relativistic part of massive and massless neutrinos, unitless.

Omega_ur(z)

Density parameter of massless neutrinos, unitless.

T_cmb(z)

The CMB temperature, in \(K\).

T_ncdm(z, species=None)

Return the ncdm temperature (massive neutrinos), in \(K\). Returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)).

property age

The current age of the Universe, in \(\mathrm{Gy}\).

angular_diameter_distance(z)

Proper angular diameter distance, in \(\mathrm{Mpc}/h\).

See eq. 18 of astro-ph/9905116 for \(D_{A}(z)\).

angular_diameter_distance_2(z1, z2)

Angular diameter distance of object at \(z_{2}\) as seen by observer at \(z_{1}\), that is, \(S_{K}((\chi(z_{2}) - \chi(z_{1})) \sqrt{|K|}) / \sqrt{|K|} / (1 + z_{2})\), where \(S_{K}\) is the identity if \(K = 0\), \(\sin\) if \(K < 0\) and \(\sinh\) if \(K > 0\). camb’s angular_diameter_distance2(z1, z2) is not used as it returns 0 when z2 < z1.

comoving_angular_distance(z)

Comoving angular distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

comoving_radial_distance(z)

Comoving radial distance, in \(mathrm{Mpc}/h\).

See eq. 15 of astro-ph/9905116 for \(D_C(z)\).

comoving_transverse_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

efunc(z)

Function giving \(E(z)\), where the Hubble parameter is defined as \(H(z) = H_{0} E(z)\), unitless.

hubble_function(z)

Hubble function, in \(\mathrm{km}/\mathrm{s}/\mathrm{Mpc}\).

luminosity_distance(z)

Luminosity distance, in \(\mathrm{Mpc}/h\).

See eq. 21 of astro-ph/9905116 for \(D_{L}(z)\).

p_ncdm(z, species=None)

Pressure of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return pressure for this species.

p_ncdm_tot(z)

Total pressure of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_Lambda(z)

Comoving density of cosmological constant \(\rho_{\Lambda}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_b(z)

Density parameter of baryons, unitless.

rho_cdm(z)

Comoving density of cold dark matter \(\rho_{cdm}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_crit(z)

Comoving critical density excluding curvature \(\rho_{c}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

This is defined as:

\[\rho_{\mathrm{crit}}(z) = \frac{3 H(z)^{2}}{8 \pi G}.\]
rho_de(z)

Comoving total density of dark energy \(\rho_{\mathrm{de}}\) (fluid + cosmological constant), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_fld(z)

Comoving density of dark energy fluid \(\rho_{\mathrm{fld}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_g(z)

Comoving density of photons \(\rho_{g}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_k(z)

Comoving density of curvature \(\rho_{k}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_m(z)

Comoving density of matter \(\rho_{m}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ncdm(z, species=None)

Comoving density of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

rho_ncdm_tot(z)

Total comoving density of non-relativistic part of massive neutrinos \(\rho_{ncdm}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_r(z)

Comoving density of radiation \(\rho_{r}\), including photons and relativistic part of massive and massless neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_tot(z)

Comoving total density \(\rho_{\mathrm{tot}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ur(z)

Comoving density of massless neutrinos \(\rho_{ur}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

time(z)

Proper time (age of universe), in \(\mathrm{Gy}\).

class cosmoprimo.camb.CambEngine(*args, **kwargs)

Bases: BaseEngine

Engine for the Boltzmann code CAMB.

Initialize engine.

Parameters:
  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Engine parameters.

compute(tasks)

The main function, which executes the desired modules.

Parameters:

tasks (list, string) – Calculation to perform, in the following list: [‘background’, ‘thermodynamics’, ‘transfer’, ‘harmonic’, ‘lensing’, ‘fourier’]

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background()

Return Background calculations.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier()

Return Fourier calculations.

get_harmonic()

Return Harmonic calculations.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations()

Return Perturbations calculations.

get_primordial()

Return Primordial calculations.

get_thermodynamics()

Return Thermodynamics calculations.

get_transfer()

Return Transfer calculations.

class cosmoprimo.camb.Fourier(engine)

Bases: BaseSection

pk_interpolator(non_linear=False, of='delta_m', **kwargs)

Return PowerSpectrumInterpolator2D instance.

Parameters:
  • non_linear (bool, default=False) – Whether to return the non_linear power spectrum (if requested in parameters, with ‘non_linear’: ‘halofit’ or ‘mead’). Computed only for of == ‘delta_m’.

  • of (string, tuple, default='delta_m') – Perturbed quantities.

  • kwargs (dict) – Arguments for PowerSpectrumInterpolator2D.

pk_kz(k, z, non_linear=False, of='delta_m')

Return power spectrum, in \((\mathrm{Mpc}/h)^{3}\).

Parameters:
  • k (array_like) – Wavenumbers, in \(h/\mathrm{Mpc}\).

  • z (array_like) – Redshifts.

  • non_linear (bool, default=False) – Whether to return the non_linear power spectrum (if requested in parameters, with ‘non_linear’: ‘halofit’ or ‘mead’).

  • of (string, default='delta_m') – Perturbed quantities.

Returns:

pk – Power spectrum array of shape (len(k),len(z)).

Return type:

array

property sigma8_m

Current r.m.s. of matter perturbations in a sphere of \(8 \mathrm{Mpc}/h\), unitless.

sigma8_z(z, of='delta_m')

Return the r.m.s. of of perturbations in sphere of \(8 \mathrm{Mpc}/h\).

sigma_rz(r, z, of='delta_m', **kwargs)

Return the r.m.s. of of perturbations in sphere of \(r \mathrm{Mpc}/h\).

table(non_linear=False, of='delta_m')

Return power spectrum table, in \((\mathrm{Mpc}/h)^{3}\).

Parameters:
  • non_linear (bool, default=False) – Whether to return the non_linear power spectrum (if requested in parameters, with ‘non_linear’: ‘halofit’ or ‘mead’). Computed only for of == ‘delta_m’ or ‘delta_cb’.

  • of (string, tuple, default='delta_m') – Perturbed quantities. ‘delta_m’ for matter perturbations, ‘delta_cb’ for cold dark matter + baryons, ‘phi’, ‘psi’ for Bardeen potentials, or ‘phi_plus_psi’ for Weyl potential. Provide a tuple, e.g. (‘delta_m’, ‘theta_cb’) for the cross matter density - cold dark matter + baryons velocity power spectra.

Returns:

  • k (numpy.ndarray) – Wavenumbers.

  • z (numpy.ndarray) – Redshifts.

  • pk (numpy.ndarray) – Power spectrum array of shape (len(k), len(z)).

class cosmoprimo.camb.Harmonic(engine)

Bases: BaseSection

lens_potential_cl(ellmax=-1)

Return potential \(C_{\ell}\) [‘pp’, ‘tp’, ‘ep’], unitless.

lensed_cl(ellmax=-1)

Return lensed \(C_{\ell}\) [‘tt’, ‘ee’, ‘bb’, ‘te’], unitless.

unlensed_cl(ellmax=-1)

Return unlensed \(C_{\ell}\) [‘tt’, ‘ee’, ‘bb’, ‘te’], unitless.

class cosmoprimo.camb.Primordial(engine)

Bases: BaseSection

property A_s

Scalar amplitude of the primordial power spectrum at \(k_\mathrm{pivot}\), unitless.

property alpha_s

Running of the scalar spectral index at \(k_\mathrm{pivot}\), unitless.

property alpha_t

Running of the tensor spectral index at \(k_\mathrm{pivot}\), unitless.

property beta_s

Running of the running of the scalar spectral index at \(k_\mathrm{pivot}\), unitless.

property k_pivot

Primordial power spectrum pivot scale, where the primordial power is equal to \(A_{s}\), in \(h/\mathrm{Mpc}\).

property ln_1e10_A_s

\(\ln(10^{10}A_s)\), unitless.

property n_s

Power-law scalar index i.e. tilt of the primordial scalar power spectrum, unitless.

property n_t

Power-law tensor index i.e. tilt of the tensor primordial power spectrum, unitless.

pk_interpolator(mode='scalar')

Return power spectrum interpolator.

Parameters:

mode (string, default='scalar') – ‘scalar’, ‘vector’ or ‘tensor’ mode.

Returns:

interp

Return type:

PowerSpectrumInterpolator1D

pk_k(k, mode='scalar')

The primordial spectrum of curvature perturbations at k, generated by inflation, in \((\mathrm{Mpc}/h)^{3}\). For scalar perturbations this is e.g. defined as:

\[\mathcal{P_R}(k) = A_s \left (\frac{k}{k_\mathrm{pivot}} \right )^{n_s - 1 + 1/2 \alpha_s \ln(k/k_\mathrm{pivot}) + 1/6 \beta_s \ln(k/k_\mathrm{pivot})^2}\]

See also: eq. 2 of this reference.

Parameters:
  • k (array_like) – Wavenumbers, in \(h/\mathrm{Mpc}\).

  • mode (string, default='scalar') – ‘scalar’, ‘vector’ or ‘tensor’ mode.

Returns:

pk – The primordial power spectrum.

Return type:

array

property r

Tensor-to-scalar power spectrum ratio at \(k_\mathrm{pivot}\), unitless.

class cosmoprimo.camb.Thermodynamics(engine)

Bases: BaseSection

rs_z(z)

Comoving sound horizon.

class cosmoprimo.camb.Transfer(engine)

Bases: BaseSection

table()

Return source functions (in array of shape (k.size, z.size)).

TODO: proper matching with CLASS, see https://github.com/lesgourg/class_public/blob/997d1ac0b64d11439948a0cc13f719ef427f87be/source/perturbations.c#L516

cosmoprimo.camb.enum(*sequential, **named)

Enumeration values to serve as ready flags.

Eisenstein & Hu

class cosmoprimo.eisenstein_hu.Background(engine)

Bases: DefaultBackground

Background quantities.

Note

Does not treat neutrinos.

Omega_Lambda(z)

Density of cosmological constant, unitless.

Omega_b(z)

Density parameter of baryons, unitless.

Omega_cdm(z)

Density parameter of cold dark matter, unitless.

Omega_de(z)

Density of total dark energy (fluid + cosmological constant), unitless.

Omega_fld(z)

Density of cosmological constant, unitless.

Omega_g(z)

Density parameter of photons, unitless.

Omega_k(z)

Density parameter of curvature, unitless.

Omega_m(z)

Density parameter of non-relativistic (matter-like) component, including non-relativistic part of massive neutrino, unitless.

Omega_ncdm(z, species=None)

Density parameter of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_ncdm_tot(z)

Total density parameter of massive neutrinos, unitless.

Omega_pncdm(z, species=None)

Density parameter of pressure of non-relativistic part of massive neutrinos, unitless. If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

Omega_pncdm_tot(z)

Total density parameter of pressure of non-relativistic part of massive neutrinos, unitless.

Omega_r(z)

Density parameter of radiation, including photons and relativistic part of massive and massless neutrinos, unitless.

Omega_ur(z)

Density parameter of massless neutrinos, unitless.

T_cmb(z)

The CMB temperature, in \(K\).

T_ncdm(z, species=None)

Return the ncdm temperature (massive neutrinos), in \(K\). Returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)).

property age

The current age of the Universe, in \(\mathrm{Gy}\).

angular_diameter_distance(z)

Proper angular diameter distance, in \(\mathrm{Mpc}/h\).

See eq. 18 of astro-ph/9905116 for \(D_{A}(z)\).

angular_diameter_distance_2(z1, z2)

Angular diameter distance of object at \(z_{2}\) as seen by observer at \(z_{1}\), that is, \(S_{K}((\chi(z_{2}) - \chi(z_{1})) \sqrt{|K|}) / \sqrt{|K|} / (1 + z_{2})\), where \(S_{K}\) is the identity if \(K = 0\), \(\sin\) if \(K < 0\) and \(\sinh\) if \(K > 0\). camb’s angular_diameter_distance2(z1, z2) is not used as it returns 0 when z2 < z1.

comoving_angular_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

comoving_radial_distance(z)

Comoving radial distance, in \(mathrm{Mpc}/h\).

See eq. 15 of astro-ph/9905116 for \(D_C(z)\).

comoving_transverse_distance(z)

Comoving transverse distance, in \(\mathrm{Mpc}/h\).

See eq. 16 of astro-ph/9905116 for \(D_{M}(z)\).

efunc(z)

Function giving \(E(z)\), where the Hubble parameter is defined as \(H(z) = H_{0} E(z)\), unitless.

growth_factor(z, znorm=None)

Approximation of growth factor.

Parameters:
  • z (array_like) – Redshifts.

  • znorm (float, default=None) – If provided, growth factor is normalized as (1 + znorm) / (1 + z) in the matter domination era. Else, growth_factor is normalized to 1 at z = 0.

References

https://arxiv.org/abs/astro-ph/9709112, eq. 4 https://ui.adsabs.harvard.edu/abs/1992ARA%26A..30..499C/abstract, eq. 29

growth_rate(z)

Approximation of growth rate.

References

https://arxiv.org/abs/astro-ph/0507263

hubble_function(z)

Hubble function ba.index_bg_H, in \(\mathrm{km}/\mathrm{s}/\mathrm{Mpc}\).

luminosity_distance(z)

Luminosity distance, in \(\mathrm{Mpc}/h\).

See eq. 21 of astro-ph/9905116 for \(D_{L}(z)\).

p_ncdm(z, species=None)

Pressure of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return pressure for this species.

p_ncdm_tot(z)

Total pressure of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_Lambda(z)

Comoving density of cosmological constant \(\rho_{\Lambda}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_b(z)

Comoving density of baryons \(\rho_{b}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_cdm(z)

Comoving density of cold dark matter \(\rho_{cdm}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_crit(z)

Comoving critical density excluding curvature \(\rho_{c}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

This is defined as:

\[\rho_{\mathrm{crit}}(z) = \frac{3 H(z)^{2}}{8 \pi G}.\]
rho_de(z)

Total comoving density of dark energy \(\rho_{\mathrm{de}}\) (fluid + cosmological constant), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_fld(z)

Comoving density of dark energy fluid \(\rho_{\mathrm{fld}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_g(z)

Comoving density of photons \(\rho_{g}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_k(z)

Comoving density of curvature \(\rho_{k}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_m(z)

Comoving density of matter \(\rho_{m}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ncdm(z, species=None)

Comoving density of non-relativistic part of massive neutrinos for each species, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\). If species is None returned shape is (N_ncdm,) if z is a scalar, else (N_ncdm, len(z)). Else if species is between 0 and N_ncdm, return density for this species.

rho_ncdm_tot(z)

Total comoving density of non-relativistic part of massive neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_r(z)

Comoving density of radiation \(\rho_{r}\), including photons and relativistic part of massive and massless neutrinos, in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_tot(z)

Comoving total density \(\rho_{\mathrm{tot}}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

rho_ur(z)

Comoving density of massless neutrinos \(\rho_{ur}\), in \(10^{10} M_{\odot}/h / (\mathrm{Mpc}/h)^{3}\).

time(z)

Proper time (age of universe), in \(\mathrm{Gy}\).

class cosmoprimo.eisenstein_hu.EisensteinHuEngine(*args, **kwargs)

Bases: BaseEngine

Implementation of Eisenstein & Hu analytic formulae.

References

https://arxiv.org/abs/astro-ph/9709112

Initialize engine.

Parameters:
  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Engine parameters.

compute()

Precompute coefficients for the transfer function.

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background()

Return Background calculations.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier()

Return Fourier calculations.

get_harmonic()

Return Harmonic calculations.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations()

Return Perturbations calculations.

get_primordial()

Return Primordial calculations.

get_thermodynamics()

Return Thermodynamics calculations.

get_transfer()

Return Transfer calculations.

class cosmoprimo.eisenstein_hu.Fourier(engine)

Bases: BaseSection

pk_interpolator(of='delta_m', **kwargs)

Return PowerSpectrumInterpolator2D instance.

Parameters:
  • of (string, tuple) – Perturbed quantities: ‘delta_m’, ‘theta_m’. No distinction is made between baryons and CDM. Requesting velocity divergence ‘theta_xx’ will rescale the power spectrum by the growth rate as a function of z.

  • kwargs (dict) – Arguments for PowerSpectrumInterpolator2D.

property sigma8_m

Current r.m.s. of matter perturbations in a sphere of \(8 \mathrm{Mpc}/h\), unitless.

sigma8_z(z, of='delta_m')

Return the r.m.s. of of perturbations in sphere of \(8 \mathrm{Mpc}/h\). No distinction is made between baryons and CDM.

sigma_rz(r, z, of='delta_m', **kwargs)

Return the r.m.s. of of perturbations in sphere of \(r \mathrm{Mpc}/h\). No distinction is made between baryons and CDM.

class cosmoprimo.eisenstein_hu.Primordial(engine)

Bases: BaseSection

Initialize Primordial.

property A_s

Scalar amplitude of the primordial power spectrum at \(k_\mathrm{pivot}\), unitless.

property ln_1e10_A_s

\(\ln(10^{10}A_s)\), unitless.

pk_interpolator(mode='scalar')

Return power spectrum interpolator.

Parameters:

mode (string, default='scalar') – ‘scalar’, ‘vector’ or ‘tensor’ mode.

Returns:

interp

Return type:

PowerSpectrumInterpolator1D

pk_k(k, mode='scalar')

The primordial spectrum of curvature perturbations at k, generated by inflation, in \((\mathrm{Mpc}/h)^{3}\). For scalar perturbations this is e.g. defined as:

\[\mathcal{P_R}(k) = A_s \left (\frac{k}{k_\mathrm{pivot}} \right )^{n_s - 1 + 1/2 \alpha_s \ln(k/k_\mathrm{pivot}) + 1/6 \beta_s \ln(k/k_\mathrm{pivot})^2}\]

See also: eq. 2 of this reference.

Parameters:
  • k (array_like) – Wavenumbers, in \(h/\mathrm{Mpc}\).

  • mode (string, default='scalar') – ‘scalar’ mode.

Returns:

pk – The primordial power spectrum.

Return type:

array

class cosmoprimo.eisenstein_hu.Thermodynamics(engine)

Bases: BaseSection

Initialize Thermodynamics.

class cosmoprimo.eisenstein_hu.Transfer(engine)

Bases: BaseSection

transfer_k(k)

Return matter transfer function.

Parameters:

k (array_like) – Wavenumbers.

Returns:

transfer

Return type:

array

Eisenstein & Hu no-wiggle

class cosmoprimo.eisenstein_hu_nowiggle.EisensteinHuNoWiggleEngine(*args, **kwargs)

Bases: EisensteinHuEngine

Implementation of Eisenstein & Hu no-wiggle analytic formulae.

References

https://arxiv.org/abs/astro-ph/9709112

Initialize engine.

Parameters:
  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Engine parameters.

compute()

Precompute coefficients for the transfer function.

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background()

Return Background calculations.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier()

Return Fourier calculations.

get_harmonic()

Return Harmonic calculations.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations()

Return Perturbations calculations.

get_primordial()

Return Primordial calculations.

get_thermodynamics()

Return Thermodynamics calculations.

get_transfer()

Return Transfer calculations.

class cosmoprimo.eisenstein_hu_nowiggle.Transfer(engine)

Bases: BaseSection

transfer_k(k)

Return matter transfer function.

Parameters:

k (array_like) – Wavenumbers.

Returns:

transfer

Return type:

array

BBKS

class cosmoprimo.bbks.BBKSEngine(*args, **kwargs)

Bases: BaseEngine

Implementation of BBKS no-wiggle analytic formulae.

References

https://ui.adsabs.harvard.edu/abs/1986ApJ…304…15B/abstract https://arxiv.org/abs/astro-ph/9412025 https://arxiv.org/abs/1812.05995

Initialize engine.

Parameters:
  • extra_params (dict, default=None) – Extra engine parameters, typically precision parameters.

  • params (dict) – Engine parameters.

compute()

Precompute coefficients for the transfer function.

get(*args, **kwargs)

Return an input (or easily derived) parameter.

get_background()

Return Background calculations.

classmethod get_default_params(of=None, include_conflicts=True)

Return default input parameters.

Parameters:
  • of (string, default=None) – One of [‘cosmology’, ‘calculation’]. If None, returns all parameters.

  • include_conflicts (bool, default=True) – Whether to include conflicting parameters (then all accepted parameters).

Returns:

params – Dictionary of default parameters.

Return type:

dict

get_fourier()

Return Fourier calculations.

get_harmonic()

Return Harmonic calculations.

get_params(of='base')

Return parameters.

Parameters:

of (string, default='base') – One of [‘cosmology’, ‘calculation’, ‘base’, ‘derived’, ‘extra’, ‘all’]. If all, returns base, derived and extra parameters.

Returns:

params – Dictionary of parameters.

Return type:

dict

get_perturbations()

Return Perturbations calculations.

get_primordial()

Return Primordial calculations.

get_thermodynamics()

Return Thermodynamics calculations.

get_transfer()

Return Transfer calculations.

class cosmoprimo.bbks.Transfer(engine)

Bases: BaseSection

transfer_k(k)

Return matter transfer function.

Parameters:

k (array_like) – Wavenumbers.

Returns:

transfer

Return type:

array

Utilities

BAO filter

Different techniques to extract the BAO signal from the power spectrum of correlation function. For the power spectrum, the most “good-looking” ones are Wallish2018PowerSpectrumBAOFilter, :class:`Brieden2022PowerSpectrumBAOFilter, PeakAveragePowerSpectrumBAOFilter. For the correlation function: Kirkby2013CorrelationFunctionBAOFilter. jax-differentiable are: - Hinton2017PowerSpectrumBAOFilter, only after initialization (to set peak maxima; should be fixable). - EHNoWigglePolyPowerSpectrumBAOFilter - PeakAveragePowerSpectrumBAOFilter, only after initialization (to set peak maxima; should be fixable). - Kirkby2013CorrelationFunctionBAOFilter

class cosmoprimo.bao_filter.BSplinePowerSpectrumBAOFilter(pk_interpolator, constraint=('sigma8',), cosmo=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

Filter BAO wiggles with B-splines.

References

https://arxiv.org/pdf/1509.02120.pdf, Appendix A (thanks to Stephen Chen for the reference and code)

Run BAO filter.

Parameters:
  • pk_interpolator (PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D) – Input power spectrum to remove BAO wiggles from.

  • constraint (str, list, tuple, default=('sigma8',)) – Quantities computed on the no-wiggle power spectrum required to match the input pk_interpolator.

  • cosmo (Cosmology, default=None) – Cosmology instance, used to compute the Eisenstein & Hu no-wiggle power spectrum.

  • kwargs (dict) – Arguments for set_k().

property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

class cosmoprimo.bao_filter.BaseCorrelationFunctionBAOFilter(xi_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BaseClass

Base BAO filter for correlation function.

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_s(ns=1024)

Set separations where to evaluate the correlation function (between xi_interpolator.extrap_smin and xi_interpolator.extrap_smax).

Parameters:

ns (int, default=1024) – Number of separations.

set_xi(xi_interpolator, cosmo=None)

Set input correlation function to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of return power spectrum interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function.

Parameters:

kwargs (dict) – Override interpolation settings of xi_interpolator.

Returns:

interp – 1D or 2D depending on xi_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

class cosmoprimo.bao_filter.BasePowerSpectrumBAOFilter(pk_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BaseClass

Base BAO filter for power spectrum.

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

class cosmoprimo.bao_filter.Brieden2022PowerSpectrumBAOFilter(pk_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

Filter BAO wiggles by averaging the minima and maxima of the wiggles.

References

https://arxiv.org/abs/2204.11868, Appendix D (thanks to Samuel Brieden for the reference)

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

cosmoprimo.bao_filter.CorrelationFunctionBAOFilter(xi_interpolator, engine='kirkby2013', **kwargs)

Run correlation function BAO filter corresponding to the provided engine.

class cosmoprimo.bao_filter.EHNoWigglePolyPowerSpectrumBAOFilter(pk_interpolator, krange=(0.001, 1.0), rescale_krange=True, cosmo=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

Remove BAO wiggles using the Eisenstein & Hu no-wiggle analytic formula, emulated with a 6-th order polynomial.

Run BAO filter.

Parameters:
  • pk_interpolator (PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D) – Input power spectrum to remove BAO wiggles from.

  • krange (tuple, default=(1e-3, 1.)) – k-range to fit the Eisenstein & Hu no-wiggle power spectrum to the input one pk_interpolator.

  • rescale_krange (bool, default=True) – Whether to rescale krange by the ratio of rs_drag relative to the fiducial cosmology (may help robustify the procedure for cosmologies far from the fiducial one).

  • cosmo (Cosmology, default=None) – Cosmology instance, used to compute the Eisenstein & Hu no-wiggle power spectrum.

  • kwargs (dict) – Arguments for set_k().

property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

class cosmoprimo.bao_filter.EHNoWiggleSavGolPowerSpectrumBAOFilter(pk_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

BAO smoothing with Savitzky-Golay filter applied on Eisenstein & Hu no-wiggle analytic formula.

References

Stephen Chen.

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

class cosmoprimo.bao_filter.Hinton2017PowerSpectrumBAOFilter(pk_interpolator, degree=12, sigma=0.5, weight=0.9, **kwargs)

Bases: BasePowerSpectrumBAOFilter

Power spectrum BAO filter consisting in fitting a high degree polynomial to the input power spectrum in log-log space.

References

https://github.com/Samreay/Barry/blob/master/barry/cosmology/power_spectrum_smoothing.py

Note

We have hand-tune parameters w.r.t. the reference.

Run BAO filter.

Parameters:
  • pk_interpolator (PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D) – Input power spectrum to remove BAO wiggles from.

  • degree (int, default=12) – Polynomial degree.

  • sigma (float, default=0.5) – Standard deviation of the Gaussian kernel that downweights the maximum of the power spectrum relative to the edges.

  • weight (float, default=0.9) – Normalisation of the Gaussian kernel.

property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

class cosmoprimo.bao_filter.Kirkby2013CorrelationFunctionBAOFilter(xi_interpolator, srange_left=(50.0, 82.0), srange_right=(150.0, 190.0), rescale_sbox=True, cosmo=None, **kwargs)

Bases: BaseCorrelationFunctionBAOFilter

Filter BAO peak by cutting the peak and interpolating with 5-order polynomial.

References

https://arxiv.org/abs/1301.3456 https://github.com/igmhub/picca/blob/master/bin/picca_compute_pk_pksb.py

Run BAO filter.

Parameters:
  • xi_interpolator (CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D) – Input correlation function to remove BAO peak from.

  • srange_left (tuple) – s-range to fit the polynomial on the left-hand side of the BAO peak.

  • srange_right (tuple) – s-range to fit the polynomial on the right-hand side of the BAO peak.

  • cosmo (Cosmology) – Cosmology instance, which may be used to tune filter settings (depending on rs_drag).

  • rescale_sbox (bool) – Whether to rescale srange_left and srange_right by the ratio of rs_drag relative to the fiducial cosmology (may help robustify the procedure for cosmologies far from the fiducial one).

  • cosmo – Cosmology instance, used to compute the Eisenstein & Hu no-wiggle power spectrum.

  • kwargs (dict) – Arguments for set_s().

property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_s(ns=1024)

Set separations where to evaluate the correlation function (between xi_interpolator.extrap_smin and xi_interpolator.extrap_smax).

Parameters:

ns (int, default=1024) – Number of separations.

set_xi(xi_interpolator, cosmo=None)

Set input correlation function to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of return power spectrum interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function.

Parameters:

kwargs (dict) – Override interpolation settings of xi_interpolator.

Returns:

interp – 1D or 2D depending on xi_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

class cosmoprimo.bao_filter.PeakAveragePowerSpectrumBAOFilter(pk_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

Filter BAO wiggles by averaging the minima and maxima of the wiggles at the fiducial positions rescaled by \(r_{\mathrm{drag}} / r_{\mathrm{drag}}^{\mathrm{fid}}\). A simpler version of Brieden2022PowerSpectrumBAOFilter.

References

https://arxiv.org/abs/2204.11868, Appendix D (thanks to Samuel Brieden for the reference)

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

cosmoprimo.bao_filter.PowerSpectrumBAOFilter(pk_interpolator, engine='wallish2018', **kwargs)

Run power spectrum BAO filter corresponding to the provided engine.

class cosmoprimo.bao_filter.RegisteredCorrelationFunctionBAOFilter(name, bases, class_dict)

Bases: type

Metaclass registering BaseCorrelationFunctionBAOFilter-derived classes.

mro()

Return a type’s method resolution order.

class cosmoprimo.bao_filter.RegisteredPowerSpectrumBAOFilter(name, bases, class_dict)

Bases: type

Metaclass registering BasePowerSpectrumBAOFilter-derived classes.

mro()

Return a type’s method resolution order.

class cosmoprimo.bao_filter.SavGolPowerSpectrumBAOFilter(pk_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

BAO smoothing with Savitzky-Golay filter.

References

Taken from https://github.com/sfschen/velocileptors/blob/master/velocileptors/EPT/ept_fullresum_fftw.py

Note

Contrary to the reference, we work in \(\log(k)\) - \(\log(k P(k))\) space.

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

class cosmoprimo.bao_filter.Wallish2018PowerSpectrumBAOFilter(pk_interpolator, cosmo=None, cosmo_fid=None, **kwargs)

Bases: BasePowerSpectrumBAOFilter

Filter BAO wiggles by sine-transforming the power spectrum to real space (where the BAO is better localized), cutting the peak and interpolating with a spline.

References

https://arxiv.org/pdf/1810.02800.pdf, Appendix D (thanks to Stephen Chen for the reference) https://arxiv.org/pdf/1003.3999.pdf

Note

We have hand-tuned parameters w.r.t. the reference.

Run BAO filter.

Parameters:
property cosmo

Cosmology.

property cosmo_fid

Reference cosmology.

rs_drag_ratio()

If cosmo is provided, return the ratio of its rs_drag to the fiducial one (from Cosmology()), else 1.

set_k(nk=1024)

Set wavenumbers where to evaluate the power spectrum (between pk_interpolator.extrap_kmin and pk_interpolator.extrap_kmax).

Parameters:

nk (int, default=1024) – Number of wavenumbers.

set_pk(pk_interpolator, cosmo=None)

Set input power spectrum to remove BAO wiggles from.

smooth_pk_interpolator(**kwargs)

Return smooth (i.e. no-wiggle) power spectrum.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of pk_interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D

smooth_xi_interpolator(**kwargs)

Return smooth (i.e. no-peak) correlation function using FFTlog.

Parameters:

kwargs (dict) – Override interpolation and extrapolation settings of returned correlation function interpolator.

Returns:

interp – 1D or 2D depending on pk_interpolator.

Return type:

CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D

property wiggles

Extracted wiggles.

FFTlog

Implementation of the FFTlog algorithm, very much inspired by mcfit (https://github.com/eelregit/mcfit) and implementation in https://github.com/sfschen/velocileptors/blob/master/velocileptors/Utils/spherical_bessel_transform_fftw.py

class cosmoprimo.fftlog.BaseBesselKernel(nu)

Bases: BaseKernel

Base Bessel kernel.

class cosmoprimo.fftlog.BaseFFTEngine(size, nparallel=1, nthreads=None)

Bases: object

Base FFT engine.

Initialize FFT engine.

Parameters:
  • size (int) – Array size.

  • nparallel (int, default=1) – Number of FFTs to be performed in parallel.

  • nthreads (int, default=None) – Number of threads.

class cosmoprimo.fftlog.BaseKernel

Bases: object

Base kernel.

class cosmoprimo.fftlog.BaseTophatKernel(ndim=1)

Bases: BaseKernel

Base tophat kernel.

class cosmoprimo.fftlog.BesselJKernel(nu)

Bases: BaseBesselKernel

(Mellin transform of) Bessel kernel.

class cosmoprimo.fftlog.CorrelationToPower(s, ell=0, q=0, complex=False, **kwargs)

Bases: FFTlog

Correlation function to power spectrum transform, defined as:

\[P_{\ell}(k) = 4 \pi i^{\ell} \int ds s^{2} \xi_{\ell}(s) j_{\ell}(ks)\]

Initialize power to correlation transform.

Parameters:
  • s (array_like) – Input log-spaced separations. If 1D, is broadcast to the number of provided ell.

  • ell (int, list of int, default=0) – Poles. If a list is provided, will perform all transforms at once.

  • q (float, list of floats, default=0) – Power-law tilt(s) to regularise integration.

  • complex (bool, default=False) – False returns the real part of even poles, and the imaginary part of odd poles.

  • kwargs (dict) – Arguments for FFTlog.

inv()

Inverse the transform.

property nparallel

Number of transforms performed in parallel.

set_fft_engine(engine='numpy', **engine_kwargs)

Set up FFT engine. See get_fft_engine()

Parameters:
  • engine (BaseEngine, string, default='numpy') – FFT engine, or one of [‘numpy’, ‘fftw’].

  • engine_kwargs (dict) – Arguments for FFT engine.

property size

Size of x-coordinates.

class cosmoprimo.fftlog.FFTWEngine(size, nparallel=1, nthreads=None, wisdom=None, plan='measure')

Bases: BaseFFTEngine

FFT engine based on pyfftw.

Initialize pyfftw engine.

Parameters:
  • size (int) – Array size.

  • nparallel (int, default=1) – Number of FFTs to be performed in parallel.

  • nthreads (int, default=None) – Number of threads.

  • wisdom (string, tuple, default=None) – pyfftw wisdom, to speed up initialization of FFTs. If a string, should be a path to the save FFT wisdom (with numpy.save()). If a tuple, directly corresponds to the wisdom.

  • plan (string, default='measure') – Choices are [‘estimate’, ‘measure’, ‘patient’, ‘exhaustive’]. The increasing amount of effort spent during the planning stage to create the fastest possible transform. Usually ‘measure’ is a good compromise.

backward(fun)

Backward transform of fun.

forward(fun)

Forward transform of fun.

class cosmoprimo.fftlog.FFTlog(x, kernel, q=0, minfolds=2, lowring=True, xy=1, check_level=0, engine='numpy', **engine_kwargs)

Bases: object

Implementation of the FFTlog algorithm presented in https://jila.colorado.edu/~ajsh/FFTLog/, which computes the generic integral:

\[G(y) = \int_{0}^{\infty} x dx F(x) K(xy)\]

with \(F(x)\) input function, \(K(xy)\) a kernel.

This transform is (mathematically) invariant under a power law transformation:

\[G_{q}(y) = \int_{0}^{\infty} x dx F_{q}(x) K_{q}(xy)\]

where \(F_{q}(x) = G(x)x^{-q}\), \(K_{q}(t) = K(t)t^{q}\) and \(G_{q}(y) = G(y)y^{q}\).

Initialize FFTlog, which can perform several transforms at once.

Parameters:
  • x (array_like) – Input log-spaced coordinates. Must be strictly increasing. If 1D, is broadcast to the number of provided kernels.

  • kernel (callable, list of callables) – Mellin transform of the kernel: .. math:: U_{K}(z) = int_{0}^{infty} t^{z-1} K(t) dt If a list of kernels is provided, will perform all transforms at once.

  • q (float, list of floats, default=0) – Power-law tilt(s) to regularise integration.

  • minfolds (int, default=2) – Padded size is 2**n, with minimum \(n\) statisfying 2**n > minfolds * x.size.

  • lowring (bool, default=True) – If True set output coordinates according to the low-ringing condition, otherwise set it with xy.

  • xy (float, list of floats, default=1) – Enforce the reciprocal product (i.e. x[0] * y[-1]) of the input x and output y coordinates.

  • check_level (int, default=0) – If non-zero run sanity checks on input.

  • engine (string, default='numpy') – FFT engine. See set_fft_engine().

  • engine_kwargs (dict) – Arguments for FFT engine.

Note

Kernel definition is different from that of https://jila.colorado.edu/~ajsh/FFTLog/, which uses (eq. 10):

\[U_{K}(z) = \int_{0}^{\infty} t^{z} K(t) dt\]

Therefore, one should use \(q = 1\) for Bessel functions to match \(q = 0\) in https://jila.colorado.edu/~ajsh/FFTLog/.

inv()

Inverse the transform.

property nparallel

Number of transforms performed in parallel.

set_fft_engine(engine='numpy', **engine_kwargs)

Set up FFT engine. See get_fft_engine()

Parameters:
  • engine (BaseEngine, string, default='numpy') – FFT engine, or one of [‘numpy’, ‘fftw’].

  • engine_kwargs (dict) – Arguments for FFT engine.

property size

Size of x-coordinates.

class cosmoprimo.fftlog.GaussianKernel

Bases: BaseKernel

(Mellin transform of) Gaussian kernel.

class cosmoprimo.fftlog.GaussianSqKernel

Bases: BaseKernel

(Mellin transform of) square of Gaussian kernel.

class cosmoprimo.fftlog.GaussianVariance(k, q=0, **kwargs)

Bases: FFTlog

Variance in Gaussian window.

It relies on Gaussian kernel.

Initialize Gaussian variance transform.

Parameters:
  • k (array_like) – Input log-spaced wavenumbers. If 1D, is broadcast to the number of provided ell.

  • q (float, list of floats, default=0) – Power-law tilt(s) to regularise integration.

  • kwargs (dict) – Arguments for FFTlog.

inv()

Inverse the transform.

property nparallel

Number of transforms performed in parallel.

set_fft_engine(engine='numpy', **engine_kwargs)

Set up FFT engine. See get_fft_engine()

Parameters:
  • engine (BaseEngine, string, default='numpy') – FFT engine, or one of [‘numpy’, ‘fftw’].

  • engine_kwargs (dict) – Arguments for FFT engine.

property size

Size of x-coordinates.

class cosmoprimo.fftlog.HankelTransform(x, nu=0, **kwargs)

Bases: FFTlog

Hankel transform implementation using FFTlog.

It relies on Bessel function kernels.

Initialize Hankel transform.

Parameters:
  • x (array_like) – Input log-spaced coordinates. If 1D, is broadcast to the number of provided nu.

  • nu (int, list of int, default=0) – Order of Bessel functions. If a list is provided, will perform all transforms at once.

  • kwargs (dict) – Arguments for FFTlog.

inv()

Inverse the transform.

property nparallel

Number of transforms performed in parallel.

set_fft_engine(engine='numpy', **engine_kwargs)

Set up FFT engine. See get_fft_engine()

Parameters:
  • engine (BaseEngine, string, default='numpy') – FFT engine, or one of [‘numpy’, ‘fftw’].

  • engine_kwargs (dict) – Arguments for FFT engine.

property size

Size of x-coordinates.

class cosmoprimo.fftlog.NumpyFFTEngine(size, nparallel=1, nthreads=None)

Bases: BaseFFTEngine

FFT engine based on numpy.fft.

Initialize FFT engine.

Parameters:
  • size (int) – Array size.

  • nparallel (int, default=1) – Number of FFTs to be performed in parallel.

  • nthreads (int, default=None) – Number of threads.

backward(fun)

Backward transform of fun.

forward(fun)

Forward transform of fun.

class cosmoprimo.fftlog.PowerToCorrelation(k, ell=0, q=0, complex=False, **kwargs)

Bases: FFTlog

Power spectrum to correlation function transform, defined as:

\[\xi_{\ell}(s) = \frac{(-i)^{\ell}}{2 \pi^{2}} \int dk k^{2} P_{\ell}(k) j_{\ell}(ks)\]

Initialize power to correlation transform.

Parameters:
  • k (array_like) – Input log-spaced wavenumbers. If 1D, is broadcast to the number of provided ell.

  • ell (int, list of int, default=0) – Poles. If a list is provided, will perform all transforms at once.

  • q (float, list of floats, default=0) – Power-law tilt(s) to regularise integration.

  • complex (bool, default=False) – False assumes the imaginary part of odd power spectrum poles is provided.

  • kwargs (dict) – Arguments for FFTlog.

inv()

Inverse the transform.

property nparallel

Number of transforms performed in parallel.

set_fft_engine(engine='numpy', **engine_kwargs)

Set up FFT engine. See get_fft_engine()

Parameters:
  • engine (BaseEngine, string, default='numpy') – FFT engine, or one of [‘numpy’, ‘fftw’].

  • engine_kwargs (dict) – Arguments for FFT engine.

property size

Size of x-coordinates.

class cosmoprimo.fftlog.SphericalBesselJKernel(nu)

Bases: BaseBesselKernel

(Mellin transform of) spherical Bessel kernel.

class cosmoprimo.fftlog.TophatKernel(ndim=1)

Bases: BaseTophatKernel

(Mellin transform of) tophat kernel.

class cosmoprimo.fftlog.TophatSqKernel(ndim=1)

Bases: BaseTophatKernel

(Mellin transform of) square of tophat kernel.

class cosmoprimo.fftlog.TophatVariance(k, q=0, **kwargs)

Bases: FFTlog

Variance in tophat window.

It relies on tophat kernel.

Initialize tophat variance transform.

Parameters:
  • k (array_like) – Input log-spaced wavenumbers. If 1D, is broadcast to the number of provided ell.

  • q (float, list of floats, default=0) – Power-law tilt(s) to regularise integration.

  • kwargs (dict) – Arguments for FFTlog.

inv()

Inverse the transform.

property nparallel

Number of transforms performed in parallel.

set_fft_engine(engine='numpy', **engine_kwargs)

Set up FFT engine. See get_fft_engine()

Parameters:
  • engine (BaseEngine, string, default='numpy') – FFT engine, or one of [‘numpy’, ‘fftw’].

  • engine_kwargs (dict) – Arguments for FFT engine.

property size

Size of x-coordinates.

cosmoprimo.fftlog.apply_along_last_axes(func, array, naxes=1, toret=None)

Apply callable func over the last naxes of array.

cosmoprimo.fftlog.get_fft_engine(engine, *args, **kwargs)

Return FFT engine.

Parameters:
  • engine (BaseFFTEngine, string) – FFT engine, or one of [‘numpy’, ‘fftw’].

  • args (tuple, dict) – Arguments for FFT engine.

  • kwargs (tuple, dict) – Arguments for FFT engine.

Returns:

engine

Return type:

BaseFFTEngine

cosmoprimo.fftlog.pad(array, pad_width, axis=-1, extrap=0)

Pad array along axis.

Parameters:
  • array (array_like) – Input array to be padded.

  • pad_width (int, tuple of ints) – Number of points to be added on both sides of the array. Pass a tuple to differentiate between left and right sides.

  • axis (int, default=-1) – Axis along which padding is to be applied.

  • extrap (string, float, default=0) – If ‘log’, performs a log-log extrapolation. If ‘edge’, pad array with its edge values. Else, pad array with the provided value. Pass a tuple to differentiate between left and right sides.

Returns:

array – Padded array.

Return type:

array

Interpolators

Utilities for power spectrum and correlation function interpolations. Useful classes are PowerSpectrumInterpolator1D, PowerSpectrumInterpolator2D, CorrelationFunctionInterpolator1D, CorrelationFunctionInterpolator2D.

class cosmoprimo.interpolator.CorrelationFunctionInterpolator1D(s, xi, interp_s='log', interp_order_s=3)

Bases: _BaseCorrelationFunctionInterpolator

1D correlation funcion interpolator.

Initialize CorrelationFunctionInterpolator1D.

Parameters:
  • s (array_like) – Seperations.

  • xi (array_like) – Correlation function to be interpolated.

  • interp_s (string, default='log') – If ‘log’, interpolation is performed in log-s coordinates.

  • interp_order_s (int, default=3) – Interpolation order, i.e. degree of smoothing spline along s.

as_dict()

Return interpolator as a dictionary.

clone(**kwargs)

Clone interpolator, i.e. return a deepcopy with (possibly) other attributes in kwargs (clone() witout arguments is the same as deepcopy()).

See deepcopy() doc for warning about interpolators built from callables.

deepcopy()

Deep copy interpolator.

If interpolator interp1 is built from callable, requires its evaluation at k, such that e.g.: >>> interp2 = interp1.clone() will not be provide exactly the same interpolated values as interp1 (due to interpolation errors).

property extrap_smax

Maximum (extrapolated) s value (same as maximum interpolated value).

property extrap_smin

Minimum (extrapolated) s value (same as minimum interpolated value).

classmethod from_callable(s=None, xi_callable=None)

Build CorrelationFunctionInterpolator1D from callable.

Parameters:
  • s (array_like, default=None) – Array of separations where the provided xi_callable can be trusted. It will be used if :attr:xi is requested. Must be strictly increasing.

  • xi_callable (callable) – Correlation function callable.

Returns:

self

Return type:

CorrelationFunctionInterpolator1D

params()

Return interpolator parameter dictionary.

rescale_sigma8(sigma8=1.0)

Rescale the correlation function to the provided sigma8 normalisation.

sigma8(**kwargs)

Return the r.m.s. of perturbations in a sphere of 8.

sigma_d(**kwargs)

Return the r.m.s. of the displacement field by transforming correlation function into power spectrum.

See PowerSpectrumInterpolator1D.sigma_d() arguments.

sigma_r(r, **kwargs)

Return the r.m.s. of perturbations in a sphere of \(r\) by transforming correlation function into power spectrum.

See PowerSpectrumInterpolator1D.sigma_r() arguments.

property smax

Maximum (interpolated) k value.

property smin

Minimum (interpolated) s value.

to_pk(ns=1024, fftlog_kwargs=None, **kwargs)

Transform correlation function into power spectrum using FFTlog.

nsint, default=1024

Number of separations used in FFTlog transform.

fftlog_kwargsdict, default=None

Arguments for FFTlog.

kwargsdict

Arguments for the new PowerSpectrumInterpolator1D instance.

Returns:

pk

Return type:

PowerSpectrumInterpolator1D

property xi

Return correlation function array (if interpolator built from callable, evaluate it), with normalisation.

class cosmoprimo.interpolator.CorrelationFunctionInterpolator2D(s, z, xi=None, interp_s='log', interp_order_s=3, interp_order_z=None, growth_factor_sq=None)

Bases: _BaseCorrelationFunctionInterpolator

2D correlation function interpolator.

Initialize CorrelationFunctionInterpolator2D.

growth_factor_sq is a callable that can be prodided to rescale the output of the base spline interpolation. Indeed, variations of \(z \rightarrow \xi(k,z)\) are (mostly) \(s\) scale independent, such that more accurate interpolation in z can be achieved by providing the z variations separately in a well-sampled growth_factor_sq.

Parameters:
  • s (array_like) – Separations.

  • z (array_like, float, default=0) – Redshifts.

  • xi (array_like) – Correlation function to be interpolated. If z is scalar, should be 1D; else 2D, with shape (s.size, z.size).

  • interp_s (string, default='log') – If ‘log’, interpolation is performed in log-s coordinates.

  • interp_order_s (int, default=3) – Interpolation order, i.e. degree of smoothing spline along s.

  • interp_order_z (int, default=None) – Interpolation order, i.e. degree of smoothing spline along z. If None, the maximum order given z size (see GenericSpline.min_spline_order()) is considered.

  • growth_factor_sq (callable, default=None) – Function that takes z as argument and returns the growth factor squared at that redshift. This will rescale the output of the base spline interpolation. Therefore, make sure that provided pk does not contain the redundant z variations.

as_dict()

Return interpolator as a dictionary.

clone(**kwargs)

Clone interpolator, i.e. return a deepcopy with (possibly) other attributes in kwargs (clone() witout arguments is the same as deepcopy()).

See deepcopy() doc for warning about interpolators built from callables.

deepcopy()

Deep copy interpolator.

If interpolator interp1 is built from callable, requires its evaluation at k, such that e.g.: >>> interp2 = interp1.clone() will not be provide exactly the same interpolated values as interp1 (due to interpolation errors).

property extrap_smax

Maximum (extrapolated) s value (same as maximum interpolated value).

property extrap_smin

Minimum (extrapolated) s value (same as minimum interpolated value).

classmethod from_callable(s=None, z=None, xi_callable=None, growth_factor_sq=None)

Build CorrelationFunctionInterpolator2D from callable.

Parameters:
  • s (array_like, default=None) – Array of separations where the provided xi_callable can be trusted. It will be used if :attr:xi is requested. Must be strictly increasing.

  • z (array_like, default=None) – Array of redshifts where the provided xi_callable can be trusted. Same remark as for s.

  • xi_callable (callable, default=None) – Correlation function callable. If growth_factor_sq is not provided, should take s, z, grid as arguments (see __call__()) else, should take s as arguments.

  • growth_factor_sq (callable, default=None) – See remark above.

Returns:

self

Return type:

CorrelationFunctionInterpolator2D

growth_rate_rz(r, z, **kwargs)

Evaluate the growth rate at the log-derivative of perturbations in a sphere of \(r\) by transforming correlation function into power spectrum.

See PowerSpectrumInterpolator2D.growth_rate_rz() arguments.

params()

Return interpolator parameter dictionary.

rescale_sigma8(sigma8=1.0)

Rescale correlation function to the provided sigma8 normalisation at \(z = 0\).

sigma8_z(z, **kwargs)

Return the r.m.s. of perturbations in a sphere of 8 by transforming correlation function into power spectrum.

See PowerSpectrumInterpolator2D.sigma8_z() arguments.

sigma_dz(z, **kwargs)

Return the r.m.s. of the displacement field by transforming correlation function into power spectrum.

See PowerSpectrumInterpolator2D.sigma_dz() arguments.

sigma_rz(r, z, **kwargs)

Return the r.m.s. of perturbations in a sphere of \(r\) by transforming correlation function into power spectrum.

See PowerSpectrumInterpolator2D.sigma_rz() arguments.

property smax

Maximum (interpolated) k value.

property smin

Minimum (interpolated) s value.

to_1d(z, **kwargs)

Return CorrelationFunctionInterpolator1D instance corresponding to interpolator at z.

Parameters:
  • z (float, default=0) – Redshift.

  • kwargs (dict) – Arguments for the new CorrelationFunctionInterpolator1D instance. By default, the new instance inherits the same s interpolattion settings.

Returns:

interp

Return type:

CorrelationFunctionInterpolator1D

to_pk(ns=1024, fftlog_kwargs=None, **kwargs)

Transform correlation function into power spectrum using FFTlog.

nsint, default=1024

Number of separations used in FFTlog transform.

fftlog_kwargsdict, default=None

Arguments for FFTlog.

kwargsdict

Arguments for the new PowerSpectrumInterpolator2D instance.

Returns:

pk

Return type:

PowerSpectrumInterpolator2D

property xi

Return xi (if interpolator built from callable, evaluate it), without growth factor, but with normalisation.

property zmax

Maximum (spline-interpolated) redshift.

property zmin

Minimum (spline-interpolated) redshift.

class cosmoprimo.interpolator.PowerSpectrumInterpolator1D(k, pk, interp_k='log', extrap_pk='log', extrap_kmin=1e-07, extrap_kmax=100.0, interp_order_k=3)

Bases: _BasePowerSpectrumInterpolator

1D power spectrum interpolator, broadly adapted from CAMB’s P(k) interpolator by Antony Lewis in https://github.com/cmbant/CAMB/blob/master/camb/results.py, providing extra useful methods, such as sigma_r() or to_xi().

as_dict()

Return interpolator as a dictionary.

clone(**kwargs)

Clone interpolator, i.e. return a deepcopy with (possibly) other attributes in kwargs (clone() witout arguments is the same as deepcopy()).

See deepcopy() doc for warning about interpolators built from callables.

deepcopy()

Deep copy interpolator.

If interpolator interp1 is built from callable, requires its evaluation at k, such that e.g.: >>> interp2 = interp1.clone() will not be provide exactly the same interpolated values as interp1 (due to interpolation errors).

classmethod from_callable(k=None, pk_callable=None, extrap_kmin=1e-07, extrap_kmax=100.0)

Build PowerSpectrumInterpolator1D from callable.

Parameters:
  • k (array_like, default=None) – Array of wavenumbers where the provided pk_callable can be trusted. It will be used if pk is requested. Must be strictly increasing.

  • pk_callable (callable, default=None) – Power spectrum callable.

Returns:

new

Return type:

PowerSpectrumInterpolator1D

property kmax

Maximum (interpolated) k value.

property kmin

Minimum (interpolated) k value.

params()

Return interpolator parameter dictionary.

property pk

Return power spectrum array (if interpolator built from callable, evaluate it), with normalisation.

rescale_sigma8(sigma8=1.0)

Rescale power spectrum to the provided sigma8 normalisation.

sigma8(**kwargs)

Return the r.m.s. of perturbations in a sphere of 8.

sigma_d(**kwargs)

Return the r.m.s. of the displacement field, i.e.:

\[\sigma_{d} = \sqrt{\frac{1}{6 \pi^{2}} \int dk P(k)}\]
Parameters:
  • nk (int, default=1024) – If not None, performs trapezoidal integration with nk points between extrap_kmin and extrap_kmax. Else, uses scipy.integrate.quad.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration).

Returns:

sigmad

Return type:

array_like

sigma_r(r, **kwargs)

Return the r.m.s. of perturbations in a sphere of \(r\), i.e.:

\[\sigma_{r} = \sqrt{\frac{1}{2 \pi^{2}} \int dk k^{2} P(k) W^{2}(kr))\]
Parameters:
  • r (array_like) – Sphere radius.

  • nk (int, default=1024) – If not None, performs trapezoidal integration with nk points between extrap_kmin and extrap_kmax. Else, uses scipy.integrate.quad.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration).

Returns:

sigmar – Array of shape (r.size,).

Return type:

array_like

to_xi(nk=1024, fftlog_kwargs=None, **kwargs)

Transform power spectrum into correlation function using FFTlog.

nkint, default=1024

Number of wavenumbers used in FFTlog transform.

fftlog_kwargsdict, default=None

Arguments for FFTlog.

kwargsdict

Arguments for the new CorrelationFunctionInterpolator1D instance.

Returns:

xi

Return type:

CorrelationFunctionInterpolator1D

class cosmoprimo.interpolator.PowerSpectrumInterpolator2D(k, z, pk, interp_k='log', extrap_pk='log', extrap_kmin=1e-07, extrap_kmax=100.0, interp_order_k=3, interp_order_z=3, growth_factor_sq=None)

Bases: _BasePowerSpectrumInterpolator

2D power spectrum interpolator, broadly adapted from CAMB’s P(k) interpolator by Antony Lewis in https://github.com/cmbant/CAMB/blob/master/camb/results.py, providing extra useful methods, such as sigma_rz() or to_xi().

as_dict()

Return interpolator as a dictionary.

clone(**kwargs)

Clone interpolator, i.e. return a deepcopy with (possibly) other attributes in kwargs (clone() witout arguments is the same as deepcopy()).

See deepcopy() doc for warning about interpolators built from callables.

deepcopy()

Deep copy interpolator.

If interpolator interp1 is built from callable, requires its evaluation at k, such that e.g.: >>> interp2 = interp1.clone() will not be provide exactly the same interpolated values as interp1 (due to interpolation errors).

classmethod from_callable(k=None, z=None, pk_callable=None, growth_factor_sq=None, extrap_kmin=1e-07, extrap_kmax=100.0)

Build PowerSpectrumInterpolator2D from callable.

Parameters:
  • k (array_like, default=None) – Array of wavenumbers where the provided pk_callable can be trusted. It will be used if pk is requested. Must be strictly increasing.

  • z (array_like, default=None) – Array of redshifts where the provided pk_callable can be trusted. Same remark as for k.

  • pk_callable (callable, default=None) – Power spectrum callable. If growth_factor_sq is not provided, should take k, z, grid as arguments (see __call__()) else, should take k as arguments.

  • growth_factor_sq (callable, default=None) – Function that takes z as argument and returns the growth factor squared at that redshift. See remark above.

Returns:

self

Return type:

PowerSpectrumInterpolator2D

growth_rate_rz(r, z, dz=0.001, **kwargs)

Evaluate the growth rate at the log-derivative of perturbations in a sphere of \(r\), i.e.:

With \(z = ln(a)\).

Parameters:
  • r (array_like) – Sphere radii.

  • z (array_like, default=0) – Redshifts.

  • dz (float, default=1e-3) – z interval used for finite differentiation.

  • nk (int, default=1024) – If not None, performs trapezoidal integration with nk points between extrap_kmin and extrap_kmax. Else, uses scipy.integrate.quad.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration).

property kmax

Maximum (interpolated) k value.

property kmin

Minimum (interpolated) k value.

params()

Return interpolator parameter dictionary.

property pk

Return power spectrum array (if interpolator built from callable, evaluate it), without growth factor, but with normalisation.

rescale_sigma8(sigma8=1.0)

Rescale power spectrum to the provided sigma8 normalisation at \(z = 0\).

sigma8_z(z=0, **kwargs)

Return the r.m.s. of perturbations in a sphere of 8.

sigma_dz(z, **kwargs)

Return the r.m.s. of the displacement field, i.e.:

\[\sigma_{d}(z) = \sqrt{\frac{1}{6 \pi^{2}} \int dk P(k,z)}\]
Parameters:
  • z (array_like, default=0) – Redshifts.

  • nk (int, default=1024) – If not None, performs trapezoidal integration with nk points between extrap_kmin and extrap_kmax. Else, uses scipy.integrate.quad.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration).

Returns:

sigmadz

Return type:

array_like

sigma_rz(r, z, **kwargs)

Return the r.m.s. of perturbations in a sphere of \(r\), i.e.:

\[\sigma_{r}(z) = \sqrt{\frac{1}{2 \pi^{2}} \int dk k^{2} P(k, z) W^{2}(kr)}\]
Parameters:
  • r (array_like) – Sphere radii.

  • z (array_like, default=0) – Redshifts.

  • nk (int, default=1024) – If not None, performs trapezoidal integration with nk points between extrap_kmin and extrap_kmax. Else, uses scipy.integrate.quad.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration).

Returns:

sigmarz – Array of shape (r.size, z.size).

Return type:

array_like

to_1d(z, **kwargs)

Return PowerSpectrumInterpolator1D instance corresponding to interpolator at z.

Parameters:
  • z (float, array) – Redshift.

  • kwargs (dict) – Arguments for the new PowerSpectrumInterpolator1D instance. By default, the new instance inherits the same k interpolattion and extrapolation settings.

Returns:

interp

Return type:

PowerSpectrumInterpolator1D

to_xi(nk=1024, fftlog_kwargs=None, **kwargs)

Transform power spectrum into correlation function using FFTlog.

nkint, default=1024

Number of wavenumbers used in FFTlog transform.

fftlog_kwargsdict, default=None

Arguments for FFTlog.

kwargsdict

Arguments for the new CorrelationFunctionInterpolator2D instance.

Returns:

xi

Return type:

CorrelationFunctionInterpolator2D

property zmax

Maximum (spline-interpolated) redshift.

property zmin

Minimum (spline-interpolated) redshift.

cosmoprimo.interpolator.integrate_sigma_d2(pk, kmin=1e-07, kmax=100.0, method='simpson', epsabs=1e-05, epsrel=1e-05, nk=None)

Return the variance of the displacement field, i.e.:

\[\sigma_{d}^{2} = \frac{1}{6 \pi^{2}} \int dk P(k)\]
Parameters:
  • pk (callable) – Power spectrum.

  • kmin (float, default=1e-7) – Minimum wavenumber.

  • kmax (float, default=1e2) – Maximum wavenumber.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration)

Returns:

sigmad – r.m.s. of the displacement field.

Return type:

float

cosmoprimo.interpolator.integrate_sigma_r2(r, pk, kmin=1e-07, kmax=100.0, method='fftlog', epsabs=1e-05, epsrel=1e-05, nk=None, kernel=<function kernel_tophat2>)

Return the variance of perturbations smoothed by a kernel \(W\) of radius \(r\), i.e.:

\[\sigma_{r}^{2} = \frac{1}{2 \pi^{2}} \int dk k^{2} P(k) W^{2}(kr)\]
Parameters:
  • r (float, array) – Smoothing radius.

  • pk (callable) – Power spectrum.

  • kmin (float, default=1e-7) – Minimum wavenumber.

  • kmax (float, default=1e2) – Maximum wavenumber.

  • epsrel (float, default=1e-5) – Relative precision (for scipy.integrate.quad() integration).

  • kernel (callable, default=kernel_tophat2) – Kernel \(W^{2}\); defaults to (square of) top-hat kernel.

Returns:

sigmar2 – Variance of perturbations.

Return type:

float

cosmoprimo.interpolator.kernel_tophat2(x)

Non-vectorized tophat function.