ler.lens_galaxy_population.sampler_functions
Module for lens galaxy parameter sampling functions.
This module provides probability density functions (PDFs) and random variable samplers (RVS) for lens galaxy parameters including redshift, velocity dispersion, and axis ratio. It also includes rejection and importance sampling algorithms for sampling lens parameters weighted by gravitational lensing cross sections.
Key Components:
Lens redshift samplers (SIS model from Haris et al. 2018)
Velocity dispersion samplers (generalized gamma distribution)
Axis ratio samplers (Rayleigh and Padilla-Strauss distributions)
Rejection and importance sampling for cross-section weighted parameters
Copyright (C) 2026 Hemantakumar Phurailatpam. Distributed under MIT License.
Module Contents
Functions
Compute lens redshift PDF for SIS model (Haris et al. 2018). |
|
Sample lens redshifts for SIS model (Haris et al. 2018). |
|
|
Calculate the lens galaxy velocity dispersion function at redshift z (Oguri et al. (2018b) + Wempe et al. (2022)). |
|
Calculate the local universe velocity dispersion function. |
|
Compute unnormalized velocity dispersion function using generalized gamma. |
|
Compute normalized velocity dispersion PDF using generalized gamma. |
|
Sample velocity dispersions from generalized gamma distribution. |
|
Sample axis ratios from velocity-dependent Rayleigh distribution. |
|
Compute truncated Rayleigh PDF for axis ratio. |
Sample axis ratios from Padilla & Strauss (2008) distribution. |
|
Compute axis ratio PDF from Padilla & Strauss (2008). |
|
|
Joint rejection sample over (zs, zl, sigma, nuisance). |
|
|
|
Importance/resampling for lens nuisance parameters at fixed (zs, zl). |
|
Full joint importance sample over (zs, zl, sigma, nuisance). |
|
Build a closure that runs one of the cross-section-weighted lens samplers. |
Attributes
- ler.lens_galaxy_population.sampler_functions.lens_redshift_strongly_lensed_sis_analytical_pdf(zl, zs, cosmo=LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0))[source]
Compute lens redshift PDF for SIS model (Haris et al. 2018).
Computes the probability density function for lens redshift between zl=0 and zl=zs using the analytical form from Haris et al. (2018) equation A7, based on the SIS (Singular Isothermal Sphere) lens model.
- Parameters:
- zl
float Redshift of the lens galaxy.
- zs
float Redshift of the source.
- cosmo
astropy.cosmology Cosmology object for distance calculations.
default: LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0)
- zl
- Returns:
- pdf
float Probability density at the given lens redshift.
- pdf
Examples
>>> from astropy.cosmology import LambdaCDM >>> cosmo = LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0) >>> pdf = lens_redshift_strongly_lensed_sis_analytical_pdf(zl=0.5, zs=1.0, cosmo=cosmo) >>> print(f"PDF at zl=0.5: {pdf:.4f}") PDF at zl=0.5: 1.8750
- ler.lens_galaxy_population.sampler_functions.lens_redshift_strongly_lensed_sis_analytical_rvs(size, zs, z_min=0.001, z_max=10.0, cosmo=LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0))[source]
Sample lens redshifts for SIS model (Haris et al. 2018).
Uses inverse transform sampling with the analytical CDF of the Haris et al. (2018) lens redshift distribution to efficiently generate samples.
- Parameters:
- size
int Number of samples to draw.
- zs
float Redshift of the source.
- z_min
float Minimum redshift for interpolation grid.
default: 0.001
- z_max
float Maximum redshift for interpolation grid.
default: 10.0
- cosmo
astropy.cosmology Cosmology object for distance calculations.
default: LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0)
- size
- Returns:
- zl
numpy.ndarray Array of sampled lens redshifts with shape (size,).
- zl
Examples
>>> from astropy.cosmology import LambdaCDM >>> cosmo = LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0) >>> zl_samples = lens_redshift_strongly_lensed_sis_analytical_rvs( ... size=1000, ... zs=1.5, ... z_min=0.001, ... z_max=10.0, ... cosmo=cosmo, ... ) >>> print(f"Mean lens redshift: {zl_samples.mean():.2f}") >>> print(f"Redshift range: [{zl_samples.min():.3f}, {zl_samples.max():.3f}]")
- ler.lens_galaxy_population.sampler_functions.velocity_dispersion_ewoud_denisty_function(sigma, z, alpha=0.94, beta=1.85, phistar=0.02099, sigmastar=113.78)[source]
Calculate the lens galaxy velocity dispersion function at redshift z (Oguri et al. (2018b) + Wempe et al. (2022)).
- Parameters:
- sigma
numpy.ndarray Velocity dispersion of the lens galaxy (km/s).
- z
float Redshift of the lens galaxy.
- alpha
float Shape parameter of the velocity dispersion function.
default: 0.94
- beta
float Slope parameter of the velocity dispersion function.
default: 1.85
- phistar
float Normalization of the velocity dispersion function (Mpc^-3).
default: 2.099e-2
- sigmastar
float Characteristic velocity dispersion (km/s).
default: 113.78
- sigma
- Returns:
- result
numpy.ndarray Velocity dispersion function values.
Negative values are clipped to 0.
- result
Notes
The function name contains the misspelling
denistyfor backward compatibility with existing serialized configuration dictionaries.
- ler.lens_galaxy_population.sampler_functions.velocity_dispersion_bernardi_denisty_function(sigma, alpha, beta, phistar, sigmastar)[source]
Calculate the local universe velocity dispersion function.
This implements the velocity dispersion function from Bernardi et al. (2010).
\[\phi(\sigma) = \phi_* \left(\frac{\sigma}{\sigma_*}\right)^\alpha \exp\left[-\left(\frac{\sigma}{\sigma_*}\right)^\beta\right] \frac{\beta}{\Gamma(\alpha/\beta)\sigma}.\]- Parameters:
- sigma
numpy.ndarray Velocity dispersion of the lens galaxy (km/s).
- alpha
float Shape parameter (alpha/beta is used in the gamma function).
For Oguri et al. (2018b): alpha=0.94
For Choi et al. (2008): alpha=2.32/2.67
- beta
float Slope parameter of the velocity dispersion function.
For Oguri et al. (2018b): beta=1.85
For Choi et al. (2008): beta=2.67
- phistar
float Normalization of the velocity dispersion function (Mpc^-3).
For Oguri et al. (2018b): phistar=2.099e-2*(h/0.7)^3
For Choi et al. (2008): phistar=8.0e-3*h^3
- sigmastar
float Characteristic velocity dispersion (km/s).
For Oguri et al. (2018b): sigmastar=113.78
For Choi et al. (2008): sigmastar=161.0
- sigma
- Returns:
- philoc
numpy.ndarray Local velocity dispersion function values.
- philoc
Notes
The function name contains the misspelling
denistyfor backward compatibility with existing serialized configuration dictionaries.
- ler.lens_galaxy_population.sampler_functions.gengamma_function(sigma, alpha=0.94, beta=1.85, phistar=0.02099, sigmastar=113.78, **kwargs)[source]
Compute unnormalized velocity dispersion function using generalized gamma.
Computes the galaxy velocity dispersion function (VDF) using the generalized gamma distribution formulation from Choi et al. (2007).
- Parameters:
- sigma
floatornumpy.ndarray Velocity dispersion in km/s.
- alpha
float Power-law index governing the low-velocity slope.
default: 0.94
- beta
float Exponential parameter for high-velocity cutoff sharpness.
default: 1.85
- phistar
float Normalization constant (comoving number density, Mpc^-3).
default: 8.0e-3
- sigmastar
float Characteristic velocity scale in km/s.
default: 113.78
- sigma
- Returns:
- density
floatornumpy.ndarray Unnormalized velocity dispersion function value.
- density
Examples
>>> import numpy as np >>> sigma = np.array([150.0, 200.0, 250.0]) >>> density = gengamma_function(sigma) >>> print(f"Density at sigma=200 km/s: {density[1]:.6f}")
- ler.lens_galaxy_population.sampler_functions.gengamma_pdf(sigma, sigma_min=100.0, sigma_max=400.0, alpha=0.94, beta=1.85, sigmastar=113.78)[source]
Compute normalized velocity dispersion PDF using generalized gamma.
Computes the probability density function for velocity dispersion using the generalized gamma distribution, normalized over the specified range.
- Parameters:
- sigma
floatornumpy.ndarray Velocity dispersion in km/s.
- sigma_min
float Minimum velocity dispersion for normalization (km/s).
default: 100.0
- sigma_max
float Maximum velocity dispersion for normalization (km/s).
default: 400.0
- alpha
float Power-law index governing the low-velocity slope.
default: 0.94
- beta
float Exponential parameter for high-velocity cutoff sharpness.
default: 1.85
- sigmastar
float Characteristic velocity scale in km/s.
default: 113.78
- sigma
- Returns:
- pdf
floatornumpy.ndarray Normalized probability density at the given velocity dispersion.
- pdf
Examples
>>> pdf = gengamma_pdf( ... sigma=200.0, ... sigma_min=100.0, ... sigma_max=400.0, ... ) >>> print(f"PDF at sigma=200 km/s: {pdf:.6f}")
- ler.lens_galaxy_population.sampler_functions.gengamma_rvs(size, sigma_min=100.0, sigma_max=400.0, alpha=0.94, beta=1.85, sigmastar=113.78)[source]
Sample velocity dispersions from generalized gamma distribution.
Uses truncated generalized gamma sampling via rejection to generate velocity dispersion samples within the specified bounds.
- Parameters:
- size
int Number of samples to draw.
- sigma_min
float Minimum velocity dispersion (km/s).
default: 100.0
- sigma_max
float Maximum velocity dispersion (km/s).
default: 400.0
- alpha
float Power-law index governing the low-velocity slope.
default: 0.94
- beta
float Exponential parameter for high-velocity cutoff sharpness.
default: 1.85
- sigmastar
float Characteristic velocity scale in km/s.
default: 113.78
- size
- Returns:
- sigma
numpy.ndarray Sampled velocity dispersions in km/s with shape (size,).
- sigma
Examples
>>> sigma_samples = gengamma_rvs( ... size=1000, ... sigma_min=100.0, ... sigma_max=400.0, ... ) >>> print(f"Mean sigma: {sigma_samples.mean():.2f} km/s") >>> print(f"Range: [{sigma_samples.min():.1f}, {sigma_samples.max():.1f}] km/s")
- ler.lens_galaxy_population.sampler_functions.rayleigh_rvs(size, sigma, q_min=0.2, q_max=1.0)[source]
Sample axis ratios from velocity-dependent Rayleigh distribution.
Generates axis ratio samples using the Rayleigh distribution with scale parameter dependent on velocity dispersion, as described in Wierda et al. (2021) Appendix C.
- Parameters:
- size
int Number of samples to draw.
- sigma
numpy.ndarray Velocity dispersions in km/s with shape (size,).
- q_min
float Minimum allowed axis ratio.
default: 0.2
- q_max
float Maximum allowed axis ratio.
default: 1.0
- size
- Returns:
- q
numpy.ndarray Sampled axis ratios with shape (size,).
- q
Examples
>>> import numpy as np >>> sigma = np.random.uniform(100, 300, 1000) >>> q_samples = rayleigh_rvs(size=1000, sigma=sigma, q_min=0.2, q_max=1.0) >>> print(f"Mean axis ratio: {q_samples.mean():.2f}") >>> print(f"Range: [{q_samples.min():.2f}, {q_samples.max():.2f}]")
- ler.lens_galaxy_population.sampler_functions.rayleigh_pdf(q, sigma, q_min=0.2, q_max=1.0)[source]
Compute truncated Rayleigh PDF for axis ratio.
Computes the probability density function for axis ratio using the truncated Rayleigh distribution with velocity-dependent scale parameter (Wierda et al. 2021 equation C16).
- Parameters:
- q
numpy.ndarray Axis ratios at which to evaluate PDF.
- sigma
numpy.ndarray Velocity dispersions in km/s (same shape as q).
- q_min
float Minimum axis ratio for truncation.
default: 0.2
- q_max
float Maximum axis ratio for truncation.
default: 1.0
- q
- Returns:
- pdf
numpy.ndarray Probability density values with same shape as q.
- pdf
Examples
>>> import numpy as np >>> q = np.array([0.5, 0.7, 0.9]) >>> sigma = np.array([150.0, 200.0, 250.0]) >>> pdf = rayleigh_pdf(q, sigma) >>> print(f"PDF values: {pdf}")
- ler.lens_galaxy_population.sampler_functions.axis_ratio_padilla_strauss_rvs(size)[source]
Sample axis ratios from Padilla & Strauss (2008) distribution.
Uses inverse transform sampling with the empirical PDF from Padilla & Strauss (2008) for early-type galaxy axis ratios.
- Parameters:
- size
int Number of samples to draw.
- size
- Returns:
- q
numpy.ndarray Sampled axis ratios with shape (size,).
- q
Examples
>>> q_samples = axis_ratio_padilla_strauss_rvs(size=1000) >>> print(f"Mean axis ratio: {q_samples.mean():.2f}") >>> print(f"Range: [{q_samples.min():.2f}, {q_samples.max():.2f}]")
- ler.lens_galaxy_population.sampler_functions.axis_ratio_padilla_strauss_pdf(q)[source]
Compute axis ratio PDF from Padilla & Strauss (2008).
Evaluates the probability density function for axis ratio using cubic spline interpolation of the Padilla & Strauss (2008) data.
- Parameters:
- q
numpy.ndarray Axis ratios at which to evaluate PDF.
- q
- Returns:
- pdf
numpy.ndarray Probability density values with same shape as q.
- pdf
Examples
>>> import numpy as np >>> q = np.array([0.3, 0.5, 0.7, 0.9]) >>> pdf = axis_ratio_padilla_strauss_pdf(q) >>> print(f"PDF at q=0.5: {pdf[1]:.4f}")
- ler.lens_galaxy_population.sampler_functions.rejection_sampler_full(size, zs_pdf, number_density, q_pdf, phi_pdf, gamma_pdf, shear1_pdf, shear2_pdf, cross_section, dVdz, ranges, threshold_factor=0.0001, n_prop=1000, lens_type='epl_shear_galaxy', npool=4)[source]
Joint rejection sample over (zs, zl, sigma, nuisance).
Intrinsic target (before cross section) uses
p(zs) * dV/dz_l * n(sigma,zl) * ...;z_lis proposed uniform in[z_min, z_s)— no separatep(zl|zs)factor (seeimportance_sampler_fulldocstring).
- ler.lens_galaxy_population.sampler_functions.rejection_sampler_partial(size, zs, zl, number_density, q_pdf, phi_pdf, gamma_pdf, shear1_pdf, shear2_pdf, cross_section, ranges, threshold_factor=0.0001, n_prop=1000, lens_type='epl_shear_galaxy', npool=4)[source]
- ler.lens_galaxy_population.sampler_functions.importance_sampler_partial(size, zs, zl, number_density, q_pdf, phi_pdf, gamma_pdf, shear1_pdf, shear2_pdf, cross_section, ranges, threshold_factor=0.0001, n_prop=100, lens_type='epl_shear_galaxy')[source]
Importance/resampling for lens nuisance parameters at fixed (zs, zl).
Redshifts are fixed (typically drawn from
zs_slandlens_redshift_sl). Conditional weights targetp(sigma,lambda|zs,zl,SL)up to a constant: proportional ton(sigma,zl), axis-ratio / shear / slope PDFs, and the strong-lensing cross section. Terms that depend only on(zs,zl)--- e.g.P(zs|SL),P(zl|zs,SL),dV/dz_l--- are omitted (constant for fixed redshifts and identical scale across sigma proposals).
- ler.lens_galaxy_population.sampler_functions.importance_sampler_full(size, zs_pdf, number_density, q_pdf, phi_pdf, gamma_pdf, shear1_pdf, shear2_pdf, cross_section, dVdz, ranges, threshold_factor=0.0001, n_prop=50, lens_type='epl_shear_galaxy')[source]
Full joint importance sample over (zs, zl, sigma, nuisance).
Target for lens+source intrinsic factors (before cross section) is taken as
p(zs) * (dV_c/dz_l) * n(sigma, zl) * (axis ratio, shear, slope PDFs),
i.e. no separate
p(zl|zs)factor:z_lis proposed uniformly in[z_min, z_s)and the comoving volume and Schechter-stylen(sigma,zl)carry the lens-redshift abundance. This matchesP(SL|...) p(zs) n(sigma,zl|...) dV/dz_l ...up to overall normalisation.
- ler.lens_galaxy_population.sampler_functions.create_sampler(number_density, q_pdf, phi_pdf, gamma_pdf, shear1_pdf, shear2_pdf, cross_section, dVdz, zs_pdf=None, use_njit_sampler=True, sampler_type='importance_sampler_partial', threshold_factor=0.0001, n_prop=50, lens_type='epl_shear_galaxy', npool=4, **range_kwargs)[source]
Build a closure that runs one of the cross-section-weighted lens samplers.
use_njit_samplerselects the implementation, not only a print message:True(default): Numbanjitkernels.cross_sectionand other callables must be Numba-compatible (e.g.@njitcross sections).False: interpreted Python samplers, which accept ordinary callables (e.g. numerical / SciPy-backed cross sections).