halomod.bias.Seljak04¶
- class halomod.bias.Seljak04(nu: ~numpy.ndarray, delta_c: float = 1.686, m: ~numpy.ndarray | None = None, mstar: float | None = None, delta_halo: float | None = 200, n: float | None = 1, sigma_8: float | None = 0.8, cosmo: ~astropy.cosmology.flrw.base.FLRW = FlatLambdaCDM(name='Planck15', H0=<Quantity 67.74 km / (Mpc s)>, Om0=0.3075, Tcmb0=<Quantity 2.7255 K>, Neff=3.046, m_nu=<Quantity [0., 0., 0.06] eV>, Ob0=0.0486), n_eff: None | ~numpy.ndarray = None, z: float = 0.0, **model_parameters)[source]¶
Bases:
Bias
Empirical bias relation from Seljak & Warren (2004), without cosmological dependence.
See documentation for
Bias
for information on input parameters. This model has no free parameters.Notes
This the form from [1] without cosmological dependence. The form is
\[a + bx^c + \frac{d}{ex+1} + fx^g\]with \(x = m/m_\star\) (and \(m_star\) the nonlinear mass – see
Bias
for details). The other parameters are all fitted, with values given [1] as(a,b,c,d,e,f,g) = (0.53, 0.39, 0.45, 0.13, 40, 5e-4, 1.5)
.- Parameters:
a (float, optional) – The fitted parameters.
b (float, optional) – The fitted parameters.
c (float, optional) – The fitted parameters.
d (float, optional) – The fitted parameters.
e (float, optional) – The fitted parameters.
f (float, optional) – The fitted parameters.
g (float, optional) – The fitted parameters.
References
[1] (1,2)Seljak, U. and Warren M. S., “Large-scale bias and stochasticity of haloes and dark matter”, https://ui.adsabs.harvard.edu/abs/2004MNRAS.355..129S, 2004.
- bias()[source]¶
Calculate the first-order, linear, deterministic halo bias.
- Returns:
b – The bias as a function of mass, as an array of values corresponding to the instance attributes m and/or nu.
- Return type:
array-like
Examples
>>> import matplotlib.pyplot as plt >>> import numpy as np >>> from halomod.bias import Mo96 >>> peak_height = np.linspace(0.1, 2, 100) >>> bias = Mo96(nu=peak_height) >>> plt.plot(peak_height, bias.bias())
- classmethod get_models() Dict[str, Type] ¶
Get a dictionary of all implemented models for this component.
- pair_hmf = ()¶
The HMF model that pairs with this bias in the peak-background split