ler.image_properties.epl_shear_njit

Module for semi-analytical image finding in EPL + external shear lens models.

Provides a numba-accelerated solver that locates multiple lensed images, computes their magnifications, arrival times, and image types for an elliptical power-law (EPL) mass profile with external shear.

Usage:

Solve for image positions of a single source:

>>> from ler.image_properties.epl_shear_njit import image_position_analytical_njit
>>> 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
... )

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

Module Contents

Functions

lensing_diagnostics_scalar(x, y, theta_E, gamma, ...)

Compute all Hessian-derived local diagnostics at one image position.

fermat_potential_scalar(x_image, y_image, x_source, ...)

Compute the Fermat potential (geometric delay minus lensing potential).

image_position_analytical_njit(x, y, q, phi, gamma, ...)

Standalone EPL + shear analytical image finder.

create_epl_shear_solver([arrival_time_sort, max_img, ...])

Create a parallel EPL + shear solver for batched lens systems.

Attributes

EPS

MAX_ROOTS

MAX_IMGS

C_LIGHT

ler.image_properties.epl_shear_njit.EPS = '1e-14'[source]
ler.image_properties.epl_shear_njit.MAX_ROOTS = '64'[source]
ler.image_properties.epl_shear_njit.MAX_IMGS = '5'[source]
ler.image_properties.epl_shear_njit.C_LIGHT = '299792458.0'[source]
ler.image_properties.epl_shear_njit.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:
xfloat

Image-plane x-coordinate.

yfloat

Image-plane y-coordinate.

theta_Efloat

Einstein radius.

gammafloat

EPL power-law slope.

gamma1float

External shear component 1.

gamma2float

External shear component 2.

qfloat

Axis ratio.

phifloat

Lens position angle in radians.

center_xfloat

Lens center x-coordinate.

center_yfloat

Lens center y-coordinate.

Returns:
f_xxfloat

Hessian component d²ψ/dx².

f_xyfloat

Hessian component d²ψ/dxdy.

f_yxfloat

Hessian component d²ψ/dydx.

f_yyfloat

Hessian component d²ψ/dy².

detAfloat

Jacobian determinant.

traceAfloat

Jacobian trace.

mufloat

Signed magnification.

image_typeint

Image classification.

1 = Type I (minimum), 2 = Type II (saddle),

3 = Type III (maximum), 0 = undefined.

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.epl_shear_njit.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_imagefloat

Image-plane x-coordinate.

y_imagefloat

Image-plane y-coordinate.

x_sourcefloat

Source-plane x-coordinate.

y_sourcefloat

Source-plane y-coordinate.

theta_Efloat

Einstein radius.

gammafloat

EPL power-law slope.

gamma1float

External shear component 1.

gamma2float

External shear component 2.

qfloat

Axis ratio.

phifloat

Lens position angle in radians.

center_xfloat

Lens center x-coordinate.

center_yfloat

Lens center y-coordinate.

ra_0float

Shear center x-coordinate.

dec_0float

Shear center y-coordinate.

Returns:
taufloat

Fermat potential value.

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.epl_shear_njit.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:
xfloat

Source-plane x-coordinate.

yfloat

Source-plane y-coordinate.

qfloat

Lens axis ratio.

phifloat

Lens position angle in radians.

gammafloat

EPL power-law slope (lenstronomy convention).

The Tessore exponent is t = gamma - 1.

gamma1float

External shear component 1.

gamma2float

External shear component 2.

theta_Efloat

Einstein radius.

default: 1.0

alpha_scalingfloat

Deflection scaling factor applied to theta_E.

default: 1.0

magnification_limitfloat

Minimum |mu| threshold; images below this are discarded.

default: 0.01

Nmeasint

Angular root-finding grid size.

default: 400

Nmeas_extraint

Extra refinement points near the source angle.

default: 80

Returns:
x_imgnumpy.ndarray

Image x-positions sorted by arrival time.

y_imgnumpy.ndarray

Image y-positions sorted by arrival time.

arrival_timenumpy.ndarray

Fermat potential values in increasing order.

magnificationnumpy.ndarray

Signed magnifications.

image_typenumpy.ndarray

Image classification codes (int64).

1 = Type I (minimum), 2 = Type II (saddle),

3 = Type III (maximum), 0 = undefined.

nimgint

Number of images returned.

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.epl_shear_njit.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_sortbool

Whether to sort images by arrival time.

default: True

max_imgint

Maximum number of images to store per system.

default: 4

num_thint

Angular samples for caustic construction.

default: 500

maginffloat

Magnification cutoff for caustic boundary.

default: -100.0

max_triesint

Maximum rejection-sampling attempts per source.

default: 100

alpha_scalingfloat

Deflection scaling factor.

default: 1.0

magnification_limitfloat

Minimum |mu| threshold for image retention.

default: 0.01

Nmeasint

Angular root-finding grid size.

default: 400

Nmeas_extraint

Extra refinement points.

default: 80

Returns:
solve_epl_shear_multithreadedcallable

Parallel solver function with signature

(theta_E, D_dt, q, phi, gamma, gamma1, gamma2)

returning a tuple of result arrays.

Examples

>>> solver = create_epl_shear_solver()
>>> results = solver(theta_E, D_dt, q, phi, gamma, gamma1, gamma2)