BiasedNFWPhaseSpace

class halotools.empirical_models.BiasedNFWPhaseSpace(profile_integration_tol=1e-05, **kwargs)[source]

Bases: NFWPhaseSpace

Model for the phase space distribution of galaxies in isotropic Jeans equilibrium in an NFW halo profile, based on Navarro, Frenk and White (1995), where the concentration of the tracers is permitted to differ from the host halo concentration.

For a review of the mathematics underlying the NFW profile, including descriptions of how the relevant equations are implemented in the Halotools code base, see Source code notes on NFWProfile and NFWPhaseSpace.

Parameters:
conc_gal_bias_logM_abscissaarray_like, optional

Numpy array of shape (num_gal_bias_bins, ) storing the values of log10(Mhalo). For each entry of conc_gal_bias_logM_abscissa, there will be a corresponding parameter in the param_dict allowing you to vary the strength of the galaxy concentration bias, using log-linear interpolation for the intermediate values between these control points.

conc_mass_modelstring or callable, optional

Specifies the function used to model the relation between NFW concentration and halo mass. Can either be a custom-built callable function, or one of the following strings: dutton_maccio14, direct_from_halo_catalog.

cosmologyobject, optional

Instance of an astropy cosmology. Default cosmology is set in sim_defaults.

redshiftfloat, optional

Default is set in sim_defaults.

mdef: str, optional

String specifying the halo mass definition, e.g., ‘vir’ or ‘200m’. Default is set in model_defaults.

halo_boundary_keystr, optional

Default behavior is to use the column associated with the input mdef.

concentration_keystring, optional

Column name of the halo catalog storing NFW concentration.

This argument is only relevant when conc_mass_model is set to direct_from_halo_catalog. In such a case, the default value is halo_nfw_conc, which is consistent with all halo catalogs provided by Halotools but may differ from the convention adopted in custom halo catalogs.

concentration_binsndarray, optional

Array storing how halo concentrations will be digitized when building a lookup table for mock-population purposes. The spacing of this array sets a limit on how accurately the concentration parameter can be recovered in a likelihood analysis.

conc_gal_bias_binsndarray, optional

Array storing how biases in galaxy concentrations will be digitized when building a lookup table for mock-population purposes. The spacing of this array sets a limit on how accurately the galaxy bias parameter can be recovered in a likelihood analysis.

profile_integration_tolfloat, optional

Default is 1e-5

Examples

>>> biased_nfw = BiasedNFWPhaseSpace()

The behavior of the BiasedNFWPhaseSpace model is controlled by the values stored in its param_dict. In the above default model, we have one parameter that controls the concentration of the satellites, conc_gal_bias. The value of \(F_{\rm gal}\) = conc_gal_bias sets the value for the concentration of satellites according to \(F_{\rm gal} * c_{\rm halo},\) where \(c_{\rm halo}\) is the value of the halo concentration specified by the concentration-mass model.

By default, \(F_{\rm gal} = 1\). Satellites with larger values of \(F_{\rm gal}\) will be more radially concentrated and have smaller radial velocity dispersions.

The BiasedNFWPhaseSpace gives you the option to allow \(F_{\rm gal}\) to vary with halo mass. You can accomplish this via the conc_gal_bias_logM_abscissa keyword:

>>> biased_nfw_mass_dep = BiasedNFWPhaseSpace(conc_gal_bias_logM_abscissa=[12., 15.])

In the above model, the values of \(F_{\rm gal} = F_{\rm gal}(M_{\rm halo})\) can be independently specified at either of the two control points, \(10^{12}M_{\odot}\) or \(10^{15}M_{\odot}\). For every element in the input conc_gal_bias_logM_abscissa, there is a correposponding value in the param_dict controlling the value of \(F_{\rm gal}\) at that mass.

For example, in the biased_nfw_mass_dep defined above, to set the value of \(F_{\rm gal}(M_{\rm halo} = 10^{15})\):

>>> biased_nfw_mass_dep.param_dict['conc_gal_bias_param1'] = 2

To triple the value of \(F_{\rm gal}(M_{\rm halo} = 10^{12})\):

>>> biased_nfw_mass_dep.param_dict['conc_gal_bias_param0'] *= 3

Values of \(F_{\rm gal}\) at masses between the control points will be determined by log-linear interpolation. When extrapolating \(F_{\rm gal}\) beyond the specified range, the values will be kept constant at the end point values.

Methods Summary

calculate_conc_gal_bias([seed])

Calculate the ratio of the galaxy concentration to the halo concentration, \(c_{\rm gal}/c_{\rm halo}\).

conc_NFWmodel(**kwargs)

NFW concentration as a function of halo mass.

cumulative_gal_PDF(scaled_radius, halo_conc, ...)

Analogous to cumulative_mass_PDF, but for the satellite galaxy distribution instead of the host halo mass distribution.

cumulative_mass_PDF(scaled_radius, halo_conc)

Analytical result for the fraction of the total mass enclosed within r/Rvir of an NFW halo,

dimensionless_radial_velocity_dispersion(...)

Analytical solution to the isotropic jeans equation for an NFW potential, rendered dimensionless via scaling by the virial velocity.

radial_velocity_dispersion(radius, ...)

Method returns the radial velocity dispersion scaled by the virial velocity as a function of the halo-centric distance.

Methods Documentation

calculate_conc_gal_bias(seed=None, **kwargs)[source]

Calculate the ratio of the galaxy concentration to the halo concentration, \(c_{\rm gal}/c_{\rm halo}\).

Parameters:
prim_haloproparray, optional

Array storing the mass-like variable, e.g., halo_mvir.

If prim_haloprop is not passed, then table keyword argument must be passed.

tableobject, optional

Table storing the halo catalog.

If your NFW model is based on the virial definition, then halo_mvir must appear in the input table, and likewise for other halo mass definitions.

If table is not passed, then prim_haloprop keyword argument must be passed.

Returns:
conc_gal_biasarray_like

Ratio of the galaxy concentration to the halo concentration, \(c_{\rm gal}/c_{\rm halo}\).

conc_NFWmodel(**kwargs)[source]

NFW concentration as a function of halo mass.

Parameters:
prim_haloproparray, optional

Array storing the mass-like variable, e.g., halo_mvir.

If prim_haloprop is not passed, then table keyword argument must be passed.

tableobject, optional

Table storing the halo catalog.

If your NFW model is based on the virial definition, then halo_mvir must appear in the input table, and likewise for other halo mass definitions.

If table is not passed, then prim_haloprop keyword argument must be passed.

Returns:
halo_concarray_like

Concentrations of the halos.

Note that concentrations will be clipped to their min/max permitted values set in the model_defaults module. For unclipped concentrations, set these variables to -np.inf, +np.inf.

cumulative_gal_PDF(scaled_radius, halo_conc, conc_gal_bias)[source]

Analogous to cumulative_mass_PDF, but for the satellite galaxy distribution instead of the host halo mass distribution.

Parameters:
scaled_radiusarray_like

Halo-centric distance r scaled by the halo boundary \(R_{\Delta}\), so that \(0 <= \tilde{r} \equiv r/R_{\Delta} <= 1\). Can be a scalar or numpy array.

halo_concarray_like

Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input scaled_radius.

conc_gal_biasarray_like

Ratio of the galaxy concentration to the halo concentration, \(c_{\rm gal}/c_{\rm halo}\).

Returns:
p: array_like

The fraction of the total mass enclosed within the input scaled_radius, in \(M_{\odot}/h\); has the same dimensions as the input scaled_radius.

Examples

>>> model = BiasedNFWPhaseSpace()
>>> scaled_radius = 0.5  # units of Rvir
>>> halo_conc = 5
>>> conc_gal_bias = 3
>>> result1 = model.cumulative_gal_PDF(scaled_radius, halo_conc, conc_gal_bias)
>>> num_halos = 50
>>> scaled_radius = np.logspace(-2, 0, num_halos)
>>> halo_conc = np.linspace(1, 25, num_halos)
>>> conc_gal_bias = np.zeros(num_halos) + 2.
>>> result2 = model.cumulative_gal_PDF(scaled_radius, halo_conc, conc_gal_bias)
cumulative_mass_PDF(scaled_radius, halo_conc)[source]

Analytical result for the fraction of the total mass enclosed within r/Rvir of an NFW halo,

\(P_{\rm NFW}(<\tilde{r}) \equiv M_{\Delta}(<\tilde{r}) / M_{\Delta} = g(c\tilde{r})/g(\tilde{r}),\)

where \(g(x) \equiv \int_{0}^{x}dy\frac{y}{(1+y)^{2}} = \log(1+x) - x / (1+x)\) is computed using g, and where \(\tilde{r} \equiv r / R_{\Delta}\).

Parameters:
scaled_radiusarray_like

Halo-centric distance r scaled by the halo boundary \(R_{\Delta}\), so that \(0 <= \tilde{r} \equiv r/R_{\Delta} <= 1\). Can be a scalar or numpy array.

halo_concarray_like

Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input scaled_radius.

Returns:
p: array_like

The fraction of the total mass enclosed within the input scaled_radius, in \(M_{\odot}/h\); has the same dimensions as the input scaled_radius.

Notes

See Spatial Profiles of Halos for derivations and implementation details.

Examples

>>> model = NFWPhaseSpace()
>>> Npts = 100
>>> scaled_radius = np.logspace(-2, 0, Npts)
>>> halo_conc = 5
>>> result = model.cumulative_mass_PDF(scaled_radius, halo_conc)
>>> halo_conc_arr = np.linspace(1, 100, Npts)
>>> result = model.cumulative_mass_PDF(scaled_radius, halo_conc_arr)
dimensionless_radial_velocity_dispersion(scaled_radius, halo_conc, conc_gal_bias)[source]

Analytical solution to the isotropic jeans equation for an NFW potential, rendered dimensionless via scaling by the virial velocity.

\(\tilde{\sigma}^{2}_{r}(\tilde{r})\equiv\sigma^{2}_{r}(\tilde{r})/V_{\rm vir}^{2} = \frac{c^{2}\tilde{r}(1 + c\tilde{r})^{2}}{g(c)}\int_{c\tilde{r}}^{\infty}{\rm d}y\frac{g(y)}{y^{3}(1 + y)^{2}}\)

See Modeling the NFW Velocity Profile for derivations and implementation details.

Parameters:
scaled_radiusarray_like

Length-Ngals numpy array storing the halo-centric distance r scaled by the halo boundary \(R_{\Delta}\), so that \(0 <= \tilde{r} \equiv r/R_{\Delta} <= 1\).

halo_concfloat

Concentration of the halo.

Returns:
resultarray_like

Radial velocity dispersion profile scaled by the virial velocity. The returned result has the same dimension as the input scaled_radius.

radial_velocity_dispersion(radius, total_mass, halo_conc, conc_gal_bias)[source]

Method returns the radial velocity dispersion scaled by the virial velocity as a function of the halo-centric distance.

Parameters:
radiusarray_like

Radius of the halo in Mpc/h units; can be a float or ndarray of shape (num_radii, )

total_massarray_like

Float or ndarray of shape (num_radii, ) storing the host halo mass

halo_concarray_like

Float or ndarray of shape (num_radii, ) storing the host halo concentration

conc_gal_biasarray_like

Ratio of the galaxy concentration to the halo concentration, \(c_{\rm gal}/c_{\rm halo}\).

Returns:
resultarray_like

Radial velocity dispersion profile as a function of the input radius, in units of km/s.