ler.image_properties.cross_section_njit
Module for analytical caustic computation and cross-section evaluation of EPL (Elliptical Power-Law) + external shear lens models.
Provides numba-accelerated routines for computing critical curves, caustic
boundaries, polygon areas, and lensing cross sections. All core functions
are decorated with @njit for high performance.
- Usage:
Compute the double-caustic boundary for a single lens:
>>> from ler.image_properties.cross_section_njit import caustics_epl_shear >>> pts = caustics_epl_shear(q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01)
Copyright (C) 2026 Phurailatpam Hemantakumar. Distributed under MIT License.
Module Contents
Functions
|
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. |
Attributes
- ler.image_properties.cross_section_njit.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.cross_section_njit.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.cross_section_njit.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.cross_section_njit.omega_scalar(phi, t, q, niter_max=200, tol=1e-16)[source]
Scalar version of
omegathat avoids temporary array allocation.
- ler.image_properties.cross_section_njit.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.cross_section_njit.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.cross_section_njit.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.cross_section_njit.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.cross_section_njit.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.cross_section_njit.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.cross_section_njit.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_njit.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