SciPy

Navarro-Frenk-White profile

This module implements the Navarro-Frenk-White form of the density profile. Please see Halo Density Profiles for a general introduction to the Colossus density profile module.

Basics

The NFW profile (Navarro et al. 1997) is defined by the density function

\[\rho(r) = \frac{\rho_s}{\left(\frac{r}{r_{\rm s}}\right) \left(1 + \frac{r}{r_s}\right)^{2}}\]

The profile class can be initialized by either passing its fundamental parameters \(\rho_{\rm s}\) and \(r_{\rm s}\), but the more convenient initialization is via mass and concentration:

from colossus.cosmology import cosmology
from colossus.halo import profile_nfw

cosmology.setCosmology('planck18')
p_nfw = profile_nfw.NFWProfile(M = 1E12, c = 10.0, z = 0.0, mdef = 'vir')

The NFW profile class is optimized by using analytical expressions instead of numerical calculations wherever possible. The radiusFromPdf() function covers the common case of drawing random radial positions given an NFW profile.

Please see the Tutorials for more code examples.

Module reference

class halo.profile_nfw.NFWProfile(**kwargs)

The Navarro-Frenk-White profile.

The constructor accepts either the free parameters in this formula, central density and scale radius, or a spherical overdensity mass and concentration (in this case the mass definition and redshift also need to be specified, and a cosmology needs to be set). The density and other commonly used routines are implemented both as class and as static routines, meaning they can be called without instantiating the class.

Parameters
rhos: float

The central density in physical \(M_{\odot} h^2 / {\rm kpc}^3\).

rs: float

The scale radius in physical kpc/h.

M: float

A spherical overdensity mass in \(M_{\odot}/h\) corresponding to the mass definition mdef at redshift z.

c: float

The concentration, \(c = R / r_{\rm s}\), corresponding to the given halo mass and mass definition.

z: float

Redshift

mdef: str

The mass definition in which M and c are given. See Halo Mass Definitions for details.

Methods

M(rhos, rs, x)

The enclosed mass in an NFW profile as a function of \(x=r/r_{\rm s}\).

M4rs()

The mass within 4 scale radii, \(M_{<4rs}\).

MDelta(z, mdef)

The spherical overdensity mass of a given mass definition.

RDelta(z, mdef)

The spherical overdensity radius of a given mass definition.

RMDelta(z, mdef)

The spherical overdensity radius and mass of a given mass definition.

Rsteepest([search_range])

The radius where the logarithmic slope of the density profile is steepest.

Vmax()

The maximum circular velocity, and the radius where it occurs.

circularVelocity(r)

The circular velocity, \(v_c \equiv \sqrt{GM(<r)/r}\).

cumulativePdf(r[, Rmax, z, mdef])

The cumulative distribution function of the profile.

deltaSigma(r[, interpolate, ...])

The excess surface density at radius r.

deltaSigmaInner(r[, interpolate, ...])

The excess surface density at radius r due to the inner profile.

deltaSigmaOuter(r[, interpolate, ...])

The excess surface density at radius r due to the outer profile.

density(r)

Density as a function of radius.

densityDerivativeLin(r)

The linear derivative of density, \(d \rho / dr\).

densityDerivativeLinInner(r)

The linear derivative of the inner density, \(d \rho_{\rm inner} / dr\).

densityDerivativeLinOuter(r)

The linear derivative of the outer density, \(d \rho_{\rm outer} / dr\).

densityDerivativeLog(r)

The logarithmic derivative of density, \(d \log(\rho) / d \log(r)\).

densityDerivativeLogInner(r)

The logarithmic derivative of the inner density, \(d \log(\rho_{\rm inner}) / d \log(r)\).

densityDerivativeLogOuter(r)

The logarithmic derivative of the outer density, \(d \log(\rho_{\rm outer}) / d \log(r)\).

densityInner(r)

Density of the inner profile as a function of radius.

densityOuter(r)

Density of the outer profile as a function of radius.

enclosedMass(r[, accuracy])

The mass enclosed within radius r.

enclosedMassInner(r[, accuracy])

The mass enclosed within radius r due to the inner profile term.

enclosedMassOuter(r[, accuracy])

The mass enclosed within radius r due to the outer profile term.

fit(r, q, quantity[, q_err, q_cov, method, ...])

Fit the density, mass, surface density, or DeltaSigma profile to a given set of data points.

getParameterArray([mask])

Returns an array of the profile parameters.

mu(x)

A function of \(x=r/r_{\rm s}\) that appears in the NFW enclosed mass.

nativeParameters(M, c, z, mdef)

The native NFW parameters, \(\rho_s\) and \(r_{\rm s}\), from mass and concentration.

rho(rhos, x)

The NFW density as a function of \(x = r/r_{\rm s}\).

setNativeParameters(M, c, z, mdef, **kwargs)

Set the native NFW parameters from mass and concentration.

setParameterArray(pars[, mask])

Set the profile parameters from an array.

surfaceDensity(r[, interpolate, accuracy, ...])

The projected surface density at radius r.

surfaceDensityInner(r, **kwargs)

The projected surface density at radius r due to the inner profile.

surfaceDensityOuter(r[, interpolate, ...])

The projected surface density at radius r due to the outer profile.

update()

Update the profile object after a change in parameters or cosmology.

updateR200m()

Update the internally stored R200m after a parameter change.

xDelta(rhos, density_threshold)

Find \(x=r/r_{\rm s}\) where the enclosed density has a particular value.

fundamentalParameters

classmethod nativeParameters(M, c, z, mdef)

The native NFW parameters, \(\rho_s\) and \(r_{\rm s}\), from mass and concentration.

This routine is called in the constructor of the NFW profile class (unless \(\rho_s\) and \(r_{\rm s}\) are passed by the user), but can also be called without instantiating an NFWProfile object.

Parameters
M: array_like

Spherical overdensity mass in \(M_{\odot}/h\); can be a number or a numpy array.

c: array_like

The concentration, \(c = R / r_{\rm s}\), corresponding to the given halo mass and mass definition; must have the same dimensions as M.

z: float

Redshift

mdef: str

The mass definition in which M and c are given. See Halo Mass Definitions for details.

Returns
rhos: array_like

The central density in physical \(M_{\odot} h^2 / {\rm kpc}^3\); has the same dimensions as M.

rs: array_like

The scale radius in physical kpc/h; has the same dimensions as M.

setNativeParameters(M, c, z, mdef, **kwargs)

Set the native NFW parameters from mass and concentration.

The NFW profile has \(\rho_s\) and \(r_{\rm s}\) as internal parameters, which are computed from a mass and concentration. This function ignores the presence of outer profiles.

Parameters
M: float

Spherical overdensity mass in \(M_{\odot}/h\).

c: float

The concentration, \(c = R / r_{\rm s}\), corresponding to the given halo mass and mass definition.

z: float

Redshift

mdef: str

The mass definition in which M and c are given. See Halo Mass Definitions for details.

static rho(rhos, x)

The NFW density as a function of \(x = r/r_{\rm s}\).

This routine can be called without instantiating an NFWProfile object. In most cases, the density() function should be used instead.

Parameters
rhos: float

The central density in physical \(M_{\odot} h^2 / {\rm kpc}^3\).

x: array_like

The radius in units of the scale radius, \(x=r/r_{\rm s}\); can be a number or a numpy array.

Returns
rho: array_like

Density in physical \(M_{\odot} h^2 / {\rm kpc}^3\); has the same dimensions as x.

See also

halo.profile_base.HaloDensityProfile.density

Density as a function of radius.

static mu(x)

A function of \(x=r/r_{\rm s}\) that appears in the NFW enclosed mass.

This routine can be called without instantiating an NFWProfile object.

Parameters
x: array_like

The radius in units of the scale radius, \(x=r/r_{\rm s}\); can be a number or a numpy array.

Returns
mu: array_like

Has the same dimensions as x.

See also

M

The enclosed mass in an NFW profile as a function of \(x=r/r_{\rm s}\).

halo.profile_base.HaloDensityProfile.enclosedMass

The mass enclosed within radius r.

classmethod M(rhos, rs, x)

The enclosed mass in an NFW profile as a function of \(x=r/r_{\rm s}\).

This routine can be called without instantiating an NFWProfile object. In most cases, the enclosedMass() function should be used instead.

Parameters
rhos: float

The central density in physical \(M_{\odot} h^2 / {\rm kpc}^3\).

rs: float

The scale radius in physical kpc/h.

x: array_like

The radius in units of the scale radius, \(x=r/r_{\rm s}\); can be a number or a numpy array.

Returns
M: array_like

The enclosed mass in \(M_{\odot}/h\); has the same dimensions as x.

See also

mu

A function of \(x=r/r_{\rm s}\) that appears in the NFW enclosed mass.

halo.profile_base.HaloDensityProfile.enclosedMass

The mass enclosed within radius r.

classmethod xDelta(rhos, density_threshold)

Find \(x=r/r_{\rm s}\) where the enclosed density has a particular value.

This function is the basis for the RDelta routine, but can be used without instantiating an NFWProfile object. This is preferable when the function needs to be evaluated many times, for example when converting a large number of mass definitions.

The function uses an interpolation table that makes it orders of magnitude faster than root finding (depending on the size of the rhos and density_threshold arrays).

Parameters
rhos: array_like

The central density in physical \(M_{\odot} h^2 / {\rm kpc}^3\); can be a number or a numpy array.

density_threshold: array_like

The desired enclosed density threshold in physical \(M_{\odot} h^2 / {\rm kpc}^3\). This number can be generated from a mass definition and redshift using the densityThreshold() function. Can be a number or a numpy array, if both density_threshold and rhos are arrays, they must have the same size.

Returns
x: array_like

The radius in units of the scale radius, \(x=r/r_{\rm s}\), where the enclosed density reaches density_threshold. Has the same dimensions as rhos and/or density_threshold.

densityInner(r)

Density of the inner profile as a function of radius.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
density: array_like

Density in physical \(M_{\odot} h^2 / {\rm kpc}^3\); has the same dimensions as r.

densityDerivativeLinInner(r)

The linear derivative of the inner density, \(d \rho_{\rm inner} / dr\).

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
derivative: array_like

The linear derivative in physical \(M_{\odot} h / {\rm kpc}^2\); has the same dimensions as r.

densityDerivativeLogInner(r)

The logarithmic derivative of the inner density, \(d \log(\rho_{\rm inner}) / d \log(r)\).

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
derivative: array_like

The dimensionless logarithmic derivative; has the same dimensions as r.

enclosedMassInner(r, accuracy=None)

The mass enclosed within radius r due to the inner profile term.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
M: array_like

The mass enclosed within radius r, in \(M_{\odot}/h\); has the same dimensions as r.

surfaceDensityInner(r, **kwargs)

The projected surface density at radius r due to the inner profile.

This function uses the analytical formula of Lokas & Mamon 2001 rather than numerical integration.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
Sigma: array_like

The surface density at radius r, in physical \(M_{\odot} h/{\rm kpc}^2\); has the same dimensions as r.

Vmax()

The maximum circular velocity, and the radius where it occurs.

Returns
vmax: float

The maximum circular velocity in km / s.

rmax: float

The radius where fmax occurs, in physical kpc/h.

See also

circularVelocity

The circular velocity, \(v_c \equiv \sqrt{GM(<r)/r}\).

RDelta(z, mdef)

The spherical overdensity radius of a given mass definition.

Parameters
z: float

Redshift

mdef: str

The mass definition for which the spherical overdensity radius is computed. See Halo Mass Definitions for details.

Returns
R: float

Spherical overdensity radius in physical kpc/h.

See also

MDelta

The spherical overdensity mass of a given mass definition.

RMDelta

The spherical overdensity radius and mass of a given mass definition.

M4rs()

The mass within 4 scale radii, \(M_{<4rs}\).

This mass definition was suggested by More et al. 2015, see the Advanced mass definition functions section for details.

Returns
M4rs: float

The mass within 4 scale radii, \(M_{<4rs}\), in \(M_{\odot} / h\).

MDelta(z, mdef)

The spherical overdensity mass of a given mass definition.

Parameters
z: float

Redshift

mdef: str

The mass definition for which the spherical overdensity mass is computed. See Halo Mass Definitions for details.

Returns
M: float

Spherical overdensity mass in \(M_{\odot} /h\).

See also

RDelta

The spherical overdensity radius of a given mass definition.

RMDelta

The spherical overdensity radius and mass of a given mass definition.

RMDelta(z, mdef)

The spherical overdensity radius and mass of a given mass definition.

This is a wrapper for the RDelta() and MDelta() functions which returns both radius and mass.

Parameters
z: float

Redshift

mdef: str

The mass definition for which the spherical overdensity mass is computed. See Halo Mass Definitions for details.

Returns
R: float

Spherical overdensity radius in physical kpc/h.

M: float

Spherical overdensity mass in \(M_{\odot} /h\).

See also

RDelta

The spherical overdensity radius of a given mass definition.

MDelta

The spherical overdensity mass of a given mass definition.

Rsteepest(search_range=10.0)

The radius where the logarithmic slope of the density profile is steepest.

This function finds the radius where the logarithmic slope of the profile is minimal, within some very generous bounds. The function makes sense only if at least one outer term has been added because the inner profile steepens with radius without ever become shallower again (for any reasonable functional form).

The radius of steepest slope is often taken as a proxy for the splashback radius, \(R_{\rm sp}\), but this correspondence is only approximate because the \(R_{\rm steep}\) is the result of a tradeoff between the orbiting and infalling profiles, whereas the splashback radius is determined by the dynamics of orbiting particles. See the Splashback radius section for a detailed description of the splashback radius.

Parameters
search_range: float

When searching for the radius of steepest slope, search within this factor of \(R_{\rm 200m}\), \(r_{\rm s}\), or another initial guess, which is determined automatically.

Returns
Rsteep: float

The radius where the slope is steepest, in physical kpc/h.

circularVelocity(r)

The circular velocity, \(v_c \equiv \sqrt{GM(<r)/r}\).

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
vc: float

The circular velocity in km / s; has the same dimensions as r.

See also

Vmax

The maximum circular velocity, and the radius where it occurs.

cumulativePdf(r, Rmax=None, z=None, mdef=None)

The cumulative distribution function of the profile.

Some density profiles do not converge to a finite mass at large radius, and the distribution thus needs to be cut off. The user can specify either a radius (in physical kpc/h) where the profile is cut off, or a mass definition and redshift to compute this radius (e.g., the virial radius \(R_{vir}\) at z = 0).

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Rmax: float

The radius where to cut off the profile in physical kpc/h.

z: float

Redshift

mdef: str

The radius definition for the cut-off radius. See Halo Mass Definitions for details.

Returns
pdf: array_like

The probability for mass to lie within radius r; has the same dimensions as r.

deltaSigma(r, interpolate=True, interpolate_surface_density=True, accuracy=0.0001, min_r_interpolate=1e-06, max_r_interpolate=100000000.0, max_r_integrate=1e+20)

The excess surface density at radius r.

This quantity is useful in weak lensing studies, and is defined as \(\Delta\Sigma(R) = \Sigma(<R)-\Sigma(R)\) where \(\Sigma(<R)\) is the averaged surface density within R weighted by area,

\[\Delta\Sigma(R) = \frac{1}{\pi R^2} \int_0^{R} 2 \pi r \Sigma(r) dr - \Sigma(R)\]
Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

interpolate: bool

Use an interpolation table for the surface density during the integration. This can speed up the evaluation significantly, as the surface density can be expensive to evaluate.

interpolate_surface_density: bool

Use an interpolation table for density during the computation of the surface density. This should make the evaluation somewhat faster, but can fail for some density terms which are negative at particular radii.

accuracy: float

The minimum accuracy of the integration (used both to compute the surface density and average it to get DeltaSigma).

min_r_interpolate: float

The minimum radius in physical kpc/h from which the surface density profile is averaged.

max_r_interpolate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using interpolating density.

max_r_integrate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using exact densities.

Returns
DeltaSigma: array_like

The excess surface density at radius r, in physical \(M_{\odot} h/{\rm kpc}^2\); has the same dimensions as r.

deltaSigmaInner(r, interpolate=True, interpolate_surface_density=True, accuracy=0.0001, min_r_interpolate=1e-06, max_r_interpolate=100000000.0, max_r_integrate=1e+20)

The excess surface density at radius r due to the inner profile.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

interpolate: bool

Use an interpolation table for the surface density during the integration. This can speed up the evaluation significantly, as the surface density can be expensive to evaluate.

interpolate_surface_density: bool

Use an interpolation table for density during the computation of the surface density. This should make the evaluation somewhat faster, but can fail for some density terms which are negative at particular radii.

accuracy: float

The minimum accuracy of the integration (used both to compute the surface density and average it to get DeltaSigma).

min_r_interpolate: float

The minimum radius in physical kpc/h from which the surface density profile is averaged.

max_r_interpolate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using interpolating density.

max_r_integrate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using exact densities.

Returns
DeltaSigma: array_like

The excess surface density at radius r, in physical \(M_{\odot} h/{\rm kpc}^2\); has the same dimensions as r.

deltaSigmaOuter(r, interpolate=True, interpolate_surface_density=True, accuracy=0.0001, min_r_interpolate=1e-06, max_r_interpolate=100000000.0, max_r_integrate=1e+20)

The excess surface density at radius r due to the outer profile.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

interpolate: bool

Use an interpolation table for the surface density during the integration. This can speed up the evaluation significantly, as the surface density can be expensive to evaluate.

interpolate_surface_density: bool

Use an interpolation table for density during the computation of the surface density. This should make the evaluation somewhat faster, but can fail for some density terms which are negative at particular radii.

accuracy: float

The minimum accuracy of the integration (used both to compute the surface density and average it to get DeltaSigma).

min_r_interpolate: float

The minimum radius in physical kpc/h from which the surface density profile is averaged.

max_r_interpolate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using interpolating density.

max_r_integrate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using exact densities.

Returns
DeltaSigma: array_like

The excess surface density at radius r, in physical \(M_{\odot} h/{\rm kpc}^2\); has the same dimensions as r.

density(r)

Density as a function of radius.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
density: array_like

Density in physical \(M_{\odot} h^2 / {\rm kpc}^3\); has the same dimensions as r.

densityDerivativeLin(r)

The linear derivative of density, \(d \rho / dr\).

This function should generally not be overwritten by child classes since it handles the general case of adding up the contributions from the inner and outer terms.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
derivative: array_like

The linear derivative in physical \(M_{\odot} h / {\rm kpc}^2\); has the same dimensions as r.

densityDerivativeLinOuter(r)

The linear derivative of the outer density, \(d \rho_{\rm outer} / dr\).

This function should generally not be overwritten by child classes since it handles the general case of adding up the contributions from all outer profile terms.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
derivative: array_like

The linear derivative in physical \(M_{\odot} h / {\rm kpc}^2\); has the same dimensions as r.

densityDerivativeLog(r)

The logarithmic derivative of density, \(d \log(\rho) / d \log(r)\).

This function should generally not be overwritten by child classes since it handles the general case of adding up the contributions from the inner and outer profile terms.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
derivative: array_like

The dimensionless logarithmic derivative; has the same dimensions as r.

densityDerivativeLogOuter(r)

The logarithmic derivative of the outer density, \(d \log(\rho_{\rm outer}) / d \log(r)\).

This function should generally not be overwritten by child classes since it handles the general case of adding up the contributions from outer profile terms.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
derivative: array_like

The dimensionless logarithmic derivative; has the same dimensions as r.

densityOuter(r)

Density of the outer profile as a function of radius.

This function should generally not be overwritten by child classes since it handles the general case of adding up the contributions from all outer profile terms.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

Returns
density: array_like

Density in physical \(M_{\odot} h^2 / {\rm kpc}^3\); has the same dimensions as r.

enclosedMass(r, accuracy=1e-06)

The mass enclosed within radius r.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

accuracy: float

The minimum accuracy of the integration.

Returns
M: array_like

The mass enclosed within radius r, in \(M_{\odot}/h\); has the same dimensions as r.

enclosedMassOuter(r, accuracy=1e-06)

The mass enclosed within radius r due to the outer profile term.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

accuracy: float

The minimum accuracy of the integration.

Returns
M: array_like

The mass enclosed within radius r, in \(M_{\odot}/h\); has the same dimensions as r.

fit(r, q, quantity, q_err=None, q_cov=None, method='leastsq', mask=None, verbose=True, tolerance=1e-05, maxfev=None, leastsq_algorithm='trf', use_legacy_leastsq=False, bounds=None, initial_step=0.01, nwalkers=100, random_seed=None, convergence_step=100, converged_GR=0.01, best_fit='median', output_every_n=100)

Fit the density, mass, surface density, or DeltaSigma profile to a given set of data points.

This function represents a general interface for finding the best-fit parameters of a halo density profile given a set of data points. These points can represent a number of different physical quantities: quantity can either be density, enclosed mass, surface density, or Delta Sigma (rho, M, Sigma, or DeltaSigma).

The data points q at radii r can optionally have error bars, and the user can pass a full covariance matrix. Please note that not passing any estimate of the uncertainty, i.e. q_err = None and q_cov = None, can lead to very poor fit results: the fitter will minimize the absolute difference between points, strongly favoring the high densities at the center.

The user can choose to vary only a sub-set of the profile parameters through the mask parameter. The current parameters of the profile instance serve as an initial guess.

By default, the parameters are transformed into log space during fitting to ensure positivity. Child classes can change this behavior by overwriting the _fitConvertParams() and _fitConvertParamsBack() functions.

There are two fundamental methods for performing the fit, a least-squares minimization (method = 'leastsq') and a Markov-Chain Monte Carlo (method = 'mcmc'). Both variants obey a few specific options (see below). The function returns a dictionary with outputs that somewhat depend on which method is chosen. After the function has completed, the profile instance represents the best-fit profile to the data points (i.e., its parameters are the best-fit parameters). Note that all output parameters are bundled into one dictionary. The explanations below refer to the entries in this dictionary.

Parameters
r: array_like

The radii of the data points, in physical kpc/h.

q: array_like

The data to fit; can either be density in physical \(M_{\odot} h^2 / {\rm kpc}^3\), enclosed mass in \(M_{\odot} /h\), or surface density in physical \(M_{\odot} h/{\rm kpc}^2\). Must have the same dimensions as r.

quantity: str

Indicates which quantity is given in as input in q, can be rho, M, Sigma, or DeltaSigma.

q_err: array_like

Optional; the uncertainty on the values in q in the same units. If method == 'mcmc', either q_err or q_cov must be passed. If method == 'leastsq' and neither q_err nor q_cov are passed, the absolute different between data points and fit is minimized. In this case, the returned chi2 is in units of absolute difference, meaning its value will depend on the units of q.

q_cov: array_like

Optional; the covariance matrix of the elements in q, as a 2-dimensional numpy array. This array must have dimensions of q**2 and be in units of the square of the units of q. If q_cov is passed, q_err is ignored since the diagonal elements of q_cov correspond to q_err**2.

method: str

The fitting method; can be leastsq for a least-squares minimization of mcmc for a Markov-Chain Monte Carlo.

mask: array_like

Optional; a numpy array of booleans that has the same length as the variables vector of the density profile class. Only variables where mask == True are varied in the fit, all others are kept constant. Important: this has to be a numpy array rather than a list.

verbose: bool / int

If true, output information about the fitting process. The flag can also be set as a number, where 1 has the same effect as True, and 2 outputs large amounts of information such as the fit parameters at each iteration.

tolerance: float

Only active when method == 'leastsq'. The accuracy to which the best-fit parameters are found.

maxfev: int

Only active when method == 'leastsq'. The maximum number of function evaluations before the fit is aborted. If None, the default value of the scipy least_squares function is used.

leastsq_algorithm: str

Only active when method == 'leastsq'. Can be any of the methods accepted by the scipy least_squares() function. Default is trf (trust region reflective), which works well for profile fits.

use_legacy_leastsq: bool

Only active when method == 'leastsq'. If True, this setting falls back to the old leastsq() in scipy rather than using the newer least_squares(). Should only be used for backward compatibility.

bounds: array_like

Only active when method == 'leastsq' and use_legacy_leastsq == False. If not None, this parameter must be an array of dimensions [2, N_par], giving two sets (lower and upper limits) of the fitted parameters (not the entire parameter array, if a mask is imposed). The limits must be given in linear space. If the parameters are fitted in log space (or some other transformation), that transformation is automatically applied to the bounds. For example, when fitting in log space (the default), the lower bounds must be positive, but the upper bounds can be np.inf.

initial_step: array_like

Only active when method == 'mcmc'. The MCMC samples (“walkers”) are initially distributed in a Gaussian around the initial guess. The width of the Gaussian is given by initial_step, either as an array of length N_par (giving the width of each Gaussian) or as a float number, in which case the width is set to initial_step times the initial value of the parameter.

nwalkers: int

Only active when method == 'mcmc'. The number of MCMC samplers that are run in parallel.

random_seed: int

Only active when method == 'mcmc'. If random_seed is not None, it is used to initialize the random number generator. This can be useful for reproducing results.

convergence_step: int

Only active when method == 'mcmc'. The convergence criteria are computed every convergence_step steps (and output is printed if verbose == True).

converged_GR: float

Only active when method == 'mcmc'. The maximum difference between different chains, according to the Gelman-Rubin criterion. Once the GR indicator is lower than this number in all parameters, the chain is ended. Setting this number too low leads to very long runtimes, but setting it too high can lead to inaccurate results.

best_fit: str

Only active when method == 'mcmc'. This parameter determines whether the mean or median value of the likelihood distribution is used as the output parameter set.

output_every_n: int

Only active when method == 'mcmc'. This parameter determines how frequently the MCMC chain outputs information. Only effective if verbose == True.

Returns
results: dict

A dictionary bundling the various fit results. Regardless of the fitting method, the dictionary always contains the following entries:

x: array_like

The best-fit result vector. If mask is passed, this vector only contains those variables that were varied in the fit.

q_fit: array_like

The fitted profile at the radii of the data points; has the same units as q and the same dimensions as r.

chi2: float

The chi^2 of the best-fit profile. If a covariance matrix was passed, the covariances are taken into account. If no uncertainty was passed at all, chi2 is in units of absolute difference, meaning its value will depend on the units of q.

ndof: int

The number of degrees of freedom, i.e. the number of fitted data points minus the number of free parameters.

chi2_ndof: float

The chi^2 per degree of freedom.

If method == 'leastsq', the dictionary additionally contains the entries returned by scipy.optimize.leastsq as well as the following:

nfev: int

The number of function calls used in the fit.

x_err: array_like

An array of dimensions [2, nparams] which contains an estimate of the lower and upper uncertainties on the fitted parameters. These uncertainties are computed from the covariance matrix estimated by the fitter. Please note that this estimate does not exactly correspond to a 68% likelihood. In order to get more statistically meaningful uncertainties, please use the MCMC samples instead of least-squares. In some cases, the fitter fails to return a covariance matrix, in which case x_err is None.

If method == 'mcmc', the dictionary contains the following entries:

x_initial: array_like

The initial positions of the walkers, in an array of dimensions [nwalkers, nparams].

chain_full: array_like

A numpy array of dimensions [n_independent_samples, nparams] with the parameters at each step in the chain. In this thin chain, only every nth step is output, where n is the auto-correlation time, meaning that the samples in this chain are truly independent.

chain_thin: array_like

Like the thin chain, but including all steps. Thus, the samples in this chain are not indepedent from each other. However, the full chain often gives better plotting results.

R: array_like

A numpy array containing the GR indicator at each step when it was saved.

x_mean: array_like

The mean of the chain for each parameter; has length nparams.

x_median: array_like

The median of the chain for each parameter; has length nparams.

x_stddev: array_like

The standard deviation of the chain for each parameter; has length nparams.

x_percentiles: array_like

The lower and upper values of each parameter that contain a certain percentile of the probability; has dimensions [n_percentages, 2, nparams] where the second dimension contains the lower/upper values.

getParameterArray(mask=None)

Returns an array of the profile parameters.

The profile parameters are internally stored in an ordered dictionary. For some applications (e.g., fitting), a simple array is more appropriate.

Parameters
mask: array_like

Optional; must be a numpy array (not a list) of booleans, with the same length as the parameter vector of the profile class (profile.N_par). Only those parameters that correspond to True values are returned.

Returns
par: array_like

A numpy array with the profile’s parameter values.

setParameterArray(pars, mask=None)

Set the profile parameters from an array.

The profile parameters are internally stored in an ordered dictionary. For some applications (e.g., fitting), setting them directly from an array might be necessary. If the profile contains values that depend on the parameters, the profile class must overwrite this function and update according to the new parameters.

Parameters
pars: array_like

The new parameter array.

mask: array_like

Optional; must be a numpy array (not a list) of booleans, with the same length as the parameter vector of the profile class (profile.N_par). If passed, only those parameters that correspond to True values are set (meaning the pars parameter must be shorter than profile.N_par).

surfaceDensity(r, interpolate=True, accuracy=0.0001, max_r_interpolate=100000000.0, max_r_integrate=1e+20)

The projected surface density at radius r.

The surface density is computed by projecting the 3D density along the line of sight,

\[\Sigma(R) = 2 \int_R^{\infty} \frac{r \rho(r)}{\sqrt{r^2-R^2}} dr\]
Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

interpolate: bool

Use an interpolation table for density during the integration. This should make the evaluation somewhat faster, depending on how large the radius array is.

accuracy: float

The minimum accuracy of the integration.

max_r_interpolate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using interpolating density.

max_r_integrate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using exact densities.

Returns
Sigma: array_like

The surface density at radius r, in physical \(M_{\odot} h/{\rm kpc}^2\); has the same dimensions as r.

surfaceDensityOuter(r, interpolate=True, accuracy=0.0001, max_r_interpolate=100000000.0, max_r_integrate=1e+20)

The projected surface density at radius r due to the outer profile.

This function checks whether there are explicit expressions for the surface density of the outer profile terms available, and uses them if possible. Note that there are some outer terms whose surface density integrates to infinity, such as the mean density of the universe which is constant to infinitely large radii.

Parameters
r: array_like

Radius in physical kpc/h; can be a number or a numpy array.

interpolate: bool

Use an interpolation table for density during the integration. This should make the evaluation somewhat faster, depending on how large the radius array is.

accuracy: float

The minimum accuracy of the integration.

max_r_interpolate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using interpolating density.

max_r_integrate: float

The maximum radius in physical kpc/h to which the density profile is integrated when using exact densities.

Returns
Sigma: array_like

The surface density at radius r, in physical \(M_{\odot} h/{\rm kpc}^2\); has the same dimensions as r.

update()

Update the profile object after a change in parameters or cosmology.

If the parameters dictionary has been changed (e.g. by the user or during fitting), this function must be called to ensure consistency within the profile object. This involves deleting any pre-computed quantities (e.g., tabulated enclosed masses) and re-computing profile properties that depend on the parameters. The corresponding functions for outer terms are automatically called as well.

updateR200m()

Update the internally stored R200m after a parameter change.

If the profile has the internal option opt['R200m'] option, that does not stay in sync with the other profile parameters if they are changed (either inside or outside the constructor). This function adjusts \(R_{\rm 200m}\), in addition to whatever action is taken in the update function of the super class. Note that this adjustment needs to be done iteratively if any outer profiles rely on \(R_{\rm 200m}\).

halo.profile_nfw.radiusFromPdf(M, c, z, mdef, cumulativePdf, interpolate=True, min_interpolate_pdf=0.01)

Get the radius where the cumulative density distribution of a halo has a certain value, assuming an NFW profile.

This function can be useful when assigning radii to satellite galaxies in mock halos, for example. The function is optimized for speed when M is a large array. The density distribution is cut off at the virial radius corresponding to the given mass definition. For example, if mdef == 'vir', the NFW profile is cut off at \(R_{\rm vir}\). The accuracy achieved is about 0.2%, unless min_interpolate_pdf is changed to a lower value; below 0.01, the accuracy of the interpolation drops.

Parameters
M: array_like

Halo mass in units of \(M_{\odot}/h\); can be a number or a numpy array.

c: array_like

Halo concentration, in the same definition as M; must have the same dimensions as M.

z: float

Redshift

mdef: str

The mass definition in which the halo mass M is given. See Halo Mass Definitions for details.

cumulativePdf: array_like

The cumulative pdf that we are seeking. If an array, this array needs to have the same dimensions as the M array.

c_model: str

The model used to evaluate concentration if c == None.

interpolate: bool

If interpolate == True, an interpolation table is built before computing the radii. This is much faster if M is a large array.

min_interpolate_pdf: float

For values of the cumulativePdf that fall below this value, the radius is computed exactly, even if interpolation == True. The reason is that the interpolation becomes unreliable for these very low pdfs.

Returns
r: array_like

The radii where the cumulative pdf(s) is/are achieved, in units of physical kpc/h; has the same dimensions as M.

Warning

If many pdf values fall below min_interpolate_pdf, this will slow the function down significantly.