ler.image_properties
Submodules
Package Contents
Classes
Class to compute image properties of strongly lensed gravitational wave events. |
Functions
|
Solve the lens equation to find image properties. |
|
Convert lens orientation and axis ratio to ellipticity components. |
|
Solve the lens equation to find image properties. |
|
Draw one source position uniformly inside the double caustic. |
|
Scalar version of |
|
Compute the real-valued dot product of two complex numbers. |
|
Convert polar coordinates to elliptical coordinates. |
|
Compute all Hessian-derived local diagnostics at one image position. |
|
Compute the Fermat potential (geometric delay minus lensing potential). |
|
Standalone EPL + shear analytical image finder. |
|
Create a parallel EPL + shear solver for batched lens systems. |
|
Convert lens orientation and axis ratio to ellipticity components. |
|
Convert polar coordinates to elliptical coordinates. |
|
Evaluate the complex angular function Omega for the EPL profile. |
|
Scalar version of |
|
Compute the real-valued dot product of two complex numbers. |
|
Convert polar coordinates to Cartesian coordinates. |
|
Convert Cartesian coordinates to polar coordinates. |
|
Analytically compute the caustics of an EPL + external shear lens model. |
|
Compute the area of a simple polygon using the Shoelace formula. |
|
Compute the area enclosed by the double caustic of an EPL + shear lens. |
|
Create a JIT-compiled cross-section evaluator for batched systems. |
|
Compute double-caustic cross sections for batched lens parameters. |
|
Analytically compute the caustics of an EPL + external shear lens model. |
|
Compute the area of a simple polygon using the Shoelace formula. |
|
Draw one uniform random point inside a polygon via rejection sampling. |
|
Draw multiple uniform random points inside a polygon. |
|
Draw one source position uniformly inside the double caustic. |
|
Draw many source positions uniformly inside the double caustic. |
Attributes
- ler.image_properties.solve_lens_equation(lens_parameters)[source]
Solve the lens equation to find image properties.
Uses the analytical solver from lenstronomy to find image positions, magnifications, time delays, and hessian properties for strongly lensed sources. Source positions are sampled from within the caustic region to ensure multiple imaging.
- Parameters:
- lens_parameters
numpy.ndarray Array of lens configuration parameters with the following structure:
[0]: e1 - ellipticity component 1
[1]: e2 - ellipticity component 2
[2]: gamma - power-law slope of mass density
[3]: gamma1 - external shear component 1
[4]: gamma2 - external shear component 2
[5]: zl - lens redshift
[6]: zs - source redshift
[7]: einstein_radius - Einstein radius (units: radians)
[8]: iteration - iteration index for tracking
- lens_parameters
- Returns:
- x_source
float Source x-position (units: radians).
- y_source
float Source y-position (units: radians).
- x0_image_position
numpy.ndarray Image x-positions (units: radians).
- x1_image_position
numpy.ndarray Image y-positions (units: radians).
- magnifications
numpy.ndarray Magnification factors for each image.
- time_delays
numpy.ndarray Time delays for each image (units: seconds).
- nImages
int Number of images formed.
- determinant
numpy.ndarray Determinant of the lensing Jacobian for each image.
- trace
numpy.ndarray Trace of the lensing Jacobian for each image.
- iteration
int Iteration index passed through for tracking.
- x_source
Examples
>>> from ler.image_properties.multiprocessing_routine import solve_lens_equation, _init_worker_multiprocessing >>> import numpy as np >>> from multiprocessing import Pool >>> lens_parameters1 = np.array([0.024, -0.016, 1.89, 0.10, 0.09, 0.25, 0.94, 2.5e-06, 0]) >>> lens_parameters2 = np.array([-0.040, -0.014, 2.00, 0.08, -0.01, 1.09, 2.55, 1.0e-06, 1]) >>> input_arguments = np.vstack((lens_parameters1, lens_parameters2)) >>> with Pool( ... processes=2, # Number of worker processes ... initializer=_init_worker_multiprocessing, # common ... initargs=( ... 2, # n_min_images ... ['EPL_NUMBA', 'SHEAR'], # lensModelList ... ), ... ) as pool: ... result = pool.map(solve_lens_equation, input_arguments)
- ler.image_properties.phi_q2_ellipticity(phi, q)[source]
Convert position angle and axis ratio to ellipticity components.
- Parameters:
- phi
numpy.ndarray Position angle of the major axis (radians).
- q
numpy.ndarray Axis ratio (0 < q <= 1).
- phi
- Returns:
- e1
numpy.ndarray First ellipticity component.
- e2
numpy.ndarray Second ellipticity component.
- e1
- class ler.image_properties.ImageProperties(npool=4, n_min_images=2, n_max_images=4, lens_model_list=['EPL_NUMBA', 'SHEAR'], image_properties_function='image_properties_epl_shear', image_properties_function_params=None, cosmology=None, time_window=365 * 24 * 3600 * 2, spin_zero=True, spin_precession=False, pdet_finder=None, include_effective_parameters=False, multiprocessing_verbose=True, include_redundant_parameters=False)[source]
Class to compute image properties of strongly lensed gravitational wave events.
This class solves the lens equation to find image positions, magnifications, time delays, and image types (morse phase) for strongly lensed sources. It uses multiprocessing for efficient computation of large samples.
Key Features:
Solves lens equations using multiprocessing for efficiency
Computes image positions, magnifications, and time delays
Classifies image types using morse phase
Calculates detection probabilities for lensed images
- Parameters:
- npool
int Number of processes for multiprocessing.
default: 4
- n_min_images
int Minimum number of images required for a valid lensing event.
default: 2
- n_max_images
int Maximum number of images to consider per event.
default: 4
- time_window
float Time window for lensed events from min(geocent_time) (units: seconds).
default: 365*24*3600*2 (2 years)
- include_effective_parameters
bool Whether to include effective parameters (effective_phase, effective_ra, effective_dec) in the output.
default: True
- lens_model_list
list List of lens models to use.
default: [‘EPL_NUMBA’, ‘SHEAR’]
- cosmology
astropy.cosmologyorNone Cosmology for distance calculations.
If None, uses default LambdaCDM.
default: LambdaCDM(H0=70, Om0=0.3, Ode0=0.7, Tcmb0=0.0, Neff=3.04, m_nu=None, Ob0=0.0)
- spin_zero
bool If True, spin parameters are set to zero (no spin sampling).
default: False
- spin_precession
bool If True (and spin_zero=False), sample precessing spin parameters.
If False (and spin_zero=False), sample aligned/anti-aligned spins.
default: False
- multiprocessing_verbose
bool If True, shows a progress bar for multiprocessing tasks.
default: True
- include_redundant_parameters
bool If True, removes redundant parameters (e.g., theta_E, n_images, mass_1, mass_2, luminosity_distance) from output to save memory.
- npool
Examples
Basic usage:
>>> from ler.image_properties import ImageProperties >>> ip = ImageProperties() >>> lens_parameters = dict( ... zs=np.array([2.0]), ... zl=np.array([0.5]), ... gamma1=np.array([0.0]), ... gamma2=np.array([0.0]), ... phi=np.array([0.0]), ... q=np.array([0.8]), ... gamma=np.array([2.0]), ... theta_E=np.array([1.0]) ... ) >>> result = ip.image_properties(lens_parameters) >>> print(result.keys())
Instance Methods
ImageProperties has the following methods:
Method
Description
Compute image properties for lensed events
Compute detection probability for lensed images
Instance Attributes
ImageProperties has the following attributes:
Attribute
Type
Unit
Description
intNumber of multiprocessing workers
boolIf True, shows a progress bar for multiprocessing tasks
intMinimum number of images required
intMaximum number of images per event
floats
Time window for lensed events
boolTo include effective parameters in output
boolIf True, removes redundant parameters from output to save memory
listList of lens models
astropy.cosmologyCosmology for calculations
boolFlag for zero spin assumption
boolFlag for spin precession
callableProbability of detection calculator
listKeys for probability of detection outputs
- property npool
Number of multiprocessing workers.
- Returns:
- npool
int Number of processes for multiprocessing.
default: 4
- npool
- property n_min_images
Minimum number of images required for a valid lensing event.
- Returns:
- n_min_images
int Minimum number of images required.
default: 2
- n_min_images
- property n_max_images
Maximum number of images per event.
- Returns:
- n_max_images
int Maximum number of images to consider per event.
default: 4
- n_max_images
- property lens_model_list
List of lens models to use.
- Returns:
- lens_model_list
list List of lens model names.
default: [‘EPL_NUMBA’, ‘SHEAR’]
- lens_model_list
- property spin_zero
Flag for zero spin assumption.
- Returns:
- spin_zero
bool Whether to assume zero spin for compact objects.
If True, spin parameters are set to zero (no spin sampling).
If False, spin parameters are sampled.
default: False
- spin_zero
- property spin_precession
Flag for spin precession.
- Returns:
- spin_precession
bool Whether to include spin precession effects.
If True (and spin_zero=False), sample precessing spin parameters.
If False (and spin_zero=False), sample aligned/anti-aligned spins.
default: False
- spin_precession
- property time_window
Time window for lensed events.
- Returns:
- time_window
float Time window for lensed events (units: s).
default: 365*24*3600*20 (20 years)
- time_window
- property cosmo
Astropy cosmology object for calculations.
- Returns:
- cosmo
astropy.cosmology Cosmology used 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)
- cosmo
- property pdet_finder
Detection probability finder function.
- Returns:
- pdet_finder
callable Function that calculates detection probability for GW events.
The function signature should be:
pdet_finder(gw_param_dict) -> dictwith key ‘pdet_net’.
- pdet_finder
- property pdet_finder_output_keys
Output keys from the detection probability finder function.
- Returns:
- pdet_finder_output_keys
list List of keys returned by the pdet_finder function.
default: None
- pdet_finder_output_keys
- property include_effective_parameters
Flag to include effective parameters in output.
- Returns:
- include_effective_parameters
bool Whether to include effective parameters in the output of get_lensed_snrs.
default: False
- include_effective_parameters
- property available_image_properties_functions
Dictionary of available functions for computing image properties.
- Returns:
- available_image_properties_functions
dict Dictionary with function names and default parameters.
- available_image_properties_functions
- image_properties_function
- multiprocessing_verbose = 'True'
- include_redundant_parameters = 'False'
- image_properties_function_params = 'None'
- image_properties_epl_shear_njit(lens_parameters)[source]
Compute image properties for strongly lensed events. This use functions similar to lenstronomy but rewritten in numba njit for speed.
Solves the lens equation using multiprocessing to find image positions, magnifications, time delays, and image types for each lensing event.
- Parameters:
- lens_parameters
dict Dictionary containing lens and source parameters shown in the table:
Parameter
Units
Description
zl
redshift of the lens
zs
redshift of the source
sigma
km s^-1
velocity dispersion
q
axis ratio
theta_E
radian
Einstein radius
phi
rad
axis rotation angle. counter-clockwise from the positive x-axis (RA-like axis) to the major axis of the projected mass distribution.
gamma
density profile slope of EPL galaxy
gamma1
external shear component in the x-direction (RA-like axis)
gamma2
external shear component in the y-direction (Dec-like axis)
- lens_parameters
- Returns:
- lens_parameters
dict Updated dictionary with additional image properties with the following description:
x0_image_positions
radian
x-coordinate (RA-like axis) of the images
x1_image_positions
radian
y-coordinate (Dec-like axis) of the images
magnifications
magnifications
time_delays
time delays
image_type
image type
n_images
number of images
x_source
radian
x-coordinate (RA-like axis) of the source
y_source
radian
y-coordinate (Dec-like axis) of the source
- lens_parameters
- image_properties_epl_shear(lens_parameters)[source]
Compute image properties for strongly lensed events. This use functions from lenstronomy.
Solves the lens equation using multiprocessing to find image positions, magnifications, time delays, and image types for each lensing event.
- Parameters:
- lens_parameters
dict Dictionary containing lens and source parameters shown in the table:
Parameter
Units
Description
zl
redshift of the lens
zs
redshift of the source
sigma
km s^-1
velocity dispersion
q
axis ratio
theta_E
radian
Einstein radius
phi
rad
axis rotation angle. counter-clockwise from the positive x-axis (RA-like axis) to the major axis of the projected mass distribution.
gamma
density profile slope of EPL galaxy
gamma1
external shear component in the x-direction (RA-like axis)
gamma2
external shear component in the y-direction (Dec-like axis)
- lens_parameters
- Returns:
- lens_parameters
dict Updated dictionary with additional image properties with the following description:
x0_image_positions
radian
x-coordinate (RA-like axis) of the images
x1_image_positions
radian
y-coordinate (Dec-like axis) of the images
magnifications
magnifications
time_delays
time delays
image_type
image type
n_images
number of images
x_source
radian
x-coordinate (RA-like axis) of the source
y_source
radian
y-coordinate (Dec-like axis) of the source
- lens_parameters
- get_lensed_snrs(lensed_param, pdet_finder=None, include_effective_parameters=False)[source]
Compute detection probability for each lensed image.
Calculates the effective luminosity distance, geocent time, and phase for each image accounting for magnification and morse phase, then computes detection probabilities using the provided calculator.
- Parameters:
- lensed_param
dict Dictionary containing lensed source and image parameters given below:
Parameter
Units
Description
geocent_time
s
geocent time
ra
rad
right ascension
dec
rad
declination
phase
rad
phase of GW at reference freq
psi
rad
polarization angle
theta_jn
rad
inclination angle
a_1
spin of the primary compact binary
a_2
spin of the secondary compact binary
luminosity_distance
Mpc
luminosity distance of the source
mass_1
Msun
mass of the primary compact binary (detector frame)
mass_2
Msun
mass of the secondary compact binary (detector frame)
x0_image_positions
radian
x-coordinate (RA-like axis) of the images
x1_image_positions
radian
y-coordinate (Dec-like axis) of the images
magnifications
magnifications
time_delays
time delays
image_type
image type
- pdet_finder
callable Function that computes detection probability given GW parameters.
- include_effective_parameters
bool If True, includes effective parameters in output lensed_param.
- lensed_param
- Returns:
- result_dict
dict Dictionary containing:
‘pdet_net’: network detection probability (shape: size x n_max_images)
Individual detector probabilities if pdet_finder outputs them
- lensed_param
dict Updated dictionary with effective parameters shown below:
Parameter
Units
Description
effective_luminosity_distance
Mpc
magnification-corrected distance luminosity_distance / sqrt(|magnifications_i|)
time-delay-corrected GPS time geocent_time + time_delays_i
morse-phase-corrected phase phi - morse_phase_i
effective_geocent_time
s
effective_phase
rad
effective_ra
rad
RA of the image ra + (x0_image_positions_i - x_source)/cos(dec)
effective_dec
rad
Dec of the image dec + (x1_image_positions_i - y_source)
- result_dict
- recover_redundant_parameters(lensed_param)[source]
Recover redundant parameters in lensed_param, i.e. theta_E, n_images, mass_1, mass_2, luminosity_distance.
- produce_effective_params(lensed_param)[source]
Produce effective parameters for each lensed image.
Calculates the effective luminosity distance, geocent time, phase, RA, and Dec for each image accounting for magnification and morse phase.
- Parameters:
- lensed_param
dict Dictionary containing lensed source and image parameters.
- lensed_param
- Returns:
- lensed_param
dict Updated dictionary with effective parameters shown below:
Parameter
Units
Description
effective_luminosity_distance
Mpc
magnification-corrected distance luminosity_distance / sqrt(|magnifications_i|)
time-delay-corrected GPS time geocent_time + time_delays_i
morse-phase-corrected phase phi - morse_phase_i
effective_geocent_time
s
effective_phase
rad
effective_ra
rad
RA of the image ra + (x0_image_positions_i - x_source)/cos(dec)
effective_dec
rad
Dec of the image dec + (x1_image_positions_i - y_source)
- lensed_param
- ler.image_properties.solve_lens_equation(lens_parameters)[source]
Solve the lens equation to find image properties.
Uses the analytical solver from lenstronomy to find image positions, magnifications, time delays, and hessian properties for strongly lensed sources. Source positions are sampled from within the caustic region to ensure multiple imaging.
- Parameters:
- lens_parameters
numpy.ndarray Array of lens configuration parameters with the following structure:
[0]: e1 - ellipticity component 1
[1]: e2 - ellipticity component 2
[2]: gamma - power-law slope of mass density
[3]: gamma1 - external shear component 1
[4]: gamma2 - external shear component 2
[5]: zl - lens redshift
[6]: zs - source redshift
[7]: einstein_radius - Einstein radius (units: radians)
[8]: iteration - iteration index for tracking
- lens_parameters
- Returns:
- x_source
float Source x-position (units: radians).
- y_source
float Source y-position (units: radians).
- x0_image_position
numpy.ndarray Image x-positions (units: radians).
- x1_image_position
numpy.ndarray Image y-positions (units: radians).
- magnifications
numpy.ndarray Magnification factors for each image.
- time_delays
numpy.ndarray Time delays for each image (units: seconds).
- nImages
int Number of images formed.
- determinant
numpy.ndarray Determinant of the lensing Jacobian for each image.
- trace
numpy.ndarray Trace of the lensing Jacobian for each image.
- iteration
int Iteration index passed through for tracking.
- x_source
Examples
>>> from ler.image_properties.multiprocessing_routine import solve_lens_equation, _init_worker_multiprocessing >>> import numpy as np >>> from multiprocessing import Pool >>> lens_parameters1 = np.array([0.024, -0.016, 1.89, 0.10, 0.09, 0.25, 0.94, 2.5e-06, 0]) >>> lens_parameters2 = np.array([-0.040, -0.014, 2.00, 0.08, -0.01, 1.09, 2.55, 1.0e-06, 1]) >>> input_arguments = np.vstack((lens_parameters1, lens_parameters2)) >>> with Pool( ... processes=2, # Number of worker processes ... initializer=_init_worker_multiprocessing, # common ... initargs=( ... 2, # n_min_images ... ['EPL_NUMBA', 'SHEAR'], # lensModelList ... ), ... ) as pool: ... result = pool.map(solve_lens_equation, input_arguments)
- ler.image_properties.sample_source_from_double_caustic(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, max_tries=10)[source]
Draw one source position uniformly inside the double caustic.
- Parameters:
- q
float Lens axis ratio.
- phi
float Lens position angle (rad).
- gamma
float EPL slope parameter.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- theta_E
float Einstein radius used to construct the caustic. default: 1.0
- num_th
int Number of angular samples for the boundary curve. default: 500
- maginf
float Magnification cutoff used in the caustic construction. default: -100.0
- max_tries
int Maximum rejection-sampling attempts. default: 100000
- q
- Returns:
- beta_x
float Sampled source x-coordinate.
- beta_y
float Sampled source y-coordinate.
- ok
int Sampling status flag where 1 indicates success.
- beta_x
Examples
>>> beta_x, beta_y, ok = sample_source_from_double_caustic( ... q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01 ... )
- ler.image_properties.omega_scalar(phi, t, q, niter_max=200, tol=1e-16)[source]
Scalar version of
omegathat avoids temporary array allocation.
- ler.image_properties.cdot(a, b)[source]
Compute the real-valued dot product of two complex numbers.
Equivalent to
Re(a) * Re(b) + Im(a) * Im(b).- Parameters:
- a
complex First complex number.
- b
complex Second complex number.
- a
- Returns:
- result
float Real-valued dot product.
- result
Examples
>>> cdot(1+2j, 3+4j) 11.0
- ler.image_properties.pol_to_ell(r, theta, q)[source]
Convert polar coordinates to elliptical coordinates.
- Parameters:
- r
floatornumpy.ndarray Radial coordinate.
- theta
floatornumpy.ndarray Polar angle in radians.
- q
float Axis ratio.
- r
- Returns:
- rell
floatornumpy.ndarray Elliptical radial coordinate.
- phi
floatornumpy.ndarray Elliptical angle in radians.
- rell
Examples
>>> rell, phi = pol_to_ell(1.0, 0.5, 0.8)
- ler.image_properties.lensing_diagnostics_scalar(x, y, theta_E, gamma, gamma1, gamma2, q, phi, center_x=0.0, center_y=0.0)[source]
Compute all Hessian-derived local diagnostics at one image position.
- Parameters:
- x
float Image-plane x-coordinate.
- y
float Image-plane y-coordinate.
- theta_E
float Einstein radius.
- gamma
float EPL power-law slope.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- q
float Axis ratio.
- phi
float Lens position angle in radians.
- center_x
float Lens center x-coordinate.
- center_y
float Lens center y-coordinate.
- x
- Returns:
- f_xx
float Hessian component d²ψ/dx².
- f_xy
float Hessian component d²ψ/dxdy.
- f_yx
float Hessian component d²ψ/dydx.
- f_yy
float Hessian component d²ψ/dy².
- detA
float Jacobian determinant.
- traceA
float Jacobian trace.
- mu
float Signed magnification.
- image_type
int Image classification.
1 = Type I (minimum), 2 = Type II (saddle),
3 = Type III (maximum), 0 = undefined.
- f_xx
Examples
>>> f_xx, f_xy, f_yx, f_yy, detA, traceA, mu, itype = lensing_diagnostics_scalar( ... x=0.5, y=0.3, theta_E=1.0, gamma=2.0, ... gamma1=0.03, gamma2=-0.01, q=0.8, phi=0.3 ... )
- ler.image_properties.fermat_potential_scalar(x_image, y_image, x_source, y_source, theta_E, gamma, gamma1, gamma2, q, phi, center_x=0.0, center_y=0.0, ra_0=0.0, dec_0=0.0)[source]
Compute the Fermat potential (geometric delay minus lensing potential).
- Parameters:
- x_image
float Image-plane x-coordinate.
- y_image
float Image-plane y-coordinate.
- x_source
float Source-plane x-coordinate.
- y_source
float Source-plane y-coordinate.
- theta_E
float Einstein radius.
- gamma
float EPL power-law slope.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- q
float Axis ratio.
- phi
float Lens position angle in radians.
- center_x
float Lens center x-coordinate.
- center_y
float Lens center y-coordinate.
- ra_0
float Shear center x-coordinate.
- dec_0
float Shear center y-coordinate.
- x_image
- Returns:
- tau
float Fermat potential value.
- tau
Examples
>>> tau = fermat_potential_scalar( ... x_image=0.5, y_image=0.3, x_source=0.1, y_source=0.05, ... theta_E=1.0, gamma=2.0, gamma1=0.03, gamma2=-0.01, q=0.8, phi=0.3 ... )
- ler.image_properties.image_position_analytical_njit(x, y, q, phi, gamma, gamma1, gamma2, theta_E=1.0, alpha_scaling=1.0, magnification_limit=0.01, Nmeas=400, Nmeas_extra=80)[source]
Standalone EPL + shear analytical image finder.
Locates lensed images for a given source position, computes their magnifications, Fermat potential (arrival time proxy), and image types. Results are sorted by ascending arrival time and filtered by a minimum magnification threshold.
- Parameters:
- x
float Source-plane x-coordinate.
- y
float Source-plane y-coordinate.
- q
float Lens axis ratio.
- phi
float Lens position angle in radians.
- gamma
float EPL power-law slope (lenstronomy convention).
The Tessore exponent is
t = gamma - 1.- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- theta_E
float Einstein radius.
default: 1.0
- alpha_scaling
float Deflection scaling factor applied to theta_E.
default: 1.0
- magnification_limit
float Minimum |mu| threshold; images below this are discarded.
default: 0.01
- Nmeas
int Angular root-finding grid size.
default: 400
- Nmeas_extra
int Extra refinement points near the source angle.
default: 80
- x
- Returns:
- x_img
numpy.ndarray Image x-positions sorted by arrival time.
- y_img
numpy.ndarray Image y-positions sorted by arrival time.
- arrival_time
numpy.ndarray Fermat potential values in increasing order.
- magnification
numpy.ndarray Signed magnifications.
- image_type
numpy.ndarray Image classification codes (int64).
1 = Type I (minimum), 2 = Type II (saddle),
3 = Type III (maximum), 0 = undefined.
- nimg
int Number of images returned.
- x_img
Examples
>>> x_img, y_img, tau, mu, itype, n = image_position_analytical_njit( ... x=0.1, y=0.05, q=0.8, phi=0.3, gamma=2.0, ... gamma1=0.03, gamma2=-0.01 ... )
- ler.image_properties.create_epl_shear_solver(arrival_time_sort=True, max_img=4, num_th=500, maginf=-100.0, max_tries=100, alpha_scaling=1.0, magnification_limit=0.01, Nmeas=400, Nmeas_extra=80)[source]
Create a parallel EPL + shear solver for batched lens systems.
Returns a JIT-compiled function that, for each system in a batch, samples a source from the double caustic and solves for image positions, magnifications, time delays, and image types.
- Parameters:
- arrival_time_sort
bool Whether to sort images by arrival time.
default: True
- max_img
int Maximum number of images to store per system.
default: 4
- num_th
int Angular samples for caustic construction.
default: 500
- maginf
float Magnification cutoff for caustic boundary.
default: -100.0
- max_tries
int Maximum rejection-sampling attempts per source.
default: 100
- alpha_scaling
float Deflection scaling factor.
default: 1.0
- magnification_limit
float Minimum |mu| threshold for image retention.
default: 0.01
- Nmeas
int Angular root-finding grid size.
default: 400
- Nmeas_extra
int Extra refinement points.
default: 80
- arrival_time_sort
- Returns:
- solve_epl_shear_multithreaded
callable Parallel solver function with signature
(theta_E, D_dt, q, phi, gamma, gamma1, gamma2)returning a tuple of result arrays.
- solve_epl_shear_multithreaded
Examples
>>> solver = create_epl_shear_solver() >>> results = solver(theta_E, D_dt, q, phi, gamma, gamma1, gamma2)
- ler.image_properties.phi_q2_ellipticity(phi, q)[source]
Convert lens orientation and axis ratio to ellipticity components.
- Parameters:
- phi
float Position angle of the lens major axis in radians.
- q
float Axis ratio (minor/major), where
0 < q <= 1.
- phi
- Returns:
- e1
float First ellipticity component.
- e2
float Second ellipticity component.
- e1
Examples
>>> e1, e2 = phi_q2_ellipticity(phi=0.25, q=0.8)
- ler.image_properties.pol_to_ell(r, theta, q)[source]
Convert polar coordinates to elliptical coordinates.
- Parameters:
- r
floatornumpy.ndarray Radial coordinate.
- theta
floatornumpy.ndarray Polar angle in radians.
- q
float Axis ratio.
- r
- Returns:
- rell
floatornumpy.ndarray Elliptical radial coordinate.
- phi
floatornumpy.ndarray Elliptical angle in radians.
- rell
Examples
>>> rell, phi = pol_to_ell(1.0, 0.5, 0.8)
- ler.image_properties.omega(phi, t, q, omegas, niter_max=200, tol=1e-16)[source]
Evaluate the complex angular function Omega for the EPL profile.
This series expansion converges geometrically with ratio
f = (1 - q)/(1 + q). Thefastmathflag provides ~4x speedup due to the reduction nature of the summation.- Parameters:
- phi
numpy.ndarray Azimuthal angles in radians.
- t
float EPL slope exponent (
t = gamma - 1).- q
float Axis ratio.
- omegas
numpy.ndarray Pre-allocated complex output buffer with the same shape as
phi. Filled in-place.- niter_max
int Maximum number of series terms.
default: 200
- tol
float Convergence tolerance.
default: 1e-16
- phi
- Returns:
- omegas
numpy.ndarray Complex Omega values at each angle (same object as input buffer).
- omegas
Examples
>>> import numpy as np >>> phi = np.linspace(0, 2 * np.pi, 100) >>> result = np.empty_like(phi, dtype=np.complex128) >>> omega(phi, t=1.0, q=0.8, omegas=result)
- ler.image_properties.omega_scalar(phi, t, q, niter_max=200, tol=1e-16)[source]
Scalar version of
omegathat avoids temporary array allocation.
- ler.image_properties.cdot(a, b)[source]
Compute the real-valued dot product of two complex numbers.
Equivalent to
Re(a) * Re(b) + Im(a) * Im(b).- Parameters:
- a
complex First complex number.
- b
complex Second complex number.
- a
- Returns:
- result
float Real-valued dot product.
- result
Examples
>>> cdot(1+2j, 3+4j) 11.0
- ler.image_properties.pol_to_cart(r, th)[source]
Convert polar coordinates to Cartesian coordinates.
- Parameters:
- r
floatornumpy.ndarray Radial coordinate.
- th
floatornumpy.ndarray Polar angle in radians.
- r
- Returns:
- x
floatornumpy.ndarray Cartesian x-coordinate.
- y
floatornumpy.ndarray Cartesian y-coordinate.
- x
Examples
>>> x, y = pol_to_cart(1.0, np.pi / 4)
- ler.image_properties.cart_to_pol(x, y)[source]
Convert Cartesian coordinates to polar coordinates.
The returned angle is wrapped to
[0, 2π).- Parameters:
- x
floatornumpy.ndarray Cartesian x-coordinate.
- y
floatornumpy.ndarray Cartesian y-coordinate.
- x
- Returns:
- r
floatornumpy.ndarray Radial coordinate.
- theta
floatornumpy.ndarray Polar angle in radians, wrapped to
[0, 2π).
- r
Examples
>>> r, theta = cart_to_pol(1.0, 1.0)
- ler.image_properties.caustics_epl_shear(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, sourceplane=True, return_which='double')[source]
Analytically compute the caustics of an EPL + external shear lens model.
For
gamma > 2the outer critical curve does not exist. In that case the routine finds the curve at a finite magnificationmaginfinstead of the true caustic.- Parameters:
- q
float Lens axis ratio.
- phi
float Lens position angle in radians.
- gamma
float EPL power-law slope.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- theta_E
float Einstein radius.
default: 1.0
- num_th
int Number of azimuthal samples for the boundary curve.
default: 500
- maginf
float Magnification cutoff for the outer curve.
default: -100.0
- sourceplane
bool If True, map critical curves to the source plane.
default: True
- return_which
str Which boundary to return.
Options:
‘double’: Double-image region boundary
‘quad’: Quad-image region boundary
‘caustic’: Inner caustic (astroid) only
‘cut’: Outer cut curve only
default: ‘double’
- q
- Returns:
- pts
numpy.ndarray 2×N array of (x, y) boundary coordinates.
- pts
Examples
>>> pts = caustics_epl_shear( ... q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01 ... )
- ler.image_properties.polygon_area(xv, yv)[source]
Compute the area of a simple polygon using the Shoelace formula.
- Parameters:
- xv
numpy.ndarray Polygon vertex x-coordinates.
- yv
numpy.ndarray Polygon vertex y-coordinates.
- xv
- Returns:
- area
float The enclosed geometric area.
- area
Examples
>>> area = polygon_area(np.array([0.0, 1.0, 0.0]), np.array([0.0, 0.0, 1.0]))
- ler.image_properties.caustic_double_area(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, sourceplane=True, return_which='double')[source]
Compute the area enclosed by the double caustic of an EPL + shear lens.
- Parameters:
- q
float Lens axis ratio.
- phi
float Lens position angle in radians.
- gamma
float Power-law slope.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- theta_E
float Einstein radius. default: 1.0
- num_th
int Number of angular samples. default: 500
- maginf
float Magnification cutoff parameter. default: -100.0
- q
- Returns:
- area
float Area of the double caustic in units of theta_E^2.
- area
Examples
>>> area = caustic_double_area(0.8, 0.0, 2.0, 0.03, -0.01)
- ler.image_properties.make_cross_section_reinit(Da_instance, num_th=500, maginf=-100.0, sourceplane=True, return_which='double')[source]
Create a JIT-compiled cross-section evaluator for batched systems.
- Parameters:
- Da_instance
callable Angular-diameter-distance function that accepts redshift arrays.
- num_th
int Number of angular samples for caustic construction. default: 500
- maginf
float Magnification cutoff used by the caustic routine. default: -100.0
- Da_instance
- Returns:
- cross_section_reinit
callable Parallel numba function with signature
(zs, zl, sigma, q, phi, gamma, gamma1, gamma2)returning the double-caustic cross section array.
- cross_section_reinit
Examples
>>> cross_section = make_cross_section_reinit(Da_instance) >>> cs = cross_section(zs, zl, sigma, q, phi, gamma, gamma1, gamma2)
- ler.image_properties.cross_section_epl_shear_unit(e1, e2, gamma, gamma1, gamma2, theta_E=None, num_th=500, maginf=-100.0, sourceplane=True, return_which='double')[source]
Compute double-caustic cross sections for batched lens parameters.
- Parameters:
- e1
numpy.ndarray Lens ellipticity component 1.
- e2
numpy.ndarray Lens ellipticity component 2.
- gamma
numpy.ndarray EPL slopes.
- gamma1
numpy.ndarray External shear component 1.
- gamma2
numpy.ndarray External shear component 2.
- theta_E
numpy.ndarrayorNone Einstein radius in angular units. If None, assumed to be 1.0 for all
- Returns
- -------
- cs
numpy.ndarray Cross section values in angular units.
- e1
- ler.image_properties.caustics_epl_shear(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, sourceplane=True, return_which='double')[source]
Analytically compute the caustics of an EPL + external shear lens model.
For
gamma > 2the outer critical curve does not exist. In that case the routine finds the curve at a finite magnificationmaginfinstead of the true caustic.- Parameters:
- q
float Lens axis ratio.
- phi
float Lens position angle in radians.
- gamma
float EPL power-law slope.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- theta_E
float Einstein radius.
default: 1.0
- num_th
int Number of azimuthal samples for the boundary curve.
default: 500
- maginf
float Magnification cutoff for the outer curve.
default: -100.0
- sourceplane
bool If True, map critical curves to the source plane.
default: True
- return_which
str Which boundary to return.
Options:
‘double’: Double-image region boundary
‘quad’: Quad-image region boundary
‘caustic’: Inner caustic (astroid) only
‘cut’: Outer cut curve only
default: ‘double’
- q
- Returns:
- pts
numpy.ndarray 2×N array of (x, y) boundary coordinates.
- pts
Examples
>>> pts = caustics_epl_shear( ... q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01 ... )
- ler.image_properties.polygon_area(xv, yv)[source]
Compute the area of a simple polygon using the Shoelace formula.
- Parameters:
- xv
numpy.ndarray Polygon vertex x-coordinates.
- yv
numpy.ndarray Polygon vertex y-coordinates.
- xv
- Returns:
- area
float The enclosed geometric area.
- area
Examples
>>> area = polygon_area(np.array([0.0, 1.0, 0.0]), np.array([0.0, 0.0, 1.0]))
- ler.image_properties.sample_uniform_in_polygon(xv, yv, max_tries=100000)[source]
Draw one uniform random point inside a polygon via rejection sampling.
- Parameters:
- xv
numpy.ndarray Polygon vertex x-coordinates.
- yv
numpy.ndarray Polygon vertex y-coordinates.
- max_tries
int Maximum number of rejection-sampling attempts. default: 100000
- xv
- Returns:
- x
float Sampled x-coordinate.
- y
float Sampled y-coordinate.
- ok
int Sampling flag where 1 indicates success and 0 indicates failure.
- x
Examples
>>> import numpy as np >>> xv = np.array([0.0, 1.0, 1.0, 0.0]) >>> yv = np.array([0.0, 0.0, 1.0, 1.0]) >>> x, y, ok = sample_uniform_in_polygon(xv, yv)
- ler.image_properties.sample_many_uniform_in_polygon(xv, yv, n_samples, max_tries_per=100000)[source]
Draw multiple uniform random points inside a polygon.
Sampling uses batched rejection from the polygon bounding box.
- Parameters:
- xv
numpy.ndarray Polygon vertex x-coordinates.
- yv
numpy.ndarray Polygon vertex y-coordinates.
- n_samples
int Number of requested samples.
- max_tries_per
int Maximum attempts per requested sample. default: 100000
- xv
- Returns:
- x_src
numpy.ndarray Sampled x-coordinates with capacity
n_samples.- y_src
numpy.ndarray Sampled y-coordinates with capacity
n_samples.- n_ok
int Number of valid samples in the leading entries.
- x_src
Examples
>>> import numpy as np >>> xv = np.array([0.0, 1.0, 1.0, 0.0]) >>> yv = np.array([0.0, 0.0, 1.0, 1.0]) >>> x_src, y_src, n_ok = sample_many_uniform_in_polygon(xv, yv, 100)
- ler.image_properties.sample_source_from_double_caustic(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, max_tries=10)[source]
Draw one source position uniformly inside the double caustic.
- Parameters:
- q
float Lens axis ratio.
- phi
float Lens position angle (rad).
- gamma
float EPL slope parameter.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- theta_E
float Einstein radius used to construct the caustic. default: 1.0
- num_th
int Number of angular samples for the boundary curve. default: 500
- maginf
float Magnification cutoff used in the caustic construction. default: -100.0
- max_tries
int Maximum rejection-sampling attempts. default: 100000
- q
- Returns:
- beta_x
float Sampled source x-coordinate.
- beta_y
float Sampled source y-coordinate.
- ok
int Sampling status flag where 1 indicates success.
- beta_x
Examples
>>> beta_x, beta_y, ok = sample_source_from_double_caustic( ... q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01 ... )
- ler.image_properties.sample_many_sources_from_double_caustic(q, phi, gamma, gamma1, gamma2, n_samples, theta_E=1.0, num_th=500, maginf=-100.0, max_tries_per=10)[source]
Draw many source positions uniformly inside the double caustic.
- Parameters:
- q
float Lens axis ratio.
- phi
float Lens position angle (rad).
- gamma
float EPL slope parameter.
- gamma1
float External shear component 1.
- gamma2
float External shear component 2.
- n_samples
int Number of source samples requested.
- theta_E
float Einstein radius used to construct the caustic. default: 1.0
- num_th
int Number of angular samples for the boundary curve. default: 500
- maginf
float Magnification cutoff used in the caustic construction. default: -100.0
- max_tries_per
int Maximum rejection attempts per requested sample. default: 100000
- q
- Returns:
- beta_x
numpy.ndarray Sampled source x-coordinates.
- beta_y
numpy.ndarray Sampled source y-coordinates.
- n_ok
int Number of valid samples in the leading entries.
- beta_x
Examples
>>> bx, by, n_ok = sample_many_sources_from_double_caustic( ... q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01, n_samples=1000 ... )