SFRBiasedNFWPhaseSpace¶
- class halotools.empirical_models.SFRBiasedNFWPhaseSpace(**kwargs)[source]¶
Bases:
BiasedNFWPhaseSpace
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, independently for red and blue galaxies.
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.
Notes
The
SFRBiasedNFWPhaseSpace
class can only be used to make mocks in concert with some other component model that is responsible for modeling anquiescent
property of thegalaxy_table
.- 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 theparam_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 insim_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
.- concentration_keystring, optional
Column name of the halo catalog storing NFW concentration.
This argument is only relevant when
conc_mass_model
is set todirect_from_halo_catalog
. In such a case, the default value ishalo_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 = SFRBiasedNFWPhaseSpace()
The behavior of the
SFRBiasedNFWPhaseSpace
model is controlled by the values stored in itsparam_dict
. In the above default model, we have two parameters that control the concentration of the satellites,quiescent_conc_gal_bias
andstar_forming_conc_gal_bias
. The value of \(F_{\rm q}\) =quiescent_conc_gal_bias
sets the value for the concentration of quiescent satellites according to \(F_{\rm q} * c_{\rm halo},\) where \(c_{\rm halo}\) is the value of the halo concentration specified by the concentration-mass model.By default, \(F_{\rm q} = F_{\rm sf} = 1\). Quiescent satellites with larger values of \(F_{\rm q}\) will be more radially concentrated and have smaller radial velocity dispersions, and likewise for star-forming galaxies.
The
SFRBiasedNFWPhaseSpace
gives you the option to allow \(F_{\rm q}\) and \(F_{\rm sf}\) to vary with halo mass. You can accomplish this via theconc_gal_bias_logM_abscissa
keyword:>>> biased_nfw_mass_dep = SFRBiasedNFWPhaseSpace(conc_gal_bias_logM_abscissa=[12., 15.])
In the above model, the values of \(F_{\rm q} = F_{\rm q}(M_{\rm halo})\) and \(F_{\rm sf} = F_{\rm sf}(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 theparam_dict
controlling the value of \(F_{\rm q}\) and \(F_{\rm sf}\) at that mass.For example, in the
biased_nfw_mass_dep
defined above, to set the value of \(F_{\rm q}(M_{\rm halo} = 10^{15})\):>>> biased_nfw_mass_dep.param_dict['quiescent_conc_gal_bias_param1'] = 2
To triple the value of \(F_{\rm sf}(M_{\rm halo} = 10^{12})\):
>>> biased_nfw_mass_dep.param_dict['star_forming_conc_gal_bias_param0'] *= 3
Values of \(F_{\rm q}\) and \(F_{\rm sf}\) at masses between the control points will be determined by log-linear interpolation. When extrapolating \(F_{\rm q}\) and \(F_{\rm sf}\) 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}\).
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, thentable
keyword argument must be passed.- quiescentarray, optional
Boolean array storing whether the galaxy is quiescent. Must be passed together with prim_haloprop argument.
- 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, thenprim_haloprop
andquiescent
keyword arguments must be passed.
- Returns:
- conc_gal_biasarray_like
Ratio of the galaxy concentration to the halo concentration, \(F_{\rm gal} = c_{\rm gal}/c_{\rm halo}\).
Examples
>>> model = SFRBiasedNFWPhaseSpace() >>> mass = np.logspace(10, 15, 100) >>> quiescent = np.zeros(100).astype(bool) >>> quiescent[::2] = True >>> conc_gal_bias = model.calculate_conc_gal_bias(prim_haloprop=mass, quiescent=quiescent)