ler.lens_galaxy_population.mp

Multiprocessing helper functions for lens galaxy population calculations.

This module provides parallelized functions for computing optical depth and cross-sections for strong gravitational lensing. These functions are designed to be used with Python’s multiprocessing module for efficient Monte Carlo integration over lens parameters.

Key Features:

  • JIT-compiled parallel function for lens redshift sampling

  • Multiprocessing-compatible wrapper functions for optical depth

  • Cross-section calculation using EPL+Shear lens models

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

Module Contents

Functions

lens_redshift_strongly_lensed_njit(zs_array, ...)

JIT-compiled parallel computation of lens redshift optical depth.

lens_redshift_strongly_lensed_mp(params)

Multiprocessing worker for lens redshift optical depth calculation.

cross_section_unit_mp(params)

Multiprocessing worker for unit Einstein radius cross-section.

cross_section_mp(params)

Multiprocessing worker for cross-section calculation.

Attributes

CS_UNIT_SLOPE

CS_UNIT_INTERCEPT

ler.lens_galaxy_population.mp.CS_UNIT_SLOPE = '0.31830988618379075'[source]
ler.lens_galaxy_population.mp.CS_UNIT_INTERCEPT = '-3.2311742677852644e-27'[source]
ler.lens_galaxy_population.mp.lens_redshift_strongly_lensed_njit(zs_array, zl_scaled, sigma_min, sigma_max, q_rvs, phi_rvs, gamma_rvs, shear_rvs, number_density, cross_section, dVcdz_function, integration_size)[source]

JIT-compiled parallel computation of lens redshift optical depth.

Computes the differential optical depth for strong lensing as a function of source and lens redshifts using Monte Carlo integration with parallel execution over source redshifts.

Parameters:
zs_arraynumpy.ndarray

1D array of source redshifts.

zl_scalednumpy.ndarray

2D array of scaled lens redshifts (zl/zs).

sigma_minfloat

Minimum velocity dispersion (km/s).

sigma_maxfloat

Maximum velocity dispersion (km/s).

q_rvscallable

Function to sample axis ratios given size and sigma.

phi_rvscallable

Function to sample axis rotation angles.

gamma_rvscallable

Function to sample density profile slopes.

shear_rvscallable

Function to sample external shear components (gamma1, gamma2).

number_densitycallable

Function to compute velocity dispersion number density.

cross_sectioncallable

Function to compute lensing cross-section.

dVcdz_functioncallable

Function to compute differential comoving volume.

integration_sizeint

Number of Monte Carlo samples per (zs, zl) pair.

Returns:
result_arraynumpy.ndarray

2D array of optical depth values with shape (len(zs_array), len(zl_scaled[0])).

ler.lens_galaxy_population.mp.lens_redshift_strongly_lensed_mp(params)[source]

Multiprocessing worker for lens redshift optical depth calculation.

Computes the differential optical depth for a single source redshift across multiple lens redshifts. Designed to be called via multiprocessing Pool.map() for parallel computation.

Parameters:
paramstuple

Packed parameters tuple containing:

  • params[0]: Source redshift (float)

  • params[1]: Scaled lens redshift array (1D array)

  • params[2]: Worker index (int)

Returns:
idxint

Worker index for result ordering.

result_arraynumpy.ndarray

1D array of optical depth values for each lens redshift.

ler.lens_galaxy_population.mp.cross_section_unit_mp(params)[source]

Multiprocessing worker for unit Einstein radius cross-section.

Computes the lensing cross-section for a lens with unit Einstein radius (theta_E = 1). Used for building cross-section interpolation grids.

Parameters:
paramstuple

Packed parameters (e1, e2, gamma, gamma1, gamma2, idx).

Returns:
idxint

Worker index for result ordering.

areafloat

Cross-section area (dimensionless, for theta_E = 1).

ler.lens_galaxy_population.mp.cross_section_mp(params)[source]

Multiprocessing worker for cross-section calculation.

Computes the lensing cross-section for given lens parameters. Designed to be called via multiprocessing Pool.map().

Parameters:
paramstuple

Packed parameters (theta_E, e1, e2, gamma, gamma1, gamma2, idx). theta_E is in radians.

Returns:
idxint

Worker index for result ordering.

areafloat

Cross-section area (radian^2).