ler.lens_galaxy_population.optical_depth

Module for optical depth and lens parameter distribution calculations.

This module provides the OpticalDepth class for computing strong lensing optical depth, cross-section, sampling velocity dispersion, axis ratio, and other lens galaxy population parameters. It supports multiple lens models including SIS, SIE, and EPL + external shear.

Key features:

  • Optical depth computation for strong gravitational lensing

  • Velocity dispersion sampling with multiple models

  • Lens redshift distribution sampling

  • Cross-section calculations for various lens models

Copyright (C) 2024 Hemantakumar Phurailatpam. Distributed under MIT License.

Module Contents

Classes

OpticalDepth

Class for computing optical depth and lens galaxy population parameters.

class ler.lens_galaxy_population.optical_depth.OpticalDepth(npool=4, z_min=0.0, z_max=10.0, cosmology=None, lens_type='epl_shear_galaxy', lens_functions=None, lens_functions_params=None, lens_param_samplers=None, lens_param_samplers_params=None, directory='./interpolator_json', create_new_interpolator=False, verbose=False)[source]

Class for computing optical depth and lens galaxy population parameters.

This class calculates strong lensing optical depth, velocity dispersion, axis ratio, and other parameters for a lens galaxy population. It supports SIS, SIE, and EPL + external shear lens models with customizable samplers and interpolators for efficient computation.

Key Features:

  • Multiple lens model support (SIS, SIE, EPL + shear)

  • Configurable velocity dispersion distributions

  • Cached interpolators for fast optical depth computation

  • Flexible parameter sampling with user-defined priors

Parameters:
npoolint

Number of processors for multiprocessing.

default: 4

z_minfloat

Minimum redshift of the lens galaxy population.

default: 0.0

z_maxfloat

Maximum redshift of the lens galaxy population.

default: 10.0

cosmologyastropy.cosmology or None

Cosmology object for distance calculations.

default: LambdaCDM(H0=70, Om0=0.3, Ode0=0.7)

lens_typestr

Type of lens galaxy model.

Options:

  • ‘epl_shear_galaxy’: Elliptical power-law with external shear

  • ‘sie_galaxy’: Singular isothermal ellipsoid

  • ‘sis_galaxy’: Singular isothermal sphere

default: ‘epl_shear_galaxy’

lens_functionsdict or None

Dictionary with lens-related functions.

default: None (uses defaults for lens_type)

lens_functions_paramsdict or None

Dictionary with parameters for lens-related functions.

default: None

lens_param_samplersdict or None

Dictionary of sampler functions for lens parameters.

default: None (uses defaults for lens_type)

lens_param_samplers_paramsdict or None

Dictionary with parameters for the samplers.

default: None

directorystr

Directory where interpolators are saved.

default: ‘./interpolator_json’

create_new_interpolatorbool or dict

Whether to create new interpolators.

default: False

verbosebool

If True, prints additional information.

default: False

Examples

Basic usage:

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> tau = od.optical_depth(zs=np.array([1.0, 2.0]))

Instance Methods

OpticalDepth has the following instance methods:

Instance Attributes

OpticalDepth has the following instance attributes:

Attribute

Type

Unit

Description

npool

int

Number of processors for multiprocessing

z_min

float

Minimum redshift

z_max

float

Maximum redshift

cosmo

astropy.cosmology

Cosmology object

lens_type

str

Type of lens galaxy model

directory

str

Directory for interpolator storage

optical_depth

FunctionConditioning

Optical depth calculator

velocity_dispersion

FunctionConditioning

km/s

Velocity dispersion sampler

axis_ratio

FunctionConditioning

Axis ratio sampler

axis_rotation_angle

FunctionConditioning

rad

Axis rotation angle sampler

lens_redshift

FunctionConditioning

Lens redshift sampler

external_shear

FunctionConditioning

External shear sampler

density_profile_slope

FunctionConditioning

Density profile slope sampler

cross_section

callable

rad²

Cross-section calculator

available_lens_samplers

dict

Available lens parameter samplers

available_lens_functions

dict

Available lens functions

property lens_type[source]

Type of lens galaxy model.

Returns:
lens_typestr

Lens type (‘epl_shear_galaxy’, ‘sie_galaxy’, or ‘sis_galaxy’).

property npool[source]

Number of processors for multiprocessing.

Returns:
npoolint

Number of parallel processors.

property z_min[source]

Minimum redshift of the lens galaxy population.

Returns:
z_minfloat

Minimum redshift.

property z_max[source]

Maximum redshift of the lens galaxy population.

Returns:
z_maxfloat

Maximum redshift.

property cosmo[source]

Cosmology object for distance calculations.

Returns:
cosmoastropy.cosmology

Cosmology object.

property directory[source]

Directory for interpolator storage.

Returns:
directorystr

Path to interpolator JSON files.

property velocity_dispersion[source]

Velocity dispersion sampler object.

Returns a FunctionConditioning object with methods:

  • rvs(size, zl): Sample velocity dispersion values

  • pdf(sigma, zl): Get probability density

  • function(sigma, zl): Get number density function

Returns:
velocity_dispersionFunctionConditioning

Sampler object for velocity dispersion (km/s).

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> sigma = od.velocity_dispersion(size=100, zl=np.ones(100)*0.5)
property axis_ratio[source]

Axis ratio sampler object.

Returns a FunctionConditioning object with methods:

  • rvs(size, sigma): Sample axis ratio values

  • pdf(q, sigma): Get probability density

  • function(q, sigma): Get distribution function

Returns:
axis_ratioFunctionConditioning

Sampler object for axis ratio.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> q = od.axis_ratio(size=100, sigma=np.ones(100)*200.)
property axis_rotation_angle[source]

Axis rotation angle sampler object.

Returns a FunctionConditioning object with methods:

  • rvs(size): Sample axis rotation angles

  • pdf(phi): Get probability density

Returns:
axis_rotation_angleFunctionConditioning

Sampler object for axis rotation angle (rad).

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> phi = od.axis_rotation_angle(size=100)
property density_profile_slope[source]

Density profile slope sampler object.

Returns a FunctionConditioning object with methods:

  • rvs(size): Sample density profile slope values

  • pdf(gamma): Get probability density

Returns:
density_profile_slopeFunctionConditioning

Sampler object for density profile slope.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> gamma = od.density_profile_slope(size=100)
property external_shear[source]

External shear sampler object.

Returns a FunctionConditioning object with methods:

  • rvs(size): Sample shear components (gamma1, gamma2)

  • pdf(gamma1, gamma2): Get probability density

Returns:
external_shearFunctionConditioning

Sampler object for external shear.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> gamma1, gamma2 = od.external_shear(size=100)
property cross_section[source]

Lensing cross-section calculator.

Returns a callable that computes lensing cross-section for individual

lensing events. Input parameters depend on lens type:

  • EPL+shear: zs, zl, sigma, q, phi, gamma, gamma1, gamma2

  • SIE: zs, zl, sigma, q

  • SIS: zs, zl, sigma

Returns:
cross_sectioncallable

Cross-section function (rad² units).

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> cs = od.cross_section(zs=zs, zl=zl, sigma=sigma, ...)
property lens_redshift[source]

Lens redshift sampler object.

Returns a FunctionConditioning object with methods:

  • rvs(size, zs): Sample lens redshifts given source redshifts

  • pdf(zl, zs): Get probability density

  • function(zl, zs): Get effective lensing cross-section

Returns:
lens_redshiftFunctionConditioning

Sampler object for lens redshift.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> zl = od.lens_redshift(size=100, zs=np.ones(100)*2.0)
property density_profile_slope_sl[source]

Density profile slope sampler object (strong lensing conditioned).

Returns a FunctionConditioning object with methods:

  • rvs(size): Sample density profile slope values

  • pdf(gamma): Get probability density

Returns:
density_profile_slope_slFunctionConditioning

Sampler object for density profile slope (strong lensing).

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> gamma = od.density_profile_slope_sl(size=100)
property external_shear_sl[source]

External shear sampler object (strong lensing conditioned).

Returns a FunctionConditioning object with methods:

  • rvs(size): Sample shear components (gamma1, gamma2)

  • pdf(gamma1, gamma2): Get probability density

Returns:
external_shear_slFunctionConditioning

Sampler object for external shear (strong lensing).

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> gamma1, gamma2 = od.external_shear_sl(size=100)
property optical_depth[source]

Strong lensing optical depth calculator.

Returns:
optical_depthFunctionConditioning

Function object with .function(zs) method that returns n optical depth for given source redshifts.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> tau = od.optical_depth.function(np.array([1.0, 2.0]))
property available_lens_samplers[source]

Dictionary of available lens parameter samplers and their default parameters.

Returns:
available_lens_samplersdict

Dictionary with sampler names and default parameters.

property available_lens_functions[source]

Dictionary of available lens functions and their default parameters.

Returns:
available_lens_functionsdict

Dictionary with function names and default parameters.

comoving_distance[source]
angular_diameter_distance[source]
angular_diameter_distance_z1z2[source]
differential_comoving_volume[source]
lens_redshift_intrinsic[source]
lens_redshift_strongly_lensed_numerical(size=1000, zs=None, get_attribute=False, **kwargs)[source]

Sample lens redshifts conditioned on strong lensing (numerical method).

This method computes the lens redshift distribution by numerically integrating over the velocity dispersion distribution (galaxy density distribution wrt), cross-section and differential comoving volume.

Parameters:
sizeint

Number of samples to generate. n default: 1000

zsnumpy.ndarray

Source redshifts.

get_attributebool

If True, returns the sampler object instead of samples. n default: False

**kwargsdict

Additional parameters.

Returns:
zlnumpy.ndarray or FunctionConditioning

Lens redshift samples or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> zl = od.lens_redshift(size=100, zs=np.ones(100)*2.0)
lens_redshift_strongly_lensed_sis_haris(size, zs, get_attribute=False, **kwargs)[source]

Sample SIS lens redshifts using Haris et al. (2018) distribution.

Parameters:
sizeint

Number of samples to generate.

zsnumpy.ndarray

Source redshifts.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters.

Returns:
zlnumpy.ndarray or FunctionConditioning

Lens redshift samples or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth(lens_type='sis_galaxy')
>>> zl = od.lens_redshift(size=100, zs=np.ones(100)*2.0)
velocity_dispersion_gengamma(size, get_attribute=False, **kwargs)[source]

Sample velocity dispersion from generalized gamma distribution.

Parameters:
sizeint

Number of samples to generate.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters.

Parameter

Unit

Description

sigma_min

km/s

minimum velocity dispersion

sigma_max

km/s

maximum velocity dispersion

alpha

dimensionless

Power-law index governing the slope of the distribution at low velocities

beta

dimensionless

Exponential parameter determining the sharpness of the high-velocity cutoff

phistar

h^3 Mpc^-3

Normalization constant representing the comoving number density of galaxy

sigmastar

km/s

Characteristic velocity scale marking the “knee” transition in the VDF

Returns:
sigmanumpy.ndarray or FunctionConditioning

Velocity dispersion samples (km/s) or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth(lens_param_samplers=dict(
...     velocity_dispersion="velocity_dispersion_gengamma"))
>>> sigma = od.velocity_dispersion(size=100)
velocity_dispersion_bernardi(size, get_attribute=False, **kwargs)[source]

Sample velocity dispersion from Bernardi et al. (2010) distribution.

Uses inverse transform sampling on the velocity dispersion function.

Parameters:
sizeint

Number of samples to generate.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters.

Parameter

Unit

Description

sigma_min

km/s

minimum velocity dispersion

sigma_max

km/s

maximum velocity dispersion

alpha

dimensionless

Power-law index governing the slope of the distribution at low velocities

beta

dimensionless

Exponential parameter determining the sharpness of the high-velocity cutoff

phistar

h^3 Mpc^-3

Normalization constant representing the comoving number density of galaxy

sigmastar

km/s

Characteristic velocity scale marking the “knee” transition in the VDF

Returns:
sigmanumpy.ndarray or FunctionConditioning

Velocity dispersion samples (km/s) or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth(lens_param_samplers=dict(
...     velocity_dispersion="velocity_dispersion_bernardi"))
>>> sigma = od.velocity_dispersion(size=100)
velocity_dispersion_ewoud(size, zl, get_attribute=False, **kwargs)[source]

Sample redshift-dependent velocity dispersion from Wempe et al. (2022).

Uses inverse transform sampling with redshift-dependent velocity dispersion function.

Parameters:
sizeint

Number of samples to generate.

zlnumpy.ndarray

Lens redshifts.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters.

Parameter

Unit

Description

sigma_min

km/s

minimum velocity dispersion

sigma_max

km/s

maximum velocity dispersion

alpha

dimensionless

Power-law index governing the slope of the distribution at low velocities

beta

dimensionless

Exponential parameter determining the sharpness of the high-velocity cutoff

phistar

h^3 Mpc^-3

Normalization constant representing the comoving number density of galaxy

sigmastar

km/s

Characteristic velocity scale marking the “knee” transition in the VDF

Returns:
sigmanumpy.ndarray or FunctionConditioning

Velocity dispersion samples (km/s) or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> sigma = od.velocity_dispersion(size=100, zl=np.ones(100)*0.5)
axis_ratio_rayleigh(size, sigma, get_attribute=False, **kwargs)[source]

Sample axis ratio from Rayleigh distribution conditioned on velocity dispersion.

Parameters:
sizeint

Number of samples to generate.

sigmanumpy.ndarray

Velocity dispersion of the lens galaxy (km/s).

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters (q_min, q_max).

Returns:
qnumpy.ndarray or FunctionConditioning

Axis ratio samples or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth(lens_param_samplers=dict(axis_ratio="axis_ratio_rayleigh"))
>>> q = od.axis_ratio(size=100, sigma=np.ones(100)*200.)
axis_ratio_padilla_strauss(size=1000, get_attribute=False, **kwargs)[source]

Sample axis ratio from Padilla & Strauss (2008) distribution.

Parameters:
sizeint

Number of samples to generate.

default: 1000

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters (q_min, q_max).

Returns:
qnumpy.ndarray or FunctionConditioning

Axis ratio samples or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth(lens_param_samplers=dict(axis_ratio="axis_ratio_padilla_strauss"))
>>> q = od.axis_ratio(size=100)
axis_rotation_angle_uniform(size, get_attribute=False, **kwargs)[source]

Sample axis rotation angle from uniform distribution.

Parameters:
sizeint

Number of samples to generate.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters (phi_min, phi_max).

Returns:
phinumpy.ndarray or FunctionConditioning

Axis rotation angle samples (rad) or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> phi = od.axis_rotation_angle(size=100)
axis_ratio_uniform(size, get_attribute=False, **kwargs)[source]

Sample axis ratio from uniform distribution.

Parameters:
sizeint

Number of samples to generate.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters (q_min, q_max).

Returns:
qnumpy.ndarray or FunctionConditioning

Axis ratio samples or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth(lens_param_samplers=dict(axis_ratio="axis_ratio_uniform"))
>>> q = od.axis_ratio(size=100)
external_shear_normal(size, get_attribute=False, **kwargs)[source]

Sample external shear parameters from 2D normal distribution.

Parameters:
sizeint

Number of samples to generate.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters (mean, std).

Returns:
shearnumpy.ndarray or FunctionConditioning

Array of shape (2, size) with gamma1, gamma2 or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> gamma1, gamma2 = od.external_shear(size=100)
density_profile_slope_normal(size, get_attribute=False, **kwargs)[source]

Sample density profile slope from normal distribution.

Parameters:
sizeint

Number of samples to generate.

get_attributebool

If True, returns the sampler object instead of samples.

default: False

**kwargsdict

Additional parameters (mean, std).

Returns:
gammanumpy.ndarray or FunctionConditioning

Density profile slope samples or sampler object.

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> gamma = od.density_profile_slope(size=100)
optical_depth_numerical(zs, get_attribute=False, **kwargs)[source]

Helper to compute optical depth numerically by integrating lens redshift.

Parameters:
zsnumpy.ndarray

Source redshifts.

get_attributebool

If True, returns the function object instead of values. n default: False

**kwargsdict

Additional parameters.

Returns:
taunumpy.ndarray or FunctionConditioning

Optical depth values or function object.

compute_einstein_radii(sigma, zl, zs)[source]

Function to compute the Einstein radii of the lens galaxies

Parameters:
sigmafloat

velocity dispersion of the lens galaxy

zlfloat

lens redshifts

zsfloat

source redshifts

Returns:
theta_Efloat

Einstein radii of the lens galaxies in radians. Multiply by

Examples

>>> from ler.lens_galaxy_population import LensGalaxyParameterDistribution
>>> lens = LensGalaxyParameterDistribution()
>>> sigma = 200.0
>>> zl = 0.5
>>> zs = 1.0
>>> lens.compute_einstein_radii(sigma, zl, zs)
optical_depth_sis_analytic(zs, get_attribute=False, **kwargs)[source]

Function to compute the strong lensing optical depth (SIS). LambdaCDM(H0=70, Om0=0.3, Ode0=0.7) was used to derive the following equation. This is the analytic version of optical depth from z=0 to z=zs.

Parameters:
zsnumpy.ndarray

Source redshifts.

get_attributebool

If True, returns the function object instead of values. n default: False

**kwargsdict

Additional parameters.

Returns:
taunumpy.ndarray or FunctionConditioning

Optical depth values or function object.

cross_section_sis(zs=None, zl=None, sigma=None, get_attribute=False, **kwargs)[source]

Function to compute the SIS cross-section

Parameters:
sigmafloat

velocity dispersion of the lens galaxy

zlfloat

redshift of the lens galaxy

zsfloat

redshift of the source galaxy

Returns:
cross_sectionfloat

SIS cross-section

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> print(self.cross_section_sis(sigma=200., zl=0.5, zs=1.0))
cross_section_sie_feixu(zs=None, zl=None, sigma=None, q=None, get_attribute=False, **kwargs)[source]

Function to compute the SIE cross-section from Fei Xu et al. (2021)

Parameters:
sigmafloat

velocity dispersion of the lens galaxy

zlfloat

redshift of the lens galaxy

zsfloat

redshift of the source galaxy

Returns:
cross_sectionfloat

SIE cross-section

Examples

>>> from ler.lens_galaxy_population import OpticalDepth
>>> od = OpticalDepth()
>>> print(self.cross_section_sie_feixu(sigma=200., zl=0.5, zs=1.0, q=1.0))
cross_section_epl_shear_numerical(zs, zl, sigma, q, phi, gamma, gamma1, gamma2)[source]

Function to compute the strong lensing cross-section numerically for EPL + external shear lenses.

Parameters:
zsnumpy.ndarray

redshift of the source galaxies

zlnumpy.ndarray

redshift of the lens galaxies

sigmanumpy.ndarray

velocity dispersion of the lens galaxies

qnumpy.ndarray

axis ratios of the lens galaxies

phinumpy.ndarray

axis rotation angles of the lens galaxies in radians

gammanumpy.ndarray

external shear magnitudes of the lens galaxies

gamma1numpy.ndarray

external shear component 1 of the lens galaxies

gamma2numpy.ndarray

external shear component 2 of the lens galaxies

Returns:
cross_sectionnumpy.ndarray

strong lensing cross-section of the lens galaxies in square radians

cross_section_epl_shear_numerical_mp(theta_E, gamma, gamma1, gamma2, q=None, phi=None, e1=None, e2=None, verbose=False, **kwargs)[source]

Function to compute the strong lensing cross-section numerically for EPL + external shear lenses.

Parameters:
theta_Enumpy.ndarray

Einstein radii of the lens galaxies in radians

gammanumpy.ndarray

external shear magnitudes of the lens galaxies

gamma1numpy.ndarray

external shear component 1 of the lens galaxies

gamma2numpy.ndarray

external shear component 2 of the lens galaxies

qnumpy.ndarray

axis ratios of the lens galaxies

phinumpy.ndarray

axis rotation angles of the lens galaxies in radians

e1numpy.ndarray

ellipticity component 1 of the lens galaxies

e2numpy.ndarray

ellipticity component 2 of the lens galaxies

Returns:
cross_sectionnumpy.ndarray

strong lensing cross-section of the lens galaxies in square radians

create_parameter_grid(size_list=[25, 25, 45, 15, 15])[source]

Create a parameter grid for lens galaxies.

Parameters:
size_listlist

List of sizes for each parameter grid.

Returns:
zlnumpy.ndarray

Lens redshifts.

sigmanumpy.ndarray

Velocity dispersions.

qnumpy.ndarray

Axis ratios.

cross_section_epl_shear_interpolation_init(file_path, size_list)[source]
cross_section_epl_shear_interpolation(zs, zl, sigma, q, phi, gamma, gamma1, gamma2, get_attribute=False, size_list=[25, 25, 45, 15, 15], **kwargs)[source]

Function to compute the cross-section correction factor