:py:mod:`ler.image_properties.cross_section_njit` ================================================= .. py:module:: ler.image_properties.cross_section_njit .. autoapi-nested-parse:: 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. .. !! processed by numpydoc !! Module Contents --------------- Functions ~~~~~~~~~ .. autoapisummary:: ler.image_properties.cross_section_njit.phi_q2_ellipticity ler.image_properties.cross_section_njit.pol_to_ell ler.image_properties.cross_section_njit.omega ler.image_properties.cross_section_njit.omega_scalar ler.image_properties.cross_section_njit.cdot ler.image_properties.cross_section_njit.pol_to_cart ler.image_properties.cross_section_njit.cart_to_pol ler.image_properties.cross_section_njit.caustics_epl_shear ler.image_properties.cross_section_njit.polygon_area ler.image_properties.cross_section_njit.caustic_double_area ler.image_properties.cross_section_njit.make_cross_section_reinit ler.image_properties.cross_section_njit.cross_section_epl_shear_unit Attributes ~~~~~~~~~~ .. autoapisummary:: ler.image_properties.cross_section_njit.C_LIGHT .. py:data:: C_LIGHT :value: '299792458.0' .. py:function:: phi_q2_ellipticity(phi, q) 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``. :Returns: **e1** : ``float`` First ellipticity component. **e2** : ``float`` Second ellipticity component. .. rubric:: Examples >>> e1, e2 = phi_q2_ellipticity(phi=0.25, q=0.8) .. !! processed by numpydoc !! .. py:function:: pol_to_ell(r, theta, q) Convert polar coordinates to elliptical coordinates. :Parameters: **r** : ``float`` or ``numpy.ndarray`` Radial coordinate. **theta** : ``float`` or ``numpy.ndarray`` Polar angle in radians. **q** : ``float`` Axis ratio. :Returns: **rell** : ``float`` or ``numpy.ndarray`` Elliptical radial coordinate. **phi** : ``float`` or ``numpy.ndarray`` Elliptical angle in radians. .. rubric:: Examples >>> rell, phi = pol_to_ell(1.0, 0.5, 0.8) .. !! processed by numpydoc !! .. py:function:: omega(phi, t, q, omegas, niter_max=200, tol=1e-16) Evaluate the complex angular function Omega for the EPL profile. This series expansion converges geometrically with ratio ``f = (1 - q)/(1 + q)``. The ``fastmath`` flag 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 :Returns: **omegas** : ``numpy.ndarray`` Complex Omega values at each angle (same object as input buffer). .. rubric:: 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) .. !! processed by numpydoc !! .. py:function:: omega_scalar(phi, t, q, niter_max=200, tol=1e-16) Scalar version of ``omega`` that avoids temporary array allocation. .. !! processed by numpydoc !! .. py:function:: cdot(a, b) 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. :Returns: **result** : ``float`` Real-valued dot product. .. rubric:: Examples >>> cdot(1+2j, 3+4j) 11.0 .. !! processed by numpydoc !! .. py:function:: pol_to_cart(r, th) Convert polar coordinates to Cartesian coordinates. :Parameters: **r** : ``float`` or ``numpy.ndarray`` Radial coordinate. **th** : ``float`` or ``numpy.ndarray`` Polar angle in radians. :Returns: **x** : ``float`` or ``numpy.ndarray`` Cartesian x-coordinate. **y** : ``float`` or ``numpy.ndarray`` Cartesian y-coordinate. .. rubric:: Examples >>> x, y = pol_to_cart(1.0, np.pi / 4) .. !! processed by numpydoc !! .. py:function:: cart_to_pol(x, y) Convert Cartesian coordinates to polar coordinates. The returned angle is wrapped to ``[0, 2π)``. :Parameters: **x** : ``float`` or ``numpy.ndarray`` Cartesian x-coordinate. **y** : ``float`` or ``numpy.ndarray`` Cartesian y-coordinate. :Returns: **r** : ``float`` or ``numpy.ndarray`` Radial coordinate. **theta** : ``float`` or ``numpy.ndarray`` Polar angle in radians, wrapped to ``[0, 2π)``. .. rubric:: Examples >>> r, theta = cart_to_pol(1.0, 1.0) .. !! processed by numpydoc !! .. py:function:: caustics_epl_shear(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, sourceplane=True, return_which='double') Analytically compute the caustics of an EPL + external shear lens model. For ``gamma > 2`` the outer critical curve does not exist. In that case the routine finds the curve at a finite magnification ``maginf`` instead 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' :Returns: **pts** : ``numpy.ndarray`` 2×N array of (x, y) boundary coordinates. .. rubric:: Examples >>> pts = caustics_epl_shear( ... q=0.8, phi=0.0, gamma=2.0, gamma1=0.03, gamma2=-0.01 ... ) .. !! processed by numpydoc !! .. py:function:: polygon_area(xv, yv) 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. :Returns: **area** : ``float`` The enclosed geometric area. .. rubric:: Examples >>> area = polygon_area(np.array([0.0, 1.0, 0.0]), np.array([0.0, 0.0, 1.0])) .. !! processed by numpydoc !! .. py:function:: caustic_double_area(q, phi, gamma, gamma1, gamma2, theta_E=1.0, num_th=500, maginf=-100.0, sourceplane=True, return_which='double') 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 :Returns: **area** : ``float`` Area of the double caustic in units of theta_E^2. .. rubric:: Examples >>> area = caustic_double_area(0.8, 0.0, 2.0, 0.03, -0.01) .. !! processed by numpydoc !! .. py:function:: make_cross_section_reinit(Da_instance, num_th=500, maginf=-100.0, sourceplane=True, return_which='double') 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 :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. .. rubric:: Examples >>> cross_section = make_cross_section_reinit(Da_instance) >>> cs = cross_section(zs, zl, sigma, q, phi, gamma, gamma1, gamma2) .. !! processed by numpydoc !! .. py:function:: cross_section_epl_shear_unit(e1, e2, gamma, gamma1, gamma2, theta_E=None, num_th=500, maginf=-100.0, sourceplane=True, return_which='double') 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.ndarray`` or ``None`` Einstein radius in angular units. If None, assumed to be 1.0 for all **Returns** .. **-------** .. **cs** : ``numpy.ndarray`` Cross section values in angular units. .. !! processed by numpydoc !!