NFWPhaseSpace

class halotools.empirical_models.NFWPhaseSpace(**kwargs)[source]

Bases: NFWProfile, MonteCarloGalProf

Model for the phase space distribution of mass and/or galaxies in isotropic Jeans equilibrium in an NFW halo profile, based on Navarro, Frenk and White (1995), where the concentration of the galaxies is the same as the concentration of the parent halo

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_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.

Examples

>>> model = NFWPhaseSpace()

Methods Summary

assign_phase_space(table[, seed])

Primary method of the NFWPhaseSpace class called during the mock-population sequence.

build_lookup_tables([logrmin, logrmax, ...])

Method used to create a lookup table of the spatial and velocity radial profiles.

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_gal_PDF(scaled_radius, conc)

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

cumulative_mass_PDF(scaled_radius, conc)

Analytical result for the fraction of the total mass enclosed within r/Rvir 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:

dimensionless_radial_velocity_dispersion(...)

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

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_phase_space_points([Ngals, ...])

Return a Monte Carlo realization of points in the phase space of an NFW halo in isotropic Jeans equilibrium.

mc_halo_centric_pos(*concentration_array, ...)

Method to generate random, three-dimensional halo-centric positions of galaxies.

mc_pos(*concentration_array, **kwargs)

Method to generate random, three-dimensional positions of galaxies.

mc_radial_velocity(scaled_radius, ...)

Method returns a Monte Carlo realization of radial velocities drawn from Gaussians with a width determined by the solution to the isotropic Jeans equation.

mc_solid_sphere(*concentration_array, **kwargs)

Method to generate random, three-dimensional, halo-centric positions of galaxies.

mc_unit_sphere(Npts, **kwargs)

Returns Npts random points on the unit sphere.

mc_vel(table[, seed])

Method assigns a Monte Carlo realization of the Jeans velocity solution to the halos in the input table.

radial_velocity_dispersion(radius, ...)

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

setup_prof_lookup_tables(*concentration_bins)

This method sets up how we will digitize halo concentrations during mock population.

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.

Methods Documentation

assign_phase_space(table, seed=None)[source]

Primary method of the NFWPhaseSpace class called during the mock-population sequence.

Parameters:
tableobject

Table storing halo catalog.

After calling the assign_phase_space method, the x, y, z, vx, vy, and vz columns of the input table will be over-written with their host-centric values.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

build_lookup_tables(logrmin=-3, logrmax=0, Npts_radius_table=101)[source]

Method used to create a lookup table of the spatial and velocity radial profiles.

Parameters:
logrminfloat, optional

Minimum radius used to build the spline table. Default is set in model_defaults.

logrmaxfloat, optional

Maximum radius used to build the spline table Default is set in model_defaults.

Npts_radius_tableint, optional

Number of control points used in the spline. Default is set in model_defaults.

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 = NFWPhaseSpace()
>>> 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, 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:
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 the conc_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 the conc_mass_model.

>>> fake_masses = np.logspace(12, 15, 10)

If you use the direct_from_halo_catalog option, you must pass a table argument storing a Table with a column name for the halo mass that is consistent with your chosen halo mass definition:

>>> from astropy.table import Table
>>> nfw = NFWPhaseSpace(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 = NFWPhaseSpace(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 a table argument, or alternatively an array of masses via the prim_haloprop argument:

>>> nfw = NFWPhaseSpace(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 = NFWPhaseSpace(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_gal_PDF(scaled_radius, conc)[source]

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

In NFWPhaseSpace there is no distinction between the two methods, but in BiasedNFWPhaseSpace these two function are different.

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 input scaled_radius.

Examples

>>> model = NFWPhaseSpace()
>>> Npts = 100
>>> scaled_radius = np.logspace(-2, 0, Npts)
>>> conc = 5
>>> result = model.cumulative_gal_PDF(scaled_radius, conc)
>>> concarr = np.linspace(1, 100, Npts)
>>> result = model.cumulative_gal_PDF(scaled_radius, concarr)
cumulative_mass_PDF(scaled_radius, 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.

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)
>>> 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)}\times\frac{1}{c\tilde{r}(1 + c\tilde{r})^{2}},\)

where \(g(x) \equiv \int_{0}^{x}dy\frac{y}{(1+y)^{2}} = \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 the mass_density function.

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 the density_threshold \(\rho_{\rm thresh}\) for the halo mass definition, cosmology, and redshift. Result is an array of the dimension as the input scaled_radius.

dimensionless_radial_velocity_dispersion(scaled_radius, *conc)[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\).

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.

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 = NFWPhaseSpace()
>>> 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 = NFWPhaseSpace()
>>> 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 input radius. Result is an array of the dimension as the input radius, reported in units of \(h^{3}/Mpc^{3}\).

Notes

See Spatial Profiles of Halos for derivations and implementation details.

Examples

>>> model = NFWPhaseSpace()
>>> 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_phase_space_points(Ngals=10000, conc=5, mass=1000000000000.0, verbose=True, seed=None)[source]

Return a Monte Carlo realization of points in the phase space of an NFW halo in isotropic Jeans equilibrium.

Parameters:
Ngalsint, optional

Number of galaxies in the Monte Carlo realization of the phase space distribution. Default is 1e4.

concfloat, optional

Concentration of the NFW profile being realized. Default is 5.

massfloat, optional

Mass of the halo whose phase space distribution is being realized in units of Msun/h. Default is 1e12.

verbosebool, optional

If True, a message prints with an estimate of the build time. Default is True.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

Returns:
ttable

Table containing the Monte Carlo realization of the phase space distribution. Keys are ‘x’, ‘y’, ‘z’, ‘vx’, ‘vy’, ‘vz’, ‘radial_position’, ‘radial_velocity’. Length units in Mpc/h, velocity units in km/s.

Sign convention on the returned radial_velocity column is such that negative (positive) values correspond to satellites moving radially inward (outward)

Examples

>>> nfw = NFWPhaseSpace()
>>> mass, conc = 1e13, 8.
>>> data = nfw.mc_generate_nfw_phase_space_points(Ngals=100, mass=mass, conc=conc, verbose=False)

Now suppose you wish to compute the radial velocity dispersion of all the returned points:

>>> vrad_disp = np.std(data['radial_velocity'])

If you wish to do the same calculation but for points in a specific range of radius:

>>> mask = data['radial_position'] < 0.1
>>> vrad_disp_inner_points = np.std(data['radial_velocity'][mask])

You may also wish to select points according to their distance to the halo center in units of the virial radius. In such as case, you can use the halo_mass_to_halo_radius method to scale the halo-centric distances. Here is an example of how to compute the velocity dispersion in the z-dimension of all points residing within \(R_{\rm vir}/2\):

>>> halo_radius = nfw.halo_mass_to_halo_radius(mass)
>>> scaled_radial_positions = data['radial_position']/halo_radius
>>> mask = scaled_radial_positions < 0.5
>>> vz_disp_inner_half = np.std(data['vz'][mask])
mc_halo_centric_pos(*concentration_array, **kwargs)[source]

Method to generate random, three-dimensional halo-centric positions of galaxies.

Parameters:
tabledata table, optional

Astropy Table storing a length-Ngals galaxy catalog. If table is not passed, concentration_array and keyword argument halo_radius must be passed.

concentration_arrayarray_like, optional

Length-Ngals numpy array storing the concentrations of the mock galaxies. If table is not passed, concentration_array and keyword argument halo_radius must be passed.

halo_radiusarray_like, optional

Length-Ngals array storing the radial boundary of the halo hosting each galaxy. Units assumed to be in Mpc/h. If concentration_array and halo_radius are not passed, table must be passed.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

Returns:
x, y, zarrays

Length-Ngals array storing a Monte Carlo realization of the galaxy positions.

mc_pos(*concentration_array, **kwargs)[source]

Method to generate random, three-dimensional positions of galaxies.

Parameters:
tabledata table, optional

Astropy Table storing a length-Ngals galaxy catalog. If table is not passed, concentration_array and halo_radius must be passed.

concentration_arrayarray_like, optional

Length-Ngals numpy array storing the concentrations of the mock galaxies. If table is not passed, concentration_array and keyword argument halo_radius must be passed. If concentration_array is passed, halo_radius must be passed as a keyword argument. The sequence must have the same order as self.gal_prof_param_keys.

halo_radiusarray_like, optional

Length-Ngals array storing the radial boundary of the halo hosting each galaxy. Units assumed to be in Mpc/h. If concentration_array and halo_radius are not passed, table must be passed.

seedint, optional

Random number seed used in Monte Carlo realization. Default is None.

Returns:
x, y, zarrays, optional

For the case where no table is passed as an argument, method will return x, y and z points distributed about the origin according to the profile model.

For the case where table is passed as an argument (this is the use case of populating halos with mock galaxies), the x, y, and z columns of the table will be over-written. When table is passed as an argument, the method assumes that the x, y, and z columns already store the position of the host halo center.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

mc_radial_velocity(scaled_radius, total_mass, *concentration_array, **kwargs)[source]

Method returns a Monte Carlo realization of radial velocities drawn from Gaussians with a width determined by the solution to the isotropic Jeans equation.

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.

total_mass: array_like

Length-Ngals numpy array storing the halo mass in \(M_{\odot}/h\).

concentration_arrayarray_like

Length-Ngals numpy array storing the concentrations of the mock galaxies.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

Returns:
radial_velocitiesarray_like

Array of radial velocities drawn from Gaussians with a width determined by the solution to the Jeans equation.

mc_solid_sphere(*concentration_array, **kwargs)[source]

Method to generate random, three-dimensional, halo-centric positions of galaxies.

Parameters:
concentration_arrayarray_like, optional

Length-Ngals numpy array storing the concentrations of the mock galaxies.

tabledata table, optional

Astropy Table storing a length-Ngals galaxy catalog. If table is not passed, concentration_array must be passed.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

Returns:
x, y, zarrays

Length-Ngals array storing a Monte Carlo realization of the galaxy positions.

mc_unit_sphere(Npts, **kwargs)[source]

Returns Npts random points on the unit sphere.

Parameters:
Nptsint

Number of 3d points to generate

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

Returns:
x, y, zarray_like

Length-Npts arrays of the coordinate positions.

mc_vel(table, seed=None)[source]

Method assigns a Monte Carlo realization of the Jeans velocity solution to the halos in the input table.

Parameters:
tableAstropy Table

astropy.table.Table object storing the halo catalog. Calling the mc_vel method will over-write the existing values of the vx, vy and vz columns.

seedint, optional

Random number seed used in the Monte Carlo realization. Default is None, which will produce stochastic results.

radial_velocity_dispersion(radius, total_mass, halo_conc)[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

Returns:
resultarray_like

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

setup_prof_lookup_tables(*concentration_bins)[source]

This method sets up how we will digitize halo concentrations during mock population.

Parameters:
concentration_binsndarray

Array storing how concentrations will be digitized during mock-population

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.

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 = NFWPhaseSpace()
>>> 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)