clmm.theory.func_layer module

@file func_layer.py Main functions to encapsule oo calls

clmm.theory.func_layer.compute_3d_density(r3d, mdelta, cdelta, z_cl, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, verbose=False, validate_input=True)[source]

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

Profiles implemented so far are:

nfw: \(\rho(r) = \frac{\rho_0}{\frac{c}{(r/R_{vir})} \left(1+\frac{c}{(r/R_{vir})}\right)^2}\) (Navarro et al. 1996)

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

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster concentration

  • z_cl (float) -- Redshift of the cluster

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition; defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

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

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

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

Return type:

numpy.ndarray, float

Notes

Need to refactor later so we only require arguments that are necessary for all profiles and use another structure to take the arguments necessary for specific models

clmm.theory.func_layer.compute_convergence(r_proj, mdelta, cdelta, z_cluster, z_src, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, z_src_info='discrete', verbose=False, validate_input=True)[source]

Computes the mass convergence

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

or

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

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster NFW concentration.

  • z_cluster (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).

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition. Defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

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

    • 'distribution' : A redshift distribution function is provided by z_src. z_src must be a one dimensional function.

    • '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.

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

kappa -- Mass convergence, kappa.

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_critical_surface_density_eff(cosmo, z_cluster, pzbins, pzpdf, validate_input=True)[source]

Computes the 'effective critical surface density'

\[\langle \Sigma_{\text{crit}}^{-1}\rangle^{-1} = \left(\int \frac{1}{\Sigma_{\text{crit}}(z)}p(z) \mathrm{d}z\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:
  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • z_cluster (float) -- Galaxy cluster redshift

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

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

  • validate_input (bool) -- Validade each input argument

Returns:

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

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_excess_surface_density(r_proj, mdelta, cdelta, z_cl, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, verbose=False, validate_input=True)[source]

Computes the excess surface density

\[\Delta\Sigma(R) = \bar{\Sigma}(<R)-\Sigma(R),\]
Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster concentration

  • z_cl (float) -- Redshift of the cluster

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition; defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

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

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

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

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_excess_surface_density_2h(r_proj, z_cl, cosmo, halobias=1.0, logkbounds=(-5, 5), ksteps=1000, loglbounds=(0, 6), lsteps=500, validate_input=True)[source]

Computes the 2-halo term excess surface density from eq.(13) of Oguri & Hamana (2011)

\[\Delta\Sigma_{\text{2h}}(R) = \frac{\rho_m(z)b(M)}{(1 + z)^3D_A(z)^2} \int\frac{ldl}{(2\pi)} P_{\text{mm}}(k_l, z)J_2(l\theta)\]

where

\[k_l = \frac{l}{D_A(z)(1 +z)}\]

and \(b(M)\) is the halo bias

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

  • z_cl (float) -- Redshift of the cluster

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • 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) -- Steps for the numerical integration

  • validate_input (bool) -- Validade each input argument

Returns:

deltasigma_2h -- 2-halo term excess surface density in units of \(M_\odot\ Mpc^{-2}\).

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_magnification(r_proj, mdelta, cdelta, z_cluster, z_src, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, z_src_info='discrete', approx=None, integ_kwargs=None, verbose=False, validate_input=True)[source]

Computes the magnification

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

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster NFW concentration.

  • z_cluster (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).

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition. Defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

  • 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 dimensional 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.

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

magnification -- Magnification \(\mu\).

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_magnification_bias(r_proj, alpha, mdelta, cdelta, z_cluster, z_src, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, z_src_info='discrete', approx=None, integ_kwargs=None, verbose=False, validate_input=True)[source]

Computes magnification bias from magnification \(\mu\) and slope parameter \(\alpha\) as :

\[\mu^{\alpha - 1}.\]

The alpha parameter depends on the source sample and is computed as the slope of the cummulative numer counts at a given magnitude :

\[\alpha \equiv \alpha(f) = - \frac{\mathrm{d}}{\mathrm{d}\log{f}} \log{n_0(>f)}\]

or,

\[\alpha \equiv \alpha(m) = 2.5 \frac{\mathrm d}{\mathrm d m} \log{n_0(<m)}\]

see e.g. Bartelmann & Schneider 2001; Umetsu 2020

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

  • alpha (array like) -- The slope of the cummulative source number counts.

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster NFW concentration.

  • z_cluster (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).

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition. Defaults to 200.

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • 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 dimensional 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)

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

magnification_bias -- magnification bias

Return type:

numpy.ndarray

clmm.theory.func_layer.compute_mean_surface_density(r_proj, mdelta, cdelta, z_cl, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, verbose=False, validate_input=True)[source]

Computes the mean value of surface density inside radius r_proj

\[\bar{\Sigma}(<R) = \frac{2}{R^2} \int^R_0 dR' R' \Sigma(R'),\]
Parameters:
  • r_proj (array_like) -- Projected radial position from the cluster center in \(M\!pc\).

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster concentration

  • z_cl (float) -- Redshift of the cluster

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition; defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo backend

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

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

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

Return type:

numpy.ndarray, float

Notes

Need to refactory so we only require arguments that are necessary for all models and use another structure to take the arguments necessary for specific models.

clmm.theory.func_layer.compute_reduced_tangential_shear(r_proj, mdelta, cdelta, z_cluster, z_src, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', z_src_info='discrete', approx=None, integ_kwargs=None, alpha_ein=None, validate_input=True, verbose=False)[source]

Computes the reduced tangential shear

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

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster NFW concentration.

  • z_cluster (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).

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition. Defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default);

    • 'critical' - not in cluster_toolkit;

    • 'virial' - not in cluster_toolkit;

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

  • 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 dimensional 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)

  • alpha_ein -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

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

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

gt -- Reduced tangential shear

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_surface_density(r_proj, mdelta, cdelta, z_cl, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, verbose=False, use_projected_quad=False, validate_input=True)[source]

Computes the surface mass density

\[\Sigma(R) = \int^\infty_{-\infty} dx\; \rho \left(\sqrt{R^2+x^2}\right),\]

where \(\rho(r)\) is the 3d density profile.

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

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster concentration

  • z_cl (float) -- Redshift of the cluster

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition; defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

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

  • 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. Default: False

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

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

Return type:

numpy.ndarray, float

Notes

Need to refactory so we only require arguments that are necessary for all models and use another structure to take the arguments necessary for specific models.

clmm.theory.func_layer.compute_surface_density_2h(r_proj, z_cl, cosmo, halobias=1, logkbounds=(-5, 5), ksteps=1000, loglbounds=(0, 6), lsteps=500, validate_input=True)[source]

Computes the 2-halo term surface density from eq.(13) of Oguri & Hamana (2011)

\[\Sigma_{\text{2h}}(R) = \frac{\rho_\text{m}(z)b(M)}{(1 + z)^3D_A(z)^2} \int\frac{ldl}{(2\pi)}P_{\text{mm}}(k_l, z)J_0(l\theta)\]

where

\[k_l = \frac{l}{D_A(z)(1 + z)}\]

and \(b(M)\) is the halo bias

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

  • z_cl (float) -- Redshift of the cluster

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • 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) -- Steps for the numerical integration

  • validate_input (bool) -- Validade each input argument

Returns:

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

Return type:

numpy.ndarray, float

clmm.theory.func_layer.compute_tangential_shear(r_proj, mdelta, cdelta, z_cluster, z_src, cosmo, delta_mdef=200, halo_profile_model='nfw', massdef='mean', alpha_ein=None, z_src_info='discrete', verbose=False, validate_input=True)[source]

Computes the tangential shear

\[\gamma_t = \frac{\Delta\Sigma}{\Sigma_{crit}} = \frac{\bar{\Sigma}-\Sigma}{\Sigma_{crit}}\]

or

\[\gamma_t = \gamma_\infty \times \beta_s\]
Parameters:
  • r_proj (array_like, float) -- The projected radial positions in \(M\!pc\).

  • mdelta (float) -- Galaxy cluster mass in \(M_\odot\).

  • cdelta (float) -- Galaxy cluster NFW concentration.

  • z_cluster (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).

  • cosmo (clmm.cosmology.Cosmology object) -- CLMM Cosmology object

  • delta_mdef (int, optional) -- Mass overdensity definition. Defaults to 200.

  • halo_profile_model (str, optional) --

    Profile model parameterization (letter case independent):

    • 'nfw' (default)

    • 'einasto' - not in cluster_toolkit

    • 'hernquist' - not in cluster_toolkit

  • massdef (str, optional) --

    Profile mass definition, with the following supported options (letter case independent):

    • 'mean' (default)

    • 'critical'

    • 'virial'

  • alpha_ein (float, None, optional) -- If halo_profile_model=='einasto', set the value of the Einasto slope. Option only available for the NumCosmo and CCL backends. If None, use the default value of the backend. (0.25 for the NumCosmo backend and a cosmology-dependent value for the CCL backend.)

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

    • 'distribution' : A redshift distribution function is provided by z_src. z_src must be a one dimensional function.

    • '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.

  • validate_input (bool, optional) -- If True (default), the types of the arguments are checked before proceeding.

Returns:

gammat -- Tangential shear

Return type:

numpy.ndarray, float