clmm.dataops package

Functions to compute polar/azimuthal averages in radial bins

clmm.dataops.compute_background_probability(z_lens, z_src=None, use_pdz=False, pzpdf=None, pzbins=None, validate_input=True)[source]

Probability for being a background galaxy, defined by:

\[P(z_s > z_l) = \int_{z_l}^{+\infty} dz_s \; p_{\text{photoz}}(z_s),\]

when the photometric probability density functions (\(p_{\text{photoz}}(z_s)\)) are provided. In the case of true redshifts, it returns 1 if \(z_s > z_l\) else returns 0.

Parameters:
  • z_lens (float) -- Redshift of the lens.

  • z_src (array, optional) -- Redshift of the source. Used only if pzpdf=pzbins=None.

Returns:

p_background -- Probability for being a background galaxy

Return type:

array

clmm.dataops.compute_galaxy_weights(use_shape_noise=False, shape_component1=None, shape_component2=None, use_shape_error=False, shape_component1_err=None, shape_component2_err=None, is_deltasigma=False, sigma_c=None, validate_input=True)[source]

Computes the individual lens-source pair weights

The weights \(w_{ls}\) express as : \(w_{ls} = w_{ls, \text{geo}} \times w_{ls, \text{shape}}\), following E. S. Sheldon et al. (2003), arXiv:astro-ph/0312036:

1. If computed for shear, the geometrical weights \(w_{ls, \text{geo}}\) are equal to 1. If computed for \(\Delta \Sigma\), it depends on lens and source redshift information via the critical surface density. This component can be expressed as:

\[w_{ls, \text{geo}} = \Sigma_\text{crit}(z_l, z_{\text{src}})^{-2}\;.\]

when only redshift point estimates are provided, or as:

\[w_{ls, \text{geo}} = \Sigma_\text{crit}^\text{eff}(z_l, z_{\text{src}})^{-2} = \left[\int_{\delta + z_l} dz_s \; p_{\text{photoz}}(z_s) \Sigma_\text{crit}(z_l, z_s)^{-1}\right]^2\]

when the redshift pdf of each source, \(p_{\text{photoz}}(z_s)\), is known.

2. The shape weight \(w_{ls,{\text{shape}}}\) depends on shapenoise and/or shape measurement errors

\[w_{ls, \text{shape}}^{-1} = \sigma_{\text{shapenoise}}^2 + \sigma_{\text{measurement}}^2\]
Parameters:
  • is_deltasigma (bool) -- If False, weights are computed for shear, else weights are computed for \(\Delta \Sigma\).

  • sigma_c (None, array_like) -- Critical (effective) surface density in units of \(M_\odot\ Mpc^{-2}\). Used only when is_deltasigma=True.

  • use_shape_noise (bool) -- If True shape noise is included in the weight computation. It then requires shape_componenet{1,2} to be provided. Default: False.

  • shape_component1 (array_like) -- The measured shear (or reduced shear or ellipticity) of the source galaxies, used if use_shapenoise=True

  • shape_component2 (array_like) -- The measured shear (or reduced shear or ellipticity) of the source galaxies, used if use_shapenoise=True

  • use_shape_error (bool) -- If True shape errors are included in the weight computation. It then requires shape_component{1,2}_err` to be provided. Default: False.

  • shape_component1_err (array_like) -- The measurement error on the 1st-component of ellipticity of the source galaxies, used if use_shape_error=True

  • shape_component2_err (array_like) -- The measurement error on the 2nd-component of ellipticity of the source galaxies, used if use_shape_error=True

  • validate_input (bool) -- Validade each input argument

Returns:

w_ls -- Individual lens source pair weights

Return type:

array_like

clmm.dataops.compute_tangential_and_cross_components(ra_lens, dec_lens, ra_source, dec_source, shear1, shear2, geometry='curve', is_deltasigma=False, sigma_c=None, validate_input=True)[source]

Computes tangential- and cross- components for shear or ellipticity

To do so, we need the right ascension and declination of the lens and of all of the sources. We also need the two shape components of all of the sources.

These quantities can be handed to compute_tangential_and_cross_components in two ways

  1. Pass in each as parameters:

    compute_tangential_and_cross_components(ra_lens, dec_lens, ra_source, dec_source,
    shape_component1, shape_component2)
    
  2. As a method of GalaxyCluster:

    cluster.compute_tangential_and_cross_components()
    

The angular separation between the source and the lens, \(\theta\), and the azimuthal position of the source relative to the lens, \(\phi\), are computed within the function and the angular separation is returned.

In the flat sky approximation, these angles are calculated using (_lens: lens, _source: source, RA is from right to left)

\[\begin{split}\theta^2 = & \left(\delta_s-\delta_l\right)^2+ \left(\alpha_l-\alpha_s\right)^2\cos^2(\delta_l)\\ \tan\phi = & \frac{\delta_s-\delta_l}{\left(\alpha_l-\alpha_s\right)\cos(\delta_l)}\end{split}\]

The tangential, \(g_t\), and cross, \(g_x\), ellipticity/shear components are calculated using the two ellipticity/shear components \(g_1\) and \(g_2\) of the source galaxies, following Eq.7 and Eq.8 in Schrabback et al. (2018), arXiv:1611:03866 which is consistent with arXiv:0509252

\[\begin{split}g_t =& -\left( g_1\cos\left(2\phi\right)+g_2\sin\left(2\phi\right)\right)\\ g_x =& g_1 \sin\left(2\phi\right)-g_2\cos\left(2\phi\right)\end{split}\]

Finally, if the critical surface density (\(\Sigma_\text{crit}\)) is provided, an estimate of the excess surface density \(\widehat{\Delta\Sigma}\) is obtained from

\[\widehat{\Delta\Sigma_{t,x}} = g_{t,x} \times \Sigma_\text{crit}(z_l, z_{\text{src}})\]

If \(g_{t,x}\) correspond to the shear, the above expression is an accurate. However, if \(g_{t,x}\) correspond to ellipticities or reduced shear, this expression only gives an estimate \(\widehat{\Delta\Sigma_{t,x}}\), valid only in the weak lensing regime.

Parameters:
  • ra_lens (float) -- Right ascension of the lensing cluster in degrees

  • dec_lens (float) -- Declination of the lensing cluster in degrees

  • ra_source (array) -- Right ascensions of each source galaxy in degrees

  • dec_source (array) -- Declinations of each source galaxy in degrees

  • shear1 (array) -- The measured shear (or reduced shear or ellipticity) of the source galaxies

  • shear2 (array) -- The measured shear (or reduced shear or ellipticity) of the source galaxies

  • geometry (str, optional) -- Sky geometry to compute angular separation. Options are curve (uses astropy) or flat.

  • is_deltasigma (bool) -- If True, the tangential and cross components returned are multiplied by Sigma_crit. It requires sigma_c argument. Results in units of \(M_\odot\ Mpc^{-2}\)

  • sigma_c (None, array_like) -- Critical (effective) surface density in units of \(M_\odot\ Mpc^{-2}\). Used only when is_deltasigma=True.

  • validate_input (bool) -- Validade each input argument

Returns:

  • angsep (array_like) -- Angular separation between lens and each source galaxy in radians

  • tangential_component (array_like) -- Tangential shear (or assimilated quantity) for each source galaxy

  • cross_component (array_like) -- Cross shear (or assimilated quantity) for each source galaxy

clmm.dataops.make_radial_profile(components, angsep, angsep_units, bin_units, bins=10, components_error=None, error_model='ste', include_empty_bins=False, return_binnumber=False, cosmo=None, z_lens=None, validate_input=True, weights=None)[source]

Compute the angular profile of given components

We assume that the cluster object contains information on the cross and tangential shears or ellipticities and angular separation of the source galaxies

This function can be called in two ways:

  1. Pass explict arguments:

    make_radial_profile([component1, component2], distances, 'radians')
    
  2. Call it as a method of a GalaxyCluster instance and specify bin_units:

    cluster.make_radial_profile('Mpc')

Parameters:
  • components (list of arrays) -- List of arrays to be binned in the radial profile

  • angsep (array) -- Transvesal distances between the sources and the lens

  • angsep_units (str) -- Units of the calculated separation of the source galaxies Allowed Options = ["radians"]

  • bin_units (str) -- Units to use for the radial bins of the radial profile Allowed Options = ["radians", deg", "arcmin", "arcsec", kpc", "Mpc"] (letter case independent)

  • bins (array_like, optional) -- User defined bins to use for the shear profile. If a list is provided, use that as the bin edges. If a scalar is provided, create that many equally spaced bins between the minimum and maximum angular separations in bin_units. If nothing is provided, default to 10 equally spaced bins.

  • components_error (list of arrays, None) -- List of errors for input arrays

  • error_model (str, optional) --

    Statistical error model to use for y uncertainties. (letter case independent)

    ste - Standard error [=std/sqrt(n) in unweighted computation] (Default). std - Standard deviation.

  • include_empty_bins (bool, optional) -- Also include empty bins in the returned table

  • return_binnumber (bool, optional) -- Also returns the indices of the bins for each object

  • cosmo (CLMM.Cosmology) -- CLMM Cosmology object to convert angular separations to physical distances

  • z_lens (array, optional) -- Redshift of the lens

  • validate_input (bool) -- Validade each input argument

  • weights (array-like, optional) -- Array of individual galaxy weights. If specified, the radial binned profile is computed using a weighted average

Returns:

  • profile (GCData) -- Output table containing the radius grid points, the profile of the components p_i, errors p_i_err and number of sources. The errors are defined as the standard errors in each bin.

  • binnumber (1-D ndarray of ints, optional) -- Indices of the bins (corresponding to xbins) in which each value of xvals belongs. Same length as yvals. A binnumber of i means the corresponding value is between (xbins[i-1], xbins[i]).

Notes

This is an example of a place where the cosmology-dependence can be sequestered to another module.

clmm.dataops.make_stacked_radial_profile(angsep, weights, components)[source]

Compute stacked profile, and mean separation distances.

Parameters:
  • angsep (2d array) -- Transvesal distances corresponding to each object with shape n_obj, n_rad_bins.

  • weights (2d array) -- Weights corresponding to each objects with shape n_obj, n_rad_bins.

  • components (list of 2d arrays) -- List of 2d properties of each array to be stacked with shape n_components, n_obj, n_rad_bins.

Returns:

  • staked_angsep (array) -- Mean transversal distance in each radial bin.

  • stacked_components (list of arrays) -- List of stacked components.