halomod.bias.Jing98

class halomod.bias.Jing98(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 of Jing (1998).

See documentation for Bias for information on input parameters. This model has no free parameters.

Notes

This is an empirical form proposed in [1], with the formula

\[(a/\nu^4 + 1)^{b - c n} \left(1 + \frac{\nu^2 - 1}{\delta_c}\right)\]

The parameters a, b and c are free parameters, with values fitted in [1] of (0.5, 0.06, 0.02), which are the defaults here.

Parameters:
  • a (float) – The fitting parameters.

  • b (float) – The fitting parameters.

  • c (float) – The fitting parameters.

References

[1] (1,2)

Jing, Y. P., “Accurate Fitting Formula for the Two-Point Correlation Function of Dark Matter Halos”, http://adsabs.harvard.edu/abs/1998ApJ…503L…9J, 1998.

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