NFWProfile¶
- class halotools.empirical_models.NFWProfile(cosmology=FlatLambdaCDM(name='WMAP5', H0=<Quantity 70.2 km / (Mpc s)>, Om0=0.277, Tcmb0=<Quantity 2.725 K>, Neff=3.04, m_nu=<Quantity [0., 0., 0.] eV>, Ob0=0.0459), redshift=0.0, mdef='vir', conc_mass_model='direct_from_halo_catalog', concentration_key='halo_nfw_conc', halo_boundary_key=None, **kwargs)[source]¶
Bases:
AnalyticDensityProf
Model for the spatial distribution of mass and/or galaxies residing in an NFW halo profile, based on Navarro, Frenk and White (1995), arXiv:9508025.
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:
- 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
.- halo_boundary_keystr, optional
Default behavior is to use the column associated with the input mdef.
- 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
.- 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.
Examples
>>> nfw = NFWProfile()
Methods Summary
circular_velocity
(radius, total_mass, conc)The circular velocity, \(V_{\rm cir} \equiv \sqrt{GM(<r)/r}\), as a function of halo-centric distance r.
conc_NFWmodel
(*args, **kwargs)NFW concentration as a function of halo mass.
cumulative_mass_PDF
(scaled_radius, conc)Analytical result for the fraction of the total mass enclosed within dimensionless radius of an NFW halo,
dimensionless_mass_density
(scaled_radius, conc)Physical density of the NFW halo scaled by the density threshold of the mass definition.
enclosed_mass
(radius, total_mass, conc)The mass enclosed within the input radius, \(M(<r) = 4\pi\int_{0}^{r}dr'r'^{2}\rho(r)\).
halo_mass_to_halo_radius
(total_mass)Spherical overdensity radius as a function of the input mass.
halo_radius_to_halo_mass
(radius)Spherical overdensity mass as a function of the input radius.
mass_density
(radius, mass, conc)Physical density of the halo at the input radius, given in units of \(h^{3}/{\rm Mpc}^{3}\).
mc_generate_nfw_radial_positions
(**kwargs)Return a Monte Carlo realization of points in an NFW profile.
rmax
(total_mass, conc)Radius at which the halo attains its maximum circular velocity, \(R_{\rm max}^{\rm NFW} = 2.16258R_{\Delta}/c\).
virial_velocity
(total_mass)The circular velocity evaluated at the halo boundary, \(V_{\rm vir} \equiv \sqrt{GM_{\rm halo}/R_{\rm halo}}\).
vmax
(total_mass, conc)Maximum circular velocity of the halo profile, \(V_{\rm max}^{\rm NFW} = V_{\rm cir}^{\rm NFW}(r = 2.16258R_{\Delta}/c)\).
Methods Documentation
- circular_velocity(radius, total_mass, conc)[source]¶
The circular velocity, \(V_{\rm cir} \equiv \sqrt{GM(<r)/r}\), as a function of halo-centric distance r.
- Parameters:
- radiusarray_like
Halo-centric distance in Mpc/h units; can be a scalar or numpy array
- total_massarray_like
Total mass of the halo; can be a scalar or numpy array of the same dimension as the input
radius
.- concarray_like
Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input
radius
.
- Returns:
- vc: array_like
The circular velocity in km/s; has the same dimensions as the input
radius
.
Notes
See Spatial Profiles of Halos for derivations and implementation details.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> radius = np.logspace(-2, -1, Npts) >>> total_mass = np.zeros(Npts) + 1e12 >>> conc = 5 >>> result = model.circular_velocity(radius, total_mass, conc) >>> concarr = np.linspace(1, 100, Npts) >>> result = model.circular_velocity(radius, total_mass, concarr)
- conc_NFWmodel(*args, **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, thentable
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, thenprim_haloprop
keyword argument must be passed.
- Returns:
- concarray_like
Concentrations of the input halos.
Note that concentrations will be clipped to their min/max permitted values set in the
model_defaults
module. The purpose of this clipping is to ensure stable results during mock galaxy population. Due to this clipping, the behavior of theconc_NFWmodel
function is different from the concentration-mass relation that underlies it.
Examples
In the examples below, we’ll demonstrate the various ways to use the
conc_NFWmodel
function, depending on the initial choice for theconc_mass_model
.>>> fake_masses = np.logspace(12, 15, 10)
If you use the
direct_from_halo_catalog
option, you must pass atable
argument storing aTable
with a column name for the halo mass that is consistent with your chosen halo mass definition:>>> from astropy.table import Table >>> nfw = NFWProfile(conc_mass_model='direct_from_halo_catalog', mdef='vir') >>> fake_conc = np.zeros_like(fake_masses) + 5. >>> fake_halo_table = Table({'halo_mvir': fake_masses, 'halo_nfw_conc': fake_conc}) >>> model_conc = nfw.conc_NFWmodel(table=fake_halo_table)
In case your halo catalog uses a different keyname from the Halotools default
halo_nfw_conc
:>>> nfw = NFWProfile(conc_mass_model='direct_from_halo_catalog', mdef='vir', concentration_key='my_conc_keyname') >>> fake_halo_table = Table({'halo_mvir': fake_masses, 'my_conc_keyname': fake_conc}) >>> model_conc = nfw.conc_NFWmodel(table=fake_halo_table)
One of the available options provided by Halotools is
dutton_maccio14
. With this option, you can either pass in atable
argument, or alternatively an array of masses via theprim_haloprop
argument:>>> nfw = NFWProfile(conc_mass_model='dutton_maccio14') >>> fake_halo_table = Table({'halo_mvir': fake_masses, 'halo_nfw_conc': fake_conc}) >>> model_conc = nfw.conc_NFWmodel(table=fake_halo_table) >>> model_conc = nfw.conc_NFWmodel(prim_haloprop=fake_masses)
Finally, you may also have chosen to define your own concentration-mass relation. If so, your function must at a minimum accept a
table
keyword argument. Below we give a trivial example of using the identity function:>>> def identity_func(*args, **kwargs): return kwargs['table']['halo_mvir'] >>> nfw = NFWProfile(conc_mass_model=identity_func, mdef='vir') >>> fake_halo_table = Table({'halo_mvir': fake_masses}) >>> model_conc = nfw.conc_NFWmodel(table=fake_halo_table)
- cumulative_mass_PDF(scaled_radius, conc)[source]¶
Analytical result for the fraction of the total mass enclosed within dimensionless radius 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}\).See Derivation of the NFW cumulative mass PDF for a derivation of this expression.
- 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.
- 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 inputscaled_radius
.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> scaled_radius = np.logspace(-2, 0, Npts) >>> conc = 5 >>> result = model.cumulative_mass_PDF(scaled_radius, conc) >>> concarr = np.linspace(1, 100, Npts) >>> result = model.cumulative_mass_PDF(scaled_radius, concarr)
- dimensionless_mass_density(scaled_radius, conc)[source]¶
Physical density of the NFW halo scaled by the density threshold of the mass definition.
The
dimensionless_mass_density
is defined as \(\tilde{\rho}_{\rm prof}(\tilde{r}) \equiv \rho_{\rm prof}(\tilde{r}) / \rho_{\rm thresh}\), where \(\tilde{r}\equiv r/R_{\Delta}\).For an NFW halo, \(\tilde{\rho}_{\rm NFW}(\tilde{r}, c) = \frac{c^{3}/3g(c)}{c\tilde{r}(1 + c\tilde{r})^{2}},\)
where \(g(x) \equiv \log(1+x) - x / (1+x)\) is computed using the
g
function.The quantity \(\rho_{\rm thresh}\) is a function of the halo mass definition, cosmology and redshift, and is computed via the
density_threshold
function. The quantity \(\rho_{\rm prof}\) is the physical mass density of the halo profile and is computed via themass_density
function. See Modeling the NFW Spatial Profile for a derivation of this expression.- 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.
- 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:
- dimensionless_density: array_like
Dimensionless density of a dark matter halo at the input
scaled_radius
, normalized by thedensity_threshold
\(\rho_{\rm thresh}\) for the halo mass definition, cosmology, and redshift. Result is an array of the dimension as the inputscaled_radius
.
- enclosed_mass(radius, total_mass, conc)[source]¶
The mass enclosed within the input radius, \(M(<r) = 4\pi\int_{0}^{r}dr'r'^{2}\rho(r)\).
- Parameters:
- radiusarray_like
Halo-centric distance in Mpc/h units; can be a scalar or numpy array
- total_massarray_like
Total mass of the halo; can be a scalar or numpy array of the same dimension as the input
radius
.- concarray_like
Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input
radius
.
- Returns:
- enclosed_mass: array_like
The mass enclosed within radius r, in \(M_{\odot}/h\); has the same dimensions as the input
radius
.
Notes
See Spatial Profiles of Halos for derivations and implementation details.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> radius = np.logspace(-2, -1, Npts) >>> total_mass = np.zeros(Npts) + 1e12 >>> conc = 5 >>> result = model.enclosed_mass(radius, total_mass, conc) >>> concarr = np.linspace(1, 100, Npts) >>> result = model.enclosed_mass(radius, total_mass, concarr)
- halo_mass_to_halo_radius(total_mass)[source]¶
Spherical overdensity radius as a function of the input mass.
Note that this function is independent of the form of the density profile.
- Parameters:
- total_mass: array_like
Total halo mass in \(M_{\odot}/h\); can be a number or a numpy array.
- Returns:
- radiusarray_like
Radius of the halo in Mpc/h units. Will have the same dimension as the input
total_mass
.
Examples
>>> model = NFWProfile() >>> halo_radius = model.halo_mass_to_halo_radius(1e13)
- halo_radius_to_halo_mass(radius)[source]¶
Spherical overdensity mass as a function of the input radius.
Note that this function is independent of the form of the density profile.
- Parameters:
- radiusarray_like
Radius of the halo in Mpc/h units; can be a number or a numpy array.
- Returns:
- total_mass: array_like
Total halo mass in \(M_{\odot}/h\). Will have the same dimension as the input
radius
.
Examples
>>> model = NFWProfile() >>> halo_mass = model.halo_mass_to_halo_radius(500.)
- mass_density(radius, mass, conc)[source]¶
Physical density of the halo at the input radius, given in units of \(h^{3}/{\rm Mpc}^{3}\).
- Parameters:
- radiusarray_like
Halo-centric distance in Mpc/h units; can be a scalar or numpy array
- massarray_like
Total mass of the halo; can be a scalar or numpy array of the same dimension as the input
radius
.- concarray_like
Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input
radius
.
- Returns:
- density: array_like
Physical density of a dark matter halo of the input
mass
at the inputradius
. Result is an array of the dimension as the inputradius
, reported in units of \(h^{3}/Mpc^{3}\).
Notes
See Spatial Profiles of Halos for derivations and implementation details.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> radius = np.logspace(-2, -1, Npts) >>> mass = np.zeros(Npts) + 1e12 >>> conc = 5 >>> result = model.mass_density(radius, mass, conc) >>> concarr = np.linspace(1, 100, Npts) >>> result = model.mass_density(radius, mass, concarr)
- mc_generate_nfw_radial_positions(**kwargs)[source]¶
Return a Monte Carlo realization of points in an NFW profile.
See Monte Carlo realizations of the NFW profile for a discussion of this technique.
- Parameters:
- num_ptsint, optional
Number of points in the Monte Carlo realization of the profile. Default is 1e4.
- concfloat, optional
Concentration of the NFW profile being realized. Default is 5.
- halo_massfloat, optional
Total mass of the halo whose profile is being realized, used to define the halo boundary for the mass definition bound to the NFWProfile instance as
mdef
.If
halo_mass
is unspecified, keyword argumenthalo_radius
must be specified.- halo_radiusfloat, optional
Physical boundary of the halo whose profile is being realized in units of Mpc/h.
If
halo_radius
is unspecified, keyword argumenthalo_mass
must be specified, in which case the outer boundary of the halo will be determined according to the mass definition bound to the NFWProfile instance asmdef
.- seedint, optional
Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.
- Returns:
- radial_positionsarray_like
Numpy array storing a Monte Carlo realization of the halo profile. All values will lie strictly between 0 and the halo boundary.
Examples
>>> nfw = NFWProfile() >>> radial_positions = nfw.mc_generate_nfw_radial_positions(halo_mass = 1e12, conc = 10) >>> radial_positions = nfw.mc_generate_nfw_radial_positions(halo_radius = 0.25)
- rmax(total_mass, conc)[source]¶
Radius at which the halo attains its maximum circular velocity, \(R_{\rm max}^{\rm NFW} = 2.16258R_{\Delta}/c\).
- Parameters:
- total_mass: array_like
Total halo mass in \(M_{\odot}/h\); can be a number or a numpy array.
- concarray_like
Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input
total_mass
.
- Returns:
- rmaxarray_like
\(R_{\rm max}\) in Mpc/h.
Notes
See Spatial Profiles of Halos for derivations and implementation details.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> total_mass = np.zeros(Npts) + 1e12 >>> conc = 5 >>> result = model.rmax(total_mass, conc) >>> concarr = np.linspace(1, 100, Npts) >>> result = model.rmax(total_mass, concarr)
- virial_velocity(total_mass)[source]¶
The circular velocity evaluated at the halo boundary, \(V_{\rm vir} \equiv \sqrt{GM_{\rm halo}/R_{\rm halo}}\).
- Parameters:
- total_massarray_like
Total mass of the halo; can be a scalar or numpy array.
- Returns:
- vvirarray_like
Virial velocity in km/s.
Notes
See Spatial Profiles of Halos for derivations and implementation details.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> mass_array = np.logspace(11, 15, Npts) >>> vvir_array = model.virial_velocity(mass_array)
- vmax(total_mass, conc)[source]¶
Maximum circular velocity of the halo profile, \(V_{\rm max}^{\rm NFW} = V_{\rm cir}^{\rm NFW}(r = 2.16258R_{\Delta}/c)\).
- Parameters:
- total_mass: array_like
Total halo mass in \(M_{\odot}/h\); can be a number or a numpy array.
- concarray_like
Value of the halo concentration. Can either be a scalar, or a numpy array of the same dimension as the input
total_mass
.
- Returns:
- vmaxarray_like
\(V_{\rm max}\) in km/s.
Notes
See Spatial Profiles of Halos for derivations and implementation details.
Examples
>>> model = NFWProfile() >>> Npts = 100 >>> total_mass = np.zeros(Npts) + 1e12 >>> conc = 5 >>> result = model.vmax(total_mass, conc) >>> concarr = np.linspace(1, 100, Npts) >>> result = model.vmax(total_mass, concarr)