Zheng07Sats

class halotools.empirical_models.Zheng07Sats(threshold=-20, prim_haloprop_key='halo_mvir', modulate_with_cenocc=False, cenocc_model=None, **kwargs)[source]

Bases: OccupationComponent

Power law model for the occupation statistics of satellite galaxies, introduced in Kravtsov et al. 2004, arXiv:0308519. This implementation uses Zheng et al. 2007, arXiv:0703457, to assign fiducial parameter values.

\(\langle N_{sat} \rangle_{M} = \left( \frac{M - M_{0}}{M_{1}} \right)^{\alpha}\).

Note

The Zheng07Sats model is part of the zheng07 prebuilt composite HOD-style model. For a tutorial on the zheng07 composite model, see Zheng et al. (2007) Composite Model.

Parameters:
thresholdfloat, optional

Luminosity threshold of the mock galaxy sample. If specified, input value must agree with one of the thresholds used in Zheng07 to fit HODs: [-18, -18.5, -19, -19.5, -20, -20.5, -21, -21.5, -22]. Default value is specified in the model_defaults module.

prim_haloprop_keystring, optional

String giving the column name of the primary halo property governing the occupation statistics of gal_type galaxies. Default value is specified in the model_defaults module.

modulate_with_cenoccbool, optional

If set to True, the Zheng07Sats.mean_occupation method will be multiplied by the the first moment of the centrals:

\(\langle N_{\mathrm{sat}}\rangle_{M}\Rightarrow\langle N_{\mathrm{sat}}\rangle_{M}\times\langle N_{\mathrm{cen}}\rangle_{M}\)

The cenocc_model keyword argument works together with the modulate_with_cenocc keyword argument to determine how the \(\langle N_{\mathrm{cen}}\rangle_{M}\) prefactor is calculated.

cenocc_modelOccupationComponent, optional

If the cenocc_model keyword argument is set to its default value of None, then the \(\langle N_{\mathrm{cen}}\rangle_{M}\) prefactor will be calculated according to Zheng07Cens.mean_occupation. However, if an instance of the OccupationComponent class is instead passed in via the cenocc_model keyword, then the first satellite moment will be multiplied by the mean_occupation function of the cenocc_model. The modulate_with_cenocc keyword must be set to True in order for the cenocc_model to be operative. See Advanced usage of the zheng07 model for further details.

Examples

>>> sat_model = Zheng07Sats()
>>> sat_model = Zheng07Sats(threshold = -21)

The param_dict attribute can be used to build an alternate model from an existing instance. This feature has a variety of uses. For example, suppose you wish to study how the choice of halo mass definition impacts HOD predictions:

>>> sat_model1 = Zheng07Sats(threshold = -19.5, prim_haloprop_key='m200b')
>>> sat_model1.param_dict['alpha_satellites'] = 1.05
>>> sat_model2 = Zheng07Sats(threshold = -19.5, prim_haloprop_key='m500c')
>>> sat_model2.param_dict = sat_model1.param_dict

After executing the above four lines of code, sat_model1 and sat_model2 are identical in every respect, excepting only for the difference in the halo mass definition.

A common convention in HOD modeling of satellite populations is for the first occupation moment of the satellites to be multiplied by the first occupation moment of the associated central population. The cenocc_model keyword arguments allows you to study the impact of this choice:

>>> sat_model1 = Zheng07Sats(threshold=-18)
>>> sat_model2 = Zheng07Sats(threshold = sat_model1.threshold, modulate_with_cenocc=True)

Now sat_model1 and sat_model2 are identical in every respect, excepting only the following difference:

\(\langle N_{\mathrm{sat}}\rangle^{\mathrm{model 2}} = \langle N_{\mathrm{cen}}\rangle\times\langle N_{\mathrm{sat}}\rangle^{\mathrm{model 1}}\)

Methods Summary

get_published_parameters(threshold[, ...])

Best-fit HOD parameters from Table 1 of Zheng et al. 2007.

mean_occupation(**kwargs)

Expected number of satellite galaxies in a halo of mass logM.

Methods Documentation

get_published_parameters(threshold, publication='Zheng07')[source]

Best-fit HOD parameters from Table 1 of Zheng et al. 2007.

Parameters:
thresholdfloat

Luminosity threshold of the mock galaxy sample. Input value must agree with one of the thresholds used in Zheng07 to fit HODs: [-18, -18.5, -19, -19.5, -20, -20.5, -21, -21.5, -22].

publicationstring, optional

String specifying the publication that will be used to set the values of param_dict. Default is Zheng et al. (2007).

Returns:
param_dictdict

Dictionary of model parameters whose values have been set to the values taken from Table 1 of Zheng et al. 2007.

Examples

>>> sat_model = Zheng07Sats()
>>> sat_model.param_dict = sat_model.get_published_parameters(sat_model.threshold)
mean_occupation(**kwargs)[source]

Expected number of satellite galaxies in a halo of mass logM. See Equation 5 of arXiv:0703457.

Parameters:
prim_haloproparray, optional

Array storing a mass-like variable that governs the occupation statistics. If prim_haloprop is not passed, then table keyword arguments must be passed.

tableobject, optional

Data table storing halo catalog. If table is not passed, then prim_haloprop keyword arguments must be passed.

Returns:
mean_nsatfloat or array

Mean number of satellite galaxies in a host halo of the specified mass.

\(\langle N_{\mathrm{sat}} \rangle_{M} = \left( \frac{M - M_{0}}{M_{1}} \right)^{\alpha} \langle N_{\mathrm{cen}} \rangle_{M}\)
or
\(\langle N_{\mathrm{sat}} \rangle_{M} = \left( \frac{M - M_{0}}{M_{1}} \right)^{\alpha}\),
depending on whether a central model was passed to the constructor.

Examples

The mean_occupation method of all OccupationComponent instances supports two different options for arguments. The first option is to directly pass the array of the primary halo property:

>>> sat_model = Zheng07Sats()
>>> testmass = np.logspace(10, 15, num=50)
>>> mean_nsat = sat_model.mean_occupation(prim_haloprop = testmass)

The second option is to pass mean_occupation a full halo catalog. In this case, the array storing the primary halo property will be selected by accessing the sat_model.prim_haloprop_key column of the input halo catalog. For illustration purposes, we’ll use a fake halo catalog rather than a (much larger) full one:

>>> from halotools.sim_manager import FakeSim
>>> fake_sim = FakeSim()
>>> mean_nsat = sat_model.mean_occupation(table=fake_sim.halo_table)