Galaxy-galaxy lensing
Excess surface mass density ΔΣ(r_p) for galaxy-galaxy lensing.
Provides:
Production measurement via the dsigma pipeline:
delta_sigma()
The delta_sigma() function follows the same steps used for DES Y3,
HSC, and KiDS DR4 measurements:
Lens-source pair precomputation with optional photo-z n(z) calibration.
Random-point subtraction to suppress additive systematics.
Adaptive jackknife resampling for covariance estimation.
Stacking with survey-specific shear corrections.
Use the survey loaders in sum_stat.lensing.surveys to prepare the
source and n(z) tables, then pass them directly to delta_sigma().
- sum_stat.lensing.esd.delta_sigma(lens, table_s, rp_bins, cosmo, *, rand=None, table_n=None, lens_source_cut=0.1, n_jackknife=None, comoving=True, corrections=None)[source]
Excess surface mass density ΔΣ(r_p) via the dsigma pipeline.
Follows the same steps as the DES Y3, HSC, and KiDS DR4 measurement scripts: precompute lens-source pair sums, subtract randoms, assign jackknife regions adaptively, and stack with survey-specific corrections.
- Parameters:
lens (GalaxyCatalogue) – Lens catalogue (ra, dec, redshift, weight required).
table_s (Table) – dsigma-formatted source shape table. Prepare with the survey-specific loaders in
sum_stat.lensing.surveys(e.g.load_des_sources()).rp_bins (ndarray, shape (n_bins+1,)) – Projected separation bin edges [Mpc].
cosmo (FlatLambdaCDM) – Cosmology for Σ_crit computation.
rand (GalaxyCatalogue, optional) – Random-point catalogue used for additive-bias subtraction. Required when
correctionscontainsrandom_subtraction=True.table_n (Table, optional) – Photo-z calibration n(z) table (columns
zandn). Used by dsigma to correct for photo-z dilution whenlens_source_cutis applied. Returned by the survey loaders alongsidetable_s.lens_source_cut (float) – Minimum redshift separation Δz = z_s − z_l required to include a lens-source pair. Typical values: 0.1 (DES, KiDS), 0.3 (HSC).
n_jackknife (int, optional) – Number of jackknife regions. When None (default), the count is chosen adaptively: 10 for N_lens < 100, 50 for N_lens < 500, 100 otherwise — matching the reference scripts.
comoving (bool) – Use comoving coordinates for Σ_crit (default True).
corrections (dict, optional) – Boolean correction flags forwarded to
dsigma.stacking.excess_surface_density. Survey-specific defaults:DES Y3
{'matrix_shear_response_correction': True, 'random_subtraction': True}
HSC
{'scalar_shear_response_correction': True, 'shear_responsivity_correction': True, 'selection_bias_correction': True, 'boost_correction': False, 'random_subtraction': True}
KiDS DR4
{'scalar_shear_response_correction': True, 'random_subtraction': True}
- Returns:
rp_centres (ndarray, shape (n_bins,)) – Projected separation bin centres [Mpc] (geometric mean of bin edges).
ds (ndarray, shape (n_bins,)) – ΔΣ(r_p) [M_sun/pc^2].
cov (ndarray, shape (n_bins, n_bins)) – Jackknife covariance [(M_sun/pc^2)^2].
- Raises:
ValueError – If
correctionsrequests random subtraction butrandis None.ValueError – If fewer than 10 lenses remain after the precomputation filter.
Performance –
----------- –
~30 s/call (N_lens=1000, N_src=3000, n_bins=10, n_jk=100, CPU x86-64) –
- Parameters:
lens (GalaxyCatalogue)
table_s (Table)
rp_bins (ndarray)
cosmo (FlatLambdaCDM)
rand (GalaxyCatalogue | None)
table_n (Table | None)
lens_source_cut (float)
n_jackknife (int | None)
comoving (bool)
corrections (dict | None)
- Return type:
Examples
DES Y3 measurement:
>>> from sum_stat.lensing.surveys import load_des_sources >>> from sum_stat import GalaxyCatalogue >>> from astropy.cosmology import Planck15 >>> import numpy as np >>> >>> table_s, table_n = load_des_sources('des_y3.hdf5') >>> lens = GalaxyCatalogue(ra=..., dec=..., redshift=...) >>> rand = GalaxyCatalogue(ra=..., dec=..., redshift=...) >>> rp_bins = 10 ** np.arange(-2.0, 1.6, 0.4) >>> rp, ds, cov = delta_sigma( ... lens, table_s, rp_bins, Planck15, ... rand=rand, table_n=table_n, ... corrections={'matrix_shear_response_correction': True, ... 'random_subtraction': True}, ... )
Shear calibration corrections for galaxy-galaxy lensing.
All correction functions are JAX-differentiable where inputs are arrays.
- sum_stat.lensing.shear_calib.apply_multiplicative_bias(delta_sigma, mean_m)[source]
Correct ΔΣ for multiplicative shear bias.
ΔΣ_corr = ΔΣ / (1 + m̄)
- Parameters:
delta_sigma (jnp.ndarray, shape (n_bins,)) – Raw excess surface density [M_sun/pc^2].
mean_m (float) – Weighted mean multiplicative bias ⟨m⟩.
- Returns:
delta_sigma_corr (jnp.ndarray, shape (n_bins,)) – Corrected ΔΣ [M_sun/pc^2].
- Parameters:
delta_sigma (Array)
mean_m (float)
- Return type:
Array
- sum_stat.lensing.shear_calib.apply_shear_response(e1, e2, R11, R22)[source]
Apply diagonal shear response correction to ellipticity components.
γ_i = e_i / R_ii
- Parameters:
e1, e2 (jnp.ndarray, shape (N,)) – Measured ellipticity components.
R11, R22 (jnp.ndarray, shape (N,)) – Diagonal elements of the per-galaxy shear response matrix.
- Returns:
gamma1, gamma2 (jnp.ndarray, shape (N,)) – Shear estimates corrected for response.
- Parameters:
e1 (Array)
e2 (Array)
R11 (Array)
R22 (Array)
- Return type:
tuple[Array, Array]
- sum_stat.lensing.shear_calib.effective_sigma_crit_inv(z_l, z_s_grid, n_z, h, omega_m, comoving=True)[source]
Effective inverse critical surface density ⟨Σ_crit^{-1}⟩.
⟨Σ_crit^{-1}⟩ = ∫ Σ_crit^{-1}(z_l, z_s) n(z_s) dz_s
- Parameters:
z_l (float) – Lens redshift.
z_s_grid (ndarray, shape (n_z,)) – Source redshift grid.
n_z (ndarray, shape (n_z,)) – Source n(z) (normalised, units sr^-1 dz^-1).
h, omega_m (float) – Cosmological parameters for distance computation.
comoving (bool) – Use comoving Σ_crit if True.
- Returns:
sigma_crit_inv_eff (float) – ⟨Σ_crit^{-1}⟩ [pc^2/M_sun].
- Parameters:
- Return type:
- sum_stat.lensing.shear_calib.nz_from_calibration_table(calib, z_bins, lens_z=None)[source]
Estimate the source n(z) from a photo-z calibration table.
Optionally selects only source galaxies behind a given lens redshift (z_true > lens_z) to avoid source contamination.
- Parameters:
calib (PhotoZCalibTable) – Photo-z calibration table.
z_bins (ndarray, shape (n_bins+1,)) – Redshift bin edges for the output n(z).
lens_z (float, optional) – If provided, only include calibration objects with z_true > lens_z.
- Returns:
z_centres (ndarray, shape (n_bins,)) – Redshift bin centres.
n_z (ndarray, shape (n_bins,)) – Normalised n(z) (integrates to 1 over z_bins).
- Parameters:
calib (PhotoZCalibTable)
z_bins (ndarray)
lens_z (float | None)
- Return type: