clmm.theory.parent_class module

@file parent_class.py CLMModeling abstract class

class clmm.theory.parent_class.CLMModeling(validate_input=True, z_inf=1000)[source]

Bases: object

Object with functions for halo mass modeling

Variables:
  • backend (str) -- Name of the backend being used

  • mdelta (float) -- Mass of the profile, in units of \(M_\odot\)

  • cdelta (float) -- Concentration of the profile

  • massdef (str) -- Profile mass definition ("mean", "critical", "virial" - letter case independent)

  • delta_mdef (int) -- Mass overdensity definition.

  • halo_profile_model (str) -- Profile model parameterization ("nfw", "einasto", "hernquist" - letter case independent)

  • cosmo (Cosmology) -- Cosmology object

  • hdpm (Object) -- Backend object with halo profiles

  • mdef_dict (dict) -- Dictionary with the definitions for mass

  • hdpm_dict (dict) -- Dictionary with the definitions for profile

  • validate_input (bool) -- Validade each input argument

  • cosmo_class (type) -- Type of used cosmology objects

  • z_inf (float) -- The value used as infinite redshift

property cdelta

Concentration of cluster

convert_mass_concentration(z_cl, massdef=None, delta_mdef=None, halo_profile_model=None, alpha=None, verbose=False)[source]

Converts current mass and concentration to the values for a different model.

Parameters:
  • z_cl (float) -- Redshift of the cluster

  • massdef (str, None) -- Profile mass definition to convert to ("mean", "critical", "virial"). If None, same value of current model is used.

  • delta_mdef (int, None) -- Mass overdensity definition to convert to. If None, same value of current model is used.

  • halo_profile_model (str, None) -- Profile model parameterization to convert to ("nfw", "einasto", "hernquist"). If None, same value of current model is used.

  • alpha (float, None) -- Einasto slope to convert to when halo_profile_model='einasto'. If None, same value of current model is used.

Returns:

  • float -- Mass of different model in units of \(M_\odot\).

  • float -- Concentration of different model.

property delta_mdef

Number of deltas in mass definition of cluster

eval_3d_density(r3d, z_cl, verbose=False)[source]

Retrieve the 3d density \(\rho(r)\).

Parameters:
  • r3d (array_like, float) -- Radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

Returns:

3-dimensional mass density in units of \(M_\odot\ Mpc^{-3}\)

Return type:

numpy.ndarray, float

eval_convergence(r_proj, z_cl, z_src, z_src_info='discrete', verbose=False)[source]

Computes the mass convergence

\[\kappa = \frac{\Sigma}{\Sigma_{crit}}\]

or

\[\kappa = \kappa_\infty \times <\beta_s>\]
Parameters:
  • r_proj (array_like) -- The projected radial positions in \(M\!pc\).

  • z_cl (float) -- Galaxy cluster redshift

  • z_src (array_like, float, function) -- Information on the background source galaxy redshift(s). Value required depends on z_src_info (see below).

  • z_src_info (str, optional) --

    Type of redshift information provided by the z_src argument. The following supported options are:

    • 'discrete' (default) : The redshift of sources is provided by z_src. It can be individual redshifts for each source galaxy when z_src is an array or all sources are at the same redshift when z_src is a float.

    • 'beta' : The averaged lensing efficiency is provided by z_src. z_src must be a tuple containing ( \(\langle \beta_s \rangle, \langle \beta_s^2 \rangle\)), the lensing efficiency and square of the lensing efficiency averaged over the galaxy redshift distribution repectively.

      \[\langle \beta_s \rangle = \left\langle \frac{D_{LS}}{D_S}\frac{D_\infty} {D_{L,\infty}}\right\rangle\]
      \[\langle \beta_s^2 \rangle = \left\langle \left(\frac{D_{LS}} {D_S}\frac{D_\infty}{D_{L,\infty}}\right)^2 \right\rangle\]

  • verbose (bool, optional) -- If True, the Einasto slope (alpha_ein) is printed out. Only availble for the NC and CCL backends.

Returns:

Mass convergence, kappa.

Return type:

numpy.ndarray, float

eval_critical_surface_density_eff(z_len, pzbins, pzpdf)[source]

Computes the 'effective critical surface density'

\[\langle \Sigma_{\rm crit}^{-1}\rangle^{-1} = \left(\int \frac{1}{\Sigma_{\rm crit}(z)} p(z) dz\right)^{-1}\]

where \(p(z)\) is the source photoz probability density function. This comes from the maximum likelihood estimator for evaluating a \(\Delta\Sigma\) profile.

For the standard \(\Sigma_{\text{crit}}(z)\) definition, use the eval_sigma_crit method of the CLMM cosmology object.

Parameters:
  • z_len (float) -- Galaxy cluster redshift

  • pzbins (array-like) -- Bins where the source redshift pdf is defined

  • pzpdf (array-like) -- Values of the source redshift pdf

Returns:

sigma_c -- Cosmology-dependent effective critical surface density in units of \(M_\odot\ Mpc^{-2}\)

Return type:

numpy.ndarray, float

eval_excess_surface_density(r_proj, z_cl, verbose=False)[source]

Computes the excess surface density

Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

Returns:

Excess surface density in units of \(M_\odot\ Mpc^{-2}\).

Return type:

numpy.ndarray, float

eval_excess_surface_density_2h(r_proj, z_cl, halobias=1.0, logkbounds=(-5, 5), ksteps=1000, loglbounds=(0, 6), lsteps=500)[source]

Computes the 2-halo term excess surface density (CCL and NC backends only)

Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

  • halobias (float, optional) -- Value of the halo bias

  • logkbounds (tuple(float, float), shape (2,), optional) -- Log10 of the upper and lower bounds for the linear matter power spectrum

  • ksteps (int, optional) -- Number of steps in k-space

  • loglbounds (tuple(float, float), shape (2,), optional) -- Log10 of the upper and lower bounds for numerical integration

  • lsteps (int, optional) -- Number of steps for numerical integration

Returns:

Excess surface density from the 2-halo term in units of \(M_\odot\ Mpc^{-2}\).

Return type:

numpy.ndarray, float

eval_magnification(r_proj, z_cl, z_src, z_src_info='discrete', approx=None, verbose=False, integ_kwargs=None)[source]

Computes the magnification

\[\mu = \frac{1}{(1-\kappa)^2-|\gamma_t|^2}\]
Parameters:
  • r_proj (array_like) -- The projected radial positions in \(M\!pc\).

  • z_cl (float) -- Galaxy cluster redshift

  • z_src (array_like, float, function) -- Information on the background source galaxy redshift(s). Value required depends on z_src_info (see below).

  • z_src_info (str, optional) --

    Type of redshift information provided by the z_src argument. The following supported options are:

    • 'discrete' (default) : The redshift of sources is provided by z_src. It can be individual redshifts for each source galaxy when z_src is an array or all sources are at the same redshift when z_src is a float (Used for approx=None).

    • 'distribution' : A redshift distribution function is provided by z_src. z_src must be a one dimentional function (Used when approx=None).

    • 'beta' : The averaged lensing efficiency is provided by z_src. z_src must be a tuple containing ( \(\langle \beta_s \rangle, \langle \beta_s^2 \rangle\)), the lensing efficiency and square of the lensing efficiency averaged over the galaxy redshift distribution repectively.

      \[\langle \beta_s \rangle = \left\langle \frac{D_{LS}}{D_S}\frac{D_\infty} {D_{L,\infty}}\right\rangle\]
      \[\langle \beta_s^2 \rangle = \left\langle \left(\frac{D_{LS}} {D_S}\frac{D_\infty}{D_{L,\infty}}\right)^2 \right\rangle\]

  • approx (str, optional) --

    Type of computation to be made for magnifications, options are:

    • None (default): Requires z_src_info to be 'discrete' or 'distribution'. If z_src_info='discrete', full computation is made for each r_proj, z_src pair individually. If z_src_info='distribution', magnification at each value of r_proj is calculated as

      \[\mu =\left<\frac{1}{\left(1-\beta_s\kappa_{\infty}\right)^2 -\left(\beta_s\gamma_{\infty}\right)^2}\right> =\frac{\int_{z_{min}}^{z_{max}}\frac{N(z)\text{d}z} {\left(1-\beta_s(z)\kappa_{\infty}\right)^2 -\left(\beta_s(z)\gamma_{\infty}\right)^2}} {\int_{z_{min}}^{z_{max}} N(z)\text{d}z}\]
    • 'order1' : Uses the weak lensing approximation of the magnification with up to first-order terms in \(\kappa_{\infty}\) or \(\gamma_{\infty}\) (z_src_info must be 'beta'):

      \[\mu \approx 1 + 2 \left<\beta_s\right>\kappa_{\infty}\]
    • 'order2' : Uses the weak lensing approximation of the magnification with up to second-order terms in \(\kappa_{\infty}\) or \(\gamma_{\infty}\) (z_src_info must be 'beta'):

      \[\mu \approx 1 + 2 \left<\beta_s\right>\kappa_{\infty} + 3 \left<\beta_s^2\right>\kappa_{\infty}^2 + \left<\beta_s^2\right>\gamma_{\infty}^2\]

  • integ_kwargs (None, dict) --

    Extra arguments for the redshift integration (when approx=None, z_src_info='distribution'). Possible keys are:

    • 'zmin' (None, float) : Minimum redshift to be set as the source of the galaxy when performing the sum. (default=None)

    • 'zmax' (float) : Maximum redshift to be set as the source of the galaxy when performing the sum. (default=10.0)

    • 'delta_z_cut' (float) : Redshift cut so that zmin = z_cl + delta_z_cut. delta_z_cut is ignored if z_min is already provided. (default=0.1)

  • verbose (bool, optional) -- If True, the Einasto slope (alpha_ein) is printed out. Only availble for the NC and CCL backends.

Returns:

mu -- magnification, mu.

Return type:

numpy.ndarray, float

eval_magnification_bias(r_proj, z_cl, z_src, alpha, z_src_info='discrete', approx=None, integ_kwargs=None, verbose=False)[source]

Computes the magnification bias

\[\mu^{\alpha - 1}\]
Parameters:
  • r_proj (array_like) -- The projected radial positions in \(M\!pc\).

  • z_cl (float) -- Galaxy cluster redshift

  • z_src (array_like, float, function) -- Information on the background source galaxy redshift(s). Value required depends on z_src_info (see below).

  • alpha (float) -- Slope of the cummulative number count of background sources at a given magnitude

  • z_src_info (str, optional) --

    Type of redshift information provided by the z_src argument. The following supported options are:

    • 'discrete' (default) : The redshift of sources is provided by z_src. It can be individual redshifts for each source galaxy when z_src is an array or all sources are at the same redshift when z_src is a float (Used for approx=None).

    • 'distribution' : A redshift distribution function is provided by z_src. z_src must be a one dimentional function (Used when approx=None).

    • 'beta' : The averaged lensing efficiency is provided by z_src. z_src must be a tuple containing ( \(\langle \beta_s \rangle, \langle \beta_s^2 \rangle\)), the lensing efficiency and square of the lensing efficiency averaged over the galaxy redshift distribution repectively.

      \[\langle \beta_s \rangle = \left\langle \frac{D_{LS}}{D_S}\frac{D_\infty} {D_{L,\infty}}\right\rangle\]
      \[\langle \beta_s^2 \rangle = \left\langle \left(\frac{D_{LS}} {D_S}\frac{D_\infty}{D_{L,\infty}}\right)^2 \right\rangle\]

  • approx (str, optional) --

    Type of computation to be made for magnification biases, options are:

    • None (default): Requires z_src_info to be 'discrete' or 'distribution'. If z_src_info='discrete', full computation is made for each r_proj, z_src pair individually. If z_src_info='distribution', magnification bias at each value of r_proj is calculated as

      \[\begin{split}\mu^{\alpha-1} &=\left(\left<\frac{1}{\left(1-\beta_s\kappa_{\infty}\right)^2 -\left(\beta_s\gamma_{\infty}\right)^2}\right>\right)^{\alpha-1} \\\\ &=\frac{\int_{z_{min}}^{z_{max}}\frac{N(z)\text{d}z} {\left(\left(1-\beta_s(z)\kappa_{\infty}\right)^2 -\left(\beta_s(z)\gamma_{\infty}\right)^2\right)^{\alpha-1}}} {\int_{z_{min}}^{z_{max}} N(z)\text{d}z}\end{split}\]
    • 'order1' : Uses the weak lensing approximation of the magnification bias with up to first-order terms in \(\kappa_{\infty}\) or \(\gamma_{\infty}\) (z_src_info must be 'beta'):

      \[\mu^{\alpha-1} \approx 1 + \left(\alpha-1\right)\left(2 \left<\beta_s\right>\kappa_{\infty}\right)\]
    • 'order2' : Uses the weak lensing approximation of the magnification bias with up to second-order terms in \(\kappa_{\infty}\) or \(\gamma_{\infty}\) (z_src_info must be 'beta'):

      \[\begin{split}\mu^{\alpha-1} \approx 1 &+ \left(\alpha-1\right)\left(2 \left<\beta_s\right>\kappa_{\infty}\right) \\\\ &+ \left(\alpha-1\right)\left(\left<\beta_s^2\right>\gamma_{\infty}^2\right) \\\\ &+ \left(2\alpha-1\right)\left(\alpha-1\right) \left(\left<\beta_s^2\right>\kappa_{\infty}^2\right)\end{split}\]

  • integ_kwargs (None, dict) --

    Extra arguments for the redshift integration (when approx=None, z_src_info='distribution'). Possible keys are:

    • 'zmin' (None, float) : Minimum redshift to be set as the source of the galaxy when performing the sum. (default=None)

    • 'zmax' (float) : Maximum redshift to be set as the source of the galaxy when performing the sum. (default=10.0)

    • 'delta_z_cut' (float) : Redshift cut so that zmin = z_cl + delta_z_cut. delta_z_cut is ignored if z_min is already provided. (default=0.1)

  • verbose (bool, optional) -- If True, the Einasto slope (alpha_ein) is printed out. Only availble for the NC and CCL backends.

Returns:

mu_bias -- magnification bias.

Return type:

numpy.ndarray, float

eval_mass_in_radius(r3d, z_cl, verbose=False)[source]

Computes the mass inside a given radius of the profile. The mass is calculated as

\[M(<\text{r3d}) = M_{\Delta}\; \frac{f\left(\frac{\text{r3d}}{r_{\Delta}/c_{\Delta}}\right)}{f(c_{\Delta})},\]

where \(f(x)\) for the different models are

NFW:

\[\quad \ln(1+x)-\frac{x}{1+x}\]

Einasto: (\(\gamma\) is the lower incomplete gamma function)

\[\gamma(\frac{3}{\alpha}, \frac{2}{\alpha}x^{\alpha})\]

Hernquist:

\[\left(\frac{x}{1+x}\right)^2\]
Parameters:
  • r3d (array_like, float) -- Radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

Returns:

Mass in units of \(M_\odot\)

Return type:

numpy.ndarray, float

eval_mean_surface_density(r_proj, z_cl, verbose=False)[source]

Computes the mean value of surface density inside radius r_proj

Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

Returns:

Excess surface density in units of \(M_\odot\ Mpc^{-2}\).

Return type:

numpy.ndarray, float

eval_rdelta(z_cl)[source]

Retrieves the radius for mdelta

\[r_\Delta=\left(\frac{3 M_\Delta}{4 \pi \Delta \rho_{bkg}(z)}\right)^{1/3}\]
Parameters:

z_cl (float) -- Redshift of the cluster

Returns:

Radius in \(M\!pc\).

Return type:

float

eval_reduced_tangential_shear(r_proj, z_cl, z_src, z_src_info='discrete', approx=None, integ_kwargs=None, verbose=False)[source]

Computes the reduced tangential shear

\[g_t = \frac{\gamma_t}{1-\kappa}\]
Parameters:
  • r_proj (array_like) -- The projected radial positions in \(M\!pc\).

  • z_cl (float) -- Galaxy cluster redshift

  • z_src (array_like, float, function) -- Information on the background source galaxy redshift(s). Value required depends on z_src_info (see below).

  • z_src_info (str, optional) --

    Type of redshift information provided by the z_src argument. The following supported options are:

    • 'discrete' (default) : The redshift of sources is provided by z_src. It can be individual redshifts for each source galaxy when z_src is an array or all sources are at the same redshift when z_src is a float (Used for approx=None).

    • 'distribution' : A redshift distribution function is provided by z_src. z_src must be a one dimentional function (Used when approx=None).

    • 'beta' : The averaged lensing efficiency is provided by z_src. z_src must be a tuple containing ( \(\langle \beta_s \rangle, \langle \beta_s^2 \rangle\)), the lensing efficiency and square of the lensing efficiency averaged over the galaxy redshift distribution repectively.

      \[\langle \beta_s \rangle = \left\langle \frac{D_{LS}}{D_S}\frac{D_\infty} {D_{L,\infty}}\right\rangle\]
      \[\langle \beta_s^2 \rangle = \left\langle \left(\frac{D_{LS}} {D_S}\frac{D_\infty}{D_{L,\infty}}\right)^2 \right\rangle\]

  • approx (str, optional) --

    Type of computation to be made for reduced tangential shears, options are:

    • None (default): Requires z_src_info to be 'discrete' or 'distribution'. If z_src_info='discrete', full computation is made for each r_proj, z_src pair individually. If z_src_info='distribution', reduced tangential shear at each value of r_proj is calculated as

      \[g_t =\left<\frac{\beta_s\gamma_{\infty}}{1-\beta_s\kappa_{\infty}}\right> =\frac{\int_{z_{min}}^{z_{max}}\frac{\beta_s(z)\gamma_{\infty}} {1-\beta_s(z)\kappa_{\infty}}N(z)\text{d}z} {\int_{z_{min}}^{z_{max}} N(z)\text{d}z}\]
    • 'order1' : Same approach as in Weighing the Giants - III (equation 6 in Applegate et al. 2014; https://arxiv.org/abs/1208.0605). z_src_info must be 'beta':

      \[g_t\approx\frac{\left<\beta_s\right>\gamma_{\infty}} {1-\left<\beta_s\right>\kappa_{\infty}}\]
    • 'order2' : Same approach as in Cluster Mass Calibration at High Redshift (equation 12 in Schrabback et al. 2017; https://arxiv.org/abs/1611.03866). z_src_info must be 'beta':

      \[g_t\approx\frac{\left<\beta_s\right>\gamma_{\infty}} {1-\left<\beta_s\right>\kappa_{\infty}} \left(1+\left(\frac{\left<\beta_s^2\right>} {\left<\beta_s\right>^2}-1\right)\left<\beta_s\right>\kappa_{\infty}\right)\]

  • integ_kwargs (None, dict) --

    Extra arguments for the redshift integration (when approx=None, z_src_info='distribution'). Possible keys are:

    • 'zmin' (None, float) : Minimum redshift to be set as the source of the galaxy when performing the sum. (default=None)

    • 'zmax' (float) : Maximum redshift to be set as the source of the galaxy when performing the sum. (default=10.0)

    • 'delta_z_cut' (float) : Redshift cut so that zmin = z_cl + delta_z_cut. delta_z_cut is ignored if z_min is already provided. (default=0.1)

  • verbose (bool, optional) -- If True, the Einasto slope (alpha_ein) is printed out. Only availble for the NC and CCL backends.

Returns:

gt -- Reduced tangential shear

Return type:

numpy.ndarray, float

Notes

Need to figure out if we want to raise exceptions rather than errors here?

eval_surface_density(r_proj, z_cl, verbose=False)[source]

Computes the surface mass density

Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

Returns:

2D projected surface density in units of \(M_\odot\ Mpc^{-2}\)

Return type:

numpy.ndarray, float

eval_surface_density_2h(r_proj, z_cl, halobias=1.0, logkbounds=(-5, 5), ksteps=1000, loglbounds=(0, 6), lsteps=500)[source]

Computes the 2-halo term surface density (CCL and NC backends only)

Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • z_cl (float) -- Redshift of the cluster

  • halobias (float, optional) -- Value of the halo bias

  • logkbounds (tuple(float,float), shape (2,), optional) -- Log10 of the upper and lower bounds for the linear matter power spectrum

  • ksteps (int, optional) -- Number of steps in k-space

  • loglbounds (tuple(float,float), shape (2,), optional) -- Log10 of the upper and lower bounds for numerical integration

  • lsteps (int, optional) -- Number of steps for numerical integration

Returns:

Excess surface density from the 2-halo term in units of \(M_\odot\ Mpc^{-2}\).

Return type:

numpy.ndarray, float

eval_tangential_shear(r_proj, z_cl, z_src, z_src_info='discrete', verbose=False)[source]

Computes the tangential shear

Parameters:
  • r_proj (array_like) -- The projected radial positions in \(M\!pc\).

  • z_cl (float) -- Galaxy cluster redshift

  • z_src (array_like, float, function) -- Information on the background source galaxy redshift(s). Value required depends on z_src_info (see below).

  • z_src_info (str, optional) --

    Type of redshift information provided by the z_src argument. The following supported options are:

    • 'discrete' (default) : The redshift of sources is provided by z_src. It can be individual redshifts for each source galaxy when z_src is an array or all sources are at the same redshift when z_src is a float.

    • 'beta' : The averaged lensing efficiency is provided by z_src. z_src must be a tuple containing ( \(\langle \beta_s \rangle, \langle \beta_s^2 \rangle\)), the lensing efficiency and square of the lensing efficiency averaged over the galaxy redshift distribution repectively.

      \[\langle \beta_s \rangle = \left\langle \frac{D_{LS}}{D_S}\frac{D_\infty} {D_{L,\infty}}\right\rangle\]
      \[\langle \beta_s^2 \rangle = \left\langle \left(\frac{D_{LS}} {D_S}\frac{D_\infty}{D_{L,\infty}}\right)^2 \right\rangle\]

  • verbose (bool, optional) -- If True, the Einasto slope (alpha_ein) is printed out. Only availble for the NC and CCL backends.

Returns:

tangential shear

Return type:

numpy.ndarray, float

get_einasto_alpha(z_cl=None)[source]

Returns the value of the \(\alpha\) parameter for the Einasto profile, if defined

Parameters:

z_cl (float) -- Cluster redshift (required for Einasto with the CCL backend, will be ignored for NC)

property halo_profile_model

Halo profile model

property massdef

Definition for the mass of cluster

property mdelta

Mass of cluster

set_concentration(cdelta)[source]

Sets the concentration

Parameters:

cdelta (float) -- Concentration

Notes

This is equivalent to doing self.cdelta = cdelta

set_cosmo(cosmo)[source]

Sets the cosmology to the internal cosmology object

Parameters:

cosmo (clmm.Comology object, None) -- CLMM Cosmology object. If is None, creates a new instance of self.cosmo_class().

set_einasto_alpha(alpha)[source]

Sets the value of the \(\alpha\) parameter for the Einasto profile

Parameters:

alpha (float, None) -- If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

set_halo_density_profile(halo_profile_model='nfw', massdef='mean', delta_mdef=200)[source]

Sets the definitions for the halo profile

Parameters:
  • halo_profile_model (str) -- Halo mass profile, supported options are 'nfw', 'einasto', 'hernquist' (letter case independent)

  • massdef (str) -- Mass definition, supported options are 'mean', 'critical', 'virial' (letter case independent)

  • delta_mdef (int) -- Overdensity number

set_mass(mdelta)[source]

Sets the value of the \(M_\Delta\).

Parameters:

mdelta (float) -- Galaxy cluster mass \(M_\Delta\) in units of \(M_\odot\)

Notes

This is equivalent to doing self.mdelta = mdelta

set_projected_quad(use_projected_quad)[source]

Control the use of quad_vec to calculate the surface density profile for CCL Einasto profile.

Parameters:

use_projected_quad (bool) -- Only available for Einasto profile with CCL as the backend. If True, CCL will use quad_vec instead of default FFTLog to calculate the surface density profile.