caustics.lenses package

Submodules

caustics.lenses.base module

class caustics.lenses.base.ThickLens(cosmology: Cosmology, name: str | None = None)[source]

Bases: Lens

Base class for modeling gravitational lenses that cannot be treated using the thin lens approximation. It is an abstract class and should be subclassed for different types of lens models.

cosmology

An instance of a Cosmology class that describes the cosmological parameters of the model.

Type:

Cosmology

effective_convergence_curl(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Use the curl of the effective reduced deflection angle vector field to compute an effective convergence which derrives specifically from the curl of the deflection field. This field is purely a result of multiplane lensing and cannot occur in single plane lensing. See: https://arxiv.org/pdf/2006.07383.pdf

effective_convergence_div(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Using the divergence of the effective reduced delfection angle we can compute the divergence component of the effective convergence field. This field produces a single plane convergence field which reproduces as much of the deflection field as possible for a single plane. See: https://arxiv.org/pdf/2006.07383.pdf see also the effective_convergence_curl method.

effective_reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

ThickLens objects do not have a reduced deflection angle since the distance D_ls is undefined. Instead we define an effective reduced deflection angle by simply assuming the relation $lpha = heta - eta$ holds, where $lpha$ is the effective reduced deflection angle, $ heta$ are the observed angular coordinates, and $eta$ are the angular coordinates to the source plane.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

jacobian_effective_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, method='autograd', pixelscale=None, **kwargs) tuple[tuple[Tensor, Tensor], tuple[Tensor, Tensor]][source]

Return the jacobian of the effective reduced deflection angle vector field. This equates to a (2,2) matrix at each (x,y) point.

method: autograd or fft

physical_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Physical deflection angles are computed with respect to a lensing plane. ThickLens objects have no unique definition of a lens plane and so cannot compute a physical_deflection_angle

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Tuple of Tensors representing the x and y components of the deflection angle, respectively.

Return type:

tuple[Tensor, Tensor]

abstract raytrace(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Performs ray tracing by computing the angular position on the source plance associated with a given input observed angular coordinate x,y.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Tuple of Tensors representing the x and y coordinates of the ray-traced light rays, respectively.

Return type:

tuple[Tensor, Tensor]

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

ThickLens objects do not have a reduced deflection angle since the distance D_ls is undefined

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Raises:

NotImplementedError

abstract surface_density(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Computes the projected mass density at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

The projected mass density at the given coordinates in units of solar masses per square Megaparsec.

Return type:

Tensor

abstract time_delay(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Computes the gravitational time delay at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor ofsource redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

The gravitational time delay at the given coordinates.

Return type:

Tensor

class caustics.lenses.base.ThinLens(cosmology: Cosmology, z_l: Tensor | float | None = None, name: str | None = None)[source]

Bases: Lens

Base class for thin gravitational lenses.

This class provides an interface for thin gravitational lenses, i.e., lenses that can be modeled using the thin lens approximation. The class provides methods to compute several lensing quantities such as the deflection angle, convergence, potential, surface mass density, and gravitational time delay.

Parameters:
  • name (str) – Name of the lens model.

  • cosmology (Cosmology) – Cosmology object that encapsulates cosmological parameters and distances.

  • z_l (Optional[Tensor], optional) – Redshift of the lens. Defaults to None.

abstract convergence(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Computes the convergence of the lens at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Convergence at the given coordinates.

Return type:

Tensor

jacobian_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, method='autograd', pixelscale=None, **kwargs) tuple[tuple[Tensor, Tensor], tuple[Tensor, Tensor]][source]

Return the jacobian of the deflection angle vector. This equates to a (2,2) matrix at each (x,y) point.

method: autograd or fft

physical_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Computes the physical deflection angle immediately after passing through this lens’s plane.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Physical deflection angle in x and y directions in arcseconds.

Return type:

tuple[Tensor, Tensor]

abstract potential(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Computes the gravitational lensing potential at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns: Tensor: Gravitational lensing potential at the given coordinates in arcsec^2.

raytrace(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Perform a ray-tracing operation by subtracting the deflection angles from the input coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Ray-traced coordinates in the x and y directions.

Return type:

tuple[Tensor, Tensor]

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Computes the reduced deflection angle of the lens at given coordinates [arcsec].

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Reduced deflection angle in x and y directions.

Return type:

tuple[Tensor, Tensor]

surface_density(x: Tensor, y: Tensor, z_s: Tensor, z_l, *args, params: Packed | None = None, **kwargs) Tensor[source]

Computes the surface mass density of the lens at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Surface mass density at the given coordinates in solar masses per Mpc^2.

Return type:

Tensor

time_delay(x: Tensor, y: Tensor, z_s: Tensor, z_l, *args, params: Packed | None = None, **kwargs)[source]

Compute the gravitational time delay for light passing through the lens at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Time delay at the given coordinates.

Return type:

Tensor

caustics.lenses.epl module

class caustics.lenses.epl.EPL(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, q: Tensor | float | None = None, phi: Tensor | float | None = None, b: Tensor | float | None = None, t: Tensor | float | None = None, s: float = 0.0, n_iter: int = 18, name: str | None = None)[source]

Bases: ThinLens

Elliptical power law (EPL, aka singular power-law ellipsoid) profile.

This class represents a thin gravitational lens model with an elliptical power law profile. The lensing equations are solved iteratively using an approach based on Tessore et al. 2015.

n_iter

Number of iterations for the iterative solver.

Type:

int

s

Softening length for the elliptical power-law profile.

Type:

float

Parameters:
  • z_l (Optional[Union[Tensor, float]]) – This is the redshift of the lens. In the context of gravitational lensing, the lens is the galaxy or other mass distribution that is bending the light from a more distant source.

  • y0 (x0 and) – These are the coordinates of the lens center in the lens plane. The lens plane is the plane perpendicular to the line of sight in which the deflection of light by the lens is considered.

  • q (Optional[Union[Tensor, float]]) – This is the axis ratio of the lens, i.e., the ratio of the minor axis to the major axis of the elliptical lens.

  • phi (Optional[Union[Tensor, float]]) – This is the orientation of the lens on the sky, typically given as an angle measured counter-clockwise from some reference direction.

  • b (Optional[Union[Tensor, float]]) – This is the scale length of the lens, which sets the overall scale of the lensing effect. In some contexts, this is referred to as the Einstein radius.

  • t (Optional[Union[Tensor, float]]) – This is the power-law slope parameter of the lens model. In the context of the EPL model, t is equivalent to the gamma parameter minus one, where gamma is the power-law index of the radial mass distribution of the lens.

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, q, phi, b, t, *args, params: Packed | None = None, **kwargs)[source]

Compute the convergence of the lens, which describes the local density of the lens.

Parameters:
  • x (Tensor) – X coordinates in the lens plane.

  • y (Tensor) – Y coordinates in the lens plane.

  • z_s (Tensor) – Source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model.

Returns:

The convergence of the lens.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, q, phi, b, t, *args, params: Packed | None = None, **kwargs)[source]

Compute the lensing potential of the lens.

Parameters:
  • x (Tensor) – X coordinates in the lens plane.

  • y (Tensor) – Y coordinates in the lens plane.

  • z_s (Tensor) – Source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, q, phi, b, t, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Compute the reduced deflection angles of the lens.

Parameters:
  • x (Tensor) – X coordinates in the lens plane.

  • y (Tensor) – Y coordinates in the lens plane.

  • z_s (Tensor) – Source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model.

Returns:

Reduced deflection angles in the x and y directions.

Return type:

tuple[Tensor, Tensor]

caustics.lenses.external_shear module

class caustics.lenses.external_shear.ExternalShear(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, gamma_1: Tensor | float | None = None, gamma_2: Tensor | float | None = None, s: float = 0.0, name: str | None = None)[source]

Bases: ThinLens

Represents an external shear effect in a gravitational lensing system.

name

Identifier for the lens instance.

Type:

str

cosmology

The cosmological model used for lensing calculations.

Type:

Cosmology

z_l

The redshift of the lens.

Type:

Optional[Union[Tensor, float]]

x0, y0

Coordinates of the shear center in the lens plane.

Type:

Optional[Union[Tensor, float]]

gamma_1, gamma_2

Shear components.

Type:

Optional[Union[Tensor, float]]

Note: The shear components gamma_1 and gamma_2 represent an external shear, a gravitational distortion that can be caused by nearby structures outside of the main lens galaxy.

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, gamma_1, gamma_2, *args, params: Packed | None = None, **kwargs) Tensor[source]

The convergence is undefined for an external shear.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Raises:
  • NotImplementedError – This method is not implemented as the convergence is not defined

  • for an external shear.

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, gamma_1, gamma_2, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculates the lensing potential.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, gamma_1, gamma_2, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Calculates the reduced deflection angle.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The reduced deflection angles in the x and y directions.

Return type:

tuple[Tensor, Tensor]

caustics.lenses.mass_sheet module

class caustics.lenses.mass_sheet.MassSheet(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, surface_density: Tensor | float | None = None, name: str | None = None)[source]

Bases: ThinLens

Represents an external shear effect in a gravitational lensing system.

name

Identifier for the lens instance.

Type:

str

cosmology

The cosmological model used for lensing calculations.

Type:

Cosmology

z_l

The redshift of the lens.

Type:

Optional[Union[Tensor, float]]

x0, y0

Coordinates of the shear center in the lens plane.

Type:

Optional[Union[Tensor, float]]

gamma_1, gamma_2

Shear components.

Type:

Optional[Union[Tensor, float]]

Note: The shear components gamma_1 and gamma_2 represent an external shear, a gravitational distortion that can be caused by nearby structures outside of the main lens galaxy.

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, surface_density, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Computes the convergence of the lens at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns:

Convergence at the given coordinates.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, surface_density, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Computes the gravitational lensing potential at given coordinates.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

Returns: Tensor: Gravitational lensing potential at the given coordinates in arcsec^2.

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, surface_density, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Calculates the reduced deflection angle.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The reduced deflection angles in the x and y directions.

Return type:

tuple[Tensor, Tensor]

caustics.lenses.multiplane module

class caustics.lenses.multiplane.Multiplane(cosmology: Cosmology, lenses: list[ThinLens], name: str | None = None)[source]

Bases: ThickLens

Class for handling gravitational lensing with multiple lens planes.

lenses

List of thin lenses.

Type:

list[ThinLens]

Parameters:
  • name (str) – Name of the lens.

  • cosmology (Cosmology) – Cosmological parameters used for calculations.

  • lenses (list[ThinLens]) – List of thin lenses.

effective_reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

ThickLens objects do not have a reduced deflection angle since the distance D_ls is undefined. Instead we define an effective reduced deflection angle by simply assuming the relation $lpha = heta - eta$ holds, where $lpha$ is the effective reduced deflection angle, $ heta$ are the observed angular coordinates, and $eta$ are the angular coordinates to the source plane.

Parameters:
  • x (Tensor) – Tensor of x coordinates in the lens plane.

  • y (Tensor) – Tensor of y coordinates in the lens plane.

  • z_s (Tensor) – Tensor of source redshifts.

  • params (Packed, optional) – Dynamic parameter container for the lens model. Defaults to None.

get_z_ls(*args, params: Packed | None = None, **kwargs) list[Tensor][source]

Get the redshifts of each lens in the multiplane.

Parameters:

params (Packed, optional) – Dynamic parameter container.

Returns:

Redshifts of the lenses.

Return type:

List[Tensor]

raytrace(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]
Calculate the angular source positions corresponding to the

observer positions x,y. See Margarita et al. 2013 for the formalism from the GLAMER -II code: https://ui.adsabs.harvard.edu/abs/2014MNRAS.445.1954P/abstract

The primary equation used here is equation 18. With a slight correction it reads:

\[\]

ec{x}^{i+1} = ec{x}^i + D_{i+1,i}left[ ec{ heta} - sum_{j=1}^{i}f{lpha}^j( ec{x}^j) ight]

As an initialization we set the physical positions at the first lensing plane to be :math:`

ec{ heta}D_{1,0}` which is just propogation through regular space to the first plane. Note that :math:` ec{lpha}` is a physical deflection angle. The equation above converts straightforwardly into a recursion formula:

\[\]

ec{x}^{i+1} = ec{x}^i + D_{i+1,i} ec{ heta}^{i}

ec{ heta}^{i+1} = ec{ heta}^{i} - lpha^i( ec{x}^{i+1})

Here we set as initialization :math:`

ec{ heta}^0 = theta` the observation angular coordinates and :math:` ec{x}^0 = 0` the initial physical coordinates (i.e. the observation rays come from a point at the observer). The indexing of :math:` ec{x}^i` and :math:` ec{ heta}^i` indicates the properties at the plane \(i\), and 0 means the observer, 1 is the first lensing plane (infinitesimally after the plane since the deflection has been applied), and so on. Note that in the actual implementation we start at :math:` ec{x}^1` and :math:` ec{ heta}^0` and begin at the second step in the recursion formula.

Args:

x (Tensor): angular x-coordinates from the observer perspective. y (Tensor): angular y-coordinates from the observer perspective. z_s (Tensor): Redshifts of the sources. params (Packed, optional): Dynamic parameter container.

Returns:

tuple[Tensor, Tensor]: The reduced deflection angle.

raytrace_z1z2(x: Tensor, y: Tensor, z_start: Tensor, z_end: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Method to do multiplane ray tracing from arbitrary start/end redshift.

surface_density(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the projected mass density.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

Projected mass density [solMass / Mpc^2].

Return type:

Tensor

Raises:

NotImplementedError – This method is not yet implemented.

time_delay(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the time delay of light caused by the lensing.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

Time delay caused by the lensing.

Return type:

Tensor

Raises:

NotImplementedError – This method is not yet implemented.

caustics.lenses.nfw module

class caustics.lenses.nfw.NFW(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, m: Tensor | float | None = None, c: Tensor | float | None = None, s: float = 0.0, use_case='batchable', name: str | None = None)[source]

Bases: ThinLens

NFW lens class. This class models a lens using the Navarro-Frenk-White (NFW) profile. The NFW profile is a spatial density profile of dark matter halo that arises in cosmological simulations.

z_l

Redshift of the lens. Default is None.

Type:

Optional[Tensor]

x0

x-coordinate of the lens center in the lens plane. Default is None.

Type:

Optional[Tensor]

y0

y-coordinate of the lens center in the lens plane. Default is None.

Type:

Optional[Tensor]

m

Mass of the lens. Default is None.

Type:

Optional[Tensor]

c

Concentration parameter of the lens. Default is None.

Type:

Optional[Tensor]

s

Softening parameter to avoid singularities at the center of the lens. Default is 0.0.

Type:

float

use_case

Due to an idyosyncratic behaviour of PyTorch, the NFW/TNFW profile specifically cant be both batchable and differentiable. You may select which version you wish to use by setting this parameter to one of: batchable, differentiable.

Type:

str

get_scale_radius()[source]

Returns the scale radius of the lens.

get_scale_density()[source]

Returns the scale density of the lens.

get_convergence_s()[source]

Returns the dimensionless surface mass density of the lens.

_f()

Helper method for computing deflection angles.

_g()

Helper method for computing lensing potential.

_h()

Helper method for computing reduced deflection angles.

deflection_angle_hat()

Computes the reduced deflection angle.

deflection_angle()

Computes the deflection angle.

convergence()[source]

Computes the convergence (dimensionless surface mass density).

potential()[source]

Computes the lensing potential.

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, m, c, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the convergence (dimensionless surface mass density).

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The convergence (dimensionless surface mass density).

Return type:

Tensor

get_convergence_s(z_s, z_l, x0, y0, m, c, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the dimensionless surface mass density of the lens.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • z_s (Tensor) – Redshift of the source.

  • m (Tensor) – Mass of the lens.

  • c (Tensor) – Concentration parameter of the lens.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The dimensionless surface mass density of the lens.

Return type:

Tensor

get_scale_density(z_l, x0, y0, m, c, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the scale density of the lens.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • c (Tensor) – Concentration parameter of the lens.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The scale density of the lens in solar masses per Mpc cubed.

Return type:

Tensor

get_scale_radius(z_l, x0, y0, m, c, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the scale radius of the lens.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • m (Tensor) – Mass of the lens.

  • c (Tensor) – Concentration parameter of the lens.

  • x (dict) – Dynamic parameter container.

Returns:

The scale radius of the lens in Mpc.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, m, c, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, m, c, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Compute the reduced deflection angle.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The reduced deflection angles in the x and y directions.

Return type:

tuple[Tensor, Tensor]

caustics.lenses.pixelated_convergence module

class caustics.lenses.pixelated_convergence.PixelatedConvergence(pixelscale: float, n_pix: int, cosmology: Cosmology, z_l: Tensor | None = None, x0: Tensor | None = tensor(0.), y0: Tensor | None = tensor(0.), convergence_map: Tensor | None = None, shape: tuple[int, ...] | None = None, convolution_mode: str = 'fft', use_next_fast_len: bool = True, padding='zero', name: str | None = None)[source]

Bases: ThinLens

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, convergence_map, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the convergence at the specified positions. This method is not implemented.

Parameters:
  • x (Tensor) – The x-coordinates of the positions to compute the convergence for.

  • y (Tensor) – The y-coordinates of the positions to compute the convergence for.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – A dictionary containing additional parameters.

Returns:

The convergence at the specified positions.

Return type:

Tensor

Raises:

NotImplementedError – This method is not implemented.

property convolution_mode

Get the convolution mode of the ConvergenceGrid object.

Returns:

The convolution mode, either “fft” or “conv2d”.

Return type:

str

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, convergence_map, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential at the specified positions using the given convergence map.

Args: x (Tensor): The x-coordinates of the positions to compute the lensing potential for. y (Tensor): The y-coordinates of the positions to compute the lensing potential for. z_s (Tensor): The source redshift. params (Packed, optional): A dictionary containing additional parameters.

Returns:

The lensing potential at the specified positions.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, convergence_map, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Compute the deflection angles at the specified positions using the given convergence map.

Parameters:
  • x (Tensor) – The x-coordinates of the positions to compute the deflection angles for.

  • y (Tensor) – The y-coordinates of the positions to compute the deflection angles for.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – A dictionary containing additional parameters.

Returns:

The x and y components of the deflection angles at the specified positions.

Return type:

tuple[Tensor, Tensor]

to(device: device | None = None, dtype: dtype | None = None)[source]

Move the ConvergenceGrid object and all its tensors to the specified device and dtype.

Parameters:
  • device (Optional[torch.device]) – The target device to move the tensors to.

  • dtype (Optional[torch.dtype]) – The target data type to cast the tensors to.

caustics.lenses.point module

class caustics.lenses.point.Point(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, th_ein: Tensor | float | None = None, s: float = 0.0, name: str | None = None)[source]

Bases: ThinLens

Class representing a point mass lens in strong gravitational lensing.

name

The name of the point lens.

Type:

str

cosmology

The cosmology used for calculations.

Type:

Cosmology

z_l

Redshift of the lens.

Type:

Optional[Union[Tensor, float]]

x0

x-coordinate of the center of the lens.

Type:

Optional[Union[Tensor, float]]

y0

y-coordinate of the center of the lens.

Type:

Optional[Union[Tensor, float]]

th_ein

Einstein radius of the lens.

Type:

Optional[Union[Tensor, float]]

s

Softening parameter to prevent numerical instabilities.

Type:

float

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, th_ein, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the convergence (dimensionless surface mass density).

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The convergence (dimensionless surface mass density).

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, th_ein, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, th_ein, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Compute the deflection angles.

Parameters:
  • x (Tensor) – x-coordinates in the lens plane.

  • y (Tensor) – y-coordinates in the lens plane.

  • z_s (Tensor) – Redshifts of the sources.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The deflection angles in the x and y directions.

Return type:

tuple[Tensor, Tensor]

caustics.lenses.pseudo_jaffe module

class caustics.lenses.pseudo_jaffe.PseudoJaffe(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, mass: Tensor | float | None = None, core_radius: Tensor | float | None = None, scale_radius: Tensor | float | None = None, s: float = 0.0, name: str | None = None)[source]

Bases: ThinLens

Class representing a Pseudo Jaffe lens in strong gravitational lensing, based on Eliasdottir et al 2007 and the lenstronomy source code.

name

The name of the Pseudo Jaffe lens.

Type:

str

cosmology

The cosmology used for calculations.

Type:

Cosmology

z_l

Redshift of the lens.

Type:

Optional[Union[Tensor, float]]

x0

x-coordinate of the center of the lens (arcsec).

Type:

Optional[Union[Tensor, float]]

y0

y-coordinate of the center of the lens (arcsec).

Type:

Optional[Union[Tensor, float]]

mass

Total mass of the lens (Msol).

Type:

Optional[Union[Tensor, float]]

core_radius

Core radius of the lens (arcsec).

Type:

Optional[Union[Tensor, float]]

scale_radius

Scaling radius of the lens (arcsec).

Type:

Optional[Union[Tensor, float]]

s

Softening parameter to prevent numerical instabilities.

Type:

float

static central_convergence(z_l, z_s, rho_0, core_radius, scale_radius, critical_surface_density)[source]

Compute the central convergence.

Parameters:
  • z_l (Tensor) – Lens redshift.

  • z_s (Tensor) – Source redshift.

  • rho_0 (Tensor) – Central mass density.

  • core_radius (Tensor) – Core radius of the lens (must be in Mpc).

  • scale_radius (Tensor) – Scaling radius of the lens (must be in Mpc).

  • cosmology (Cosmology) – The cosmology used for calculations.

Returns:

The central convergence.

Return type:

Tensor

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, mass, core_radius, scale_radius, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the projected mass density, based on equation A6.

Parameters:
  • x (Tensor) – x-coordinate of the lens.

  • y (Tensor) – y-coordinate of the lens.

  • z_s (Tensor) – Source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The projected mass density.

Return type:

Tensor

get_convergence_0(z_s, z_l, x0, y0, mass, core_radius, scale_radius, *args, params: Packed | None = None, **kwargs)[source]
mass_enclosed_2d(theta, z_s, z_l, x0, y0, mass, core_radius, scale_radius, *args, params: Packed | None = None, **kwargs)[source]

Calculate the mass enclosed within a two-dimensional radius.

Parameters:
  • theta (Tensor) – Radius at which to calculate enclosed mass (arcsec).

  • z_s (Tensor) – Source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The mass enclosed within the given radius.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, mass, core_radius, scale_radius, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential. This calculation is based on equation A18.

Parameters:
  • x (Tensor) – x-coordinate of the lens.

  • y (Tensor) – y-coordinate of the lens.

  • z_s (Tensor) – Source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, mass, core_radius, scale_radius, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Calculate the deflection angle.

Parameters:
  • x (Tensor) – x-coordinate of the lens.

  • y (Tensor) – y-coordinate of the lens.

  • z_s (Tensor) – Source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The deflection angle in the x and y directions.

Return type:

Tuple[Tensor, Tensor]

caustics.lenses.sie module

class caustics.lenses.sie.SIE(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, q: Tensor | float | None = None, phi: Tensor | float | None = None, b: Tensor | float | None = None, s: float = 0.0, name: str | None = None)[source]

Bases: ThinLens

A class representing a Singular Isothermal Ellipsoid (SIE) strong gravitational lens model. This model is based on Keeton 2001, which can be found at https://arxiv.org/abs/astro-ph/0102341.

name

The name of the lens.

Type:

str

cosmology

An instance of the Cosmology class.

Type:

Cosmology

z_l

The redshift of the lens.

Type:

Optional[Union[Tensor, float]]

x0

The x-coordinate of the lens center.

Type:

Optional[Union[Tensor, float]]

y0

The y-coordinate of the lens center.

Type:

Optional[Union[Tensor, float]]

q

The axis ratio of the lens.

Type:

Optional[Union[Tensor, float]]

phi

The orientation angle of the lens (position angle).

Type:

Optional[Union[Tensor, float]]

b

The Einstein radius of the lens.

Type:

Optional[Union[Tensor, float]]

s

The core radius of the lens (defaults to 0.0).

Type:

float

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, q, phi, b, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the projected mass density.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The projected mass.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, q, phi, b, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, q, phi, b, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Calculate the physical deflection angle.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The deflection angle in the x and y directions.

Return type:

Tuple[Tensor, Tensor]

caustics.lenses.singleplane module

class caustics.lenses.singleplane.SinglePlane(cosmology: Cosmology, lenses: list[ThinLens], name: str | None = None, **kwargs)[source]

Bases: ThinLens

A class for combining multiple thin lenses into a single lensing plane. This model inherits from the base ThinLens class.

name

The name of the single plane lens.

Type:

str

cosmology

An instance of the Cosmology class.

Type:

Cosmology

lenses

A list of ThinLens objects that are being combined into a single lensing plane.

Type:

List[ThinLens]

convergence(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the total projected mass density by summing the mass densities of all individual lenses.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The total projected mass density.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the total lensing potential by summing the lensing potentials of all individual lenses.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The total lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Calculate the total deflection angle by summing the deflection angles of all individual lenses.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The total deflection angle in the x and y directions.

Return type:

Tuple[Tensor, Tensor]

caustics.lenses.sis module

class caustics.lenses.sis.SIS(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, th_ein: Tensor | float | None = None, s: float = 0.0, name: str | None = None)[source]

Bases: ThinLens

A class representing the Singular Isothermal Sphere (SIS) model. This model inherits from the base ThinLens class.

name

The name of the SIS lens.

Type:

str

cosmology

An instance of the Cosmology class.

Type:

Cosmology

z_l

The lens redshift.

Type:

Optional[Union[Tensor, float]]

x0

The x-coordinate of the lens center.

Type:

Optional[Union[Tensor, float]]

y0

The y-coordinate of the lens center.

Type:

Optional[Union[Tensor, float]]

th_ein

The Einstein radius of the lens.

Type:

Optional[Union[Tensor, float]]

s

A smoothing factor, default is 0.0.

Type:

float

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, th_ein, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the projected mass density of the SIS lens.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The projected mass density.

Return type:

Tensor

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, th_ein, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential of the SIS lens.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

reduced_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, th_ein, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Calculate the deflection angle of the SIS lens.

Parameters:
  • x (Tensor) – The x-coordinate of the lens.

  • y (Tensor) – The y-coordinate of the lens.

  • z_s (Tensor) – The source redshift.

  • params (Packed, optional) – Dynamic parameter container.

Returns:

The deflection angle in the x and y directions.

Return type:

Tuple[Tensor, Tensor]

caustics.lenses.tnfw module

class caustics.lenses.tnfw.TNFW(cosmology: Cosmology, z_l: Tensor | float | None = None, x0: Tensor | float | None = None, y0: Tensor | float | None = None, mass: Tensor | float | None = None, scale_radius: Tensor | float | None = None, tau: Tensor | float | None = None, s: float = 0.0, interpret_m_total_mass: bool = True, use_case='batchable', name: str | None = None)[source]

Bases: ThinLens

Truncated Navaro-Frenk-White profile

TNFW lens class. This class models a lens using the truncated Navarro-Frenk-White (NFW) profile. The NFW profile is a spatial density profile of dark matter halo that arises in cosmological simulations. It is truncated with an extra scaling term which smoothly reduces the density such that it does not diverge to infinity. This is based off the paper by Baltz et al. 2009:

https://arxiv.org/abs/0705.0682

https://ui.adsabs.harvard.edu/abs/2009JCAP…01..015B/abstract

Note

The mass m in the TNFW profile corresponds to the total mass of the lens. This is different from the NFW profile where the mass m parameter corresponds to the mass within R200. If you prefer the “mass within R200” version you can set: interpret_m_total_mass = False on initialization of the object. However, the mass within R200 will be computed for an NFW profile, not a TNFW profile. This is in line with how lenstronomy inteprets the mass parameter.

Parameters:
  • name (str) – Name of the lens instance.

  • cosmology (Cosmology) – An instance of the Cosmology class which contains information about the cosmological model and parameters.

  • z_l (Optional[Tensor]) – Redshift of the lens.

  • x0 (Optional[Tensor]) – Center of lens position on x-axis (arcsec).

  • y0 (Optional[Tensor]) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • s (float) – Softening parameter to avoid singularities at the center of the lens. Default is 0.0.

  • interpret_m_total_mass (bool) – Indicates how to interpret the mass variable “m”. If true the mass is intepreted as the total mass of the halo (good because it makes sense). If false it is intepreted as what the mass would have been within R200 of a an NFW that isn’t truncated (good because it is easily compared with an NFW).

  • use_case (str) – Due to an idyosyncratic behaviour of PyTorch, the NFW/TNFW profile specifically cant be both batchable and differentiable. You may select which version you wish to use by setting this parameter to one of: batchable, differentiable.

convergence(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

TNFW convergence as given in Baltz et al. 2009. This is unitless since it is Sigma(x) / Sigma_crit.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

unitless convergence at requested position

Return type:

Tensor

get_M0(z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the reference mass. This is an abstract reference mass used internally in the equations from Baltz et al. 2009.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

The reference mass of the lens in Msol.

Return type:

Tensor

get_concentration(z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the scale radius of the lens. This is the same formula used for the classic NFW profile.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

The scale radius of the lens in Mpc.

Return type:

Tensor

get_scale_density(z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the scale density of the lens.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

The scale density of the lens in solar masses per Mpc cubed.

Return type:

Tensor

get_truncation_radius(z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

Calculate the truncation radius of the TNFW lens.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

The truncation radius of the lens in arcsec.

Return type:

Tensor

mass_enclosed_2d(r: Tensor, z_s: Tensor, z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

Total projected mass (Msol) within a radius r (arcsec).

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

Integrated mass projected in infinite cylinder within radius r.

Return type:

Tensor

physical_deflection_angle(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) tuple[Tensor, Tensor][source]

Compute the physical deflection angle (arcsec) for this lens at the requested position. Note that the NFW/TNFW profile is more naturally represented as a physical deflection angle, this is easily internally converted to a reduced deflection angle.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

The physical deflection angles in the x and y directions (arcsec).

Return type:

tuple[Tensor, Tensor]

potential(x: Tensor, y: Tensor, z_s: Tensor, z_l, x0, y0, mass, scale_radius, tau, *args, params: Packed | None = None, **kwargs) Tensor[source]

Compute the lensing potential. Note that this is not a unitless potential! This is the potential as given in Baltz et al. 2009.

TODO: convert to dimensionless potential.

Parameters:
  • z_l (Tensor) – Redshift of the lens.

  • x0 (Tensor) – Center of lens position on x-axis (arcsec).

  • y0 (Tensor) – Center of lens position on y-axis (arcsec).

  • mass (Optional[Tensor]) – Mass of the lens (Msol).

  • scale_radius (Optional[Tensor]) – Scale radius of the TNFW lens (arcsec).

  • tau (Optional[Tensor]) – Truncation scale. Ratio of truncation radius to scale radius (rt/rs).

  • params (dict) – Dynamic parameter container.

Returns:

The lensing potential.

Return type:

Tensor

caustics.lenses.utils module

caustics.lenses.utils.get_magnification(raytrace, x, y, z_s) Tensor[source]

Computes the magnification over a grid on the lensing plane. This is done by calling get_pix_magnification for each point on the grid.

Parameters:
  • raytrace – A function that maps the lensing plane coordinates to the source plane coordinates.

  • x (Tensor) – The x-coordinates on the lensing plane.

  • y (Tensor) – The y-coordinates on the lensing plane.

  • z_s (Tensor) – The redshift of the source.

Returns:

A tensor representing the magnification at each point on the grid.

caustics.lenses.utils.get_pix_jacobian(raytrace, x, y, z_s) Tuple[Tuple[Tensor, Tensor], Tuple[Tensor, Tensor]][source]

Computes the Jacobian matrix of the partial derivatives of the image position with respect to the source position (\(\partial eta / \partial heta\)). This is done at a single point on the lensing plane.

Parameters:
  • raytrace – A function that maps the lensing plane coordinates to the source plane coordinates.

  • x (Tensor) – The x-coordinate on the lensing plane.

  • y (Tensor) – The y-coordinate on the lensing plane.

  • z_s (Tensor) – The redshift of the source.

Returns:

The Jacobian matrix of the image position with respect to the source position at the given point.

caustics.lenses.utils.get_pix_magnification(raytrace, x, y, z_s) Tensor[source]

Computes the magnification at a single point on the lensing plane. The magnification is derived from the determinant of the Jacobian matrix of the image position with respect to the source position.

Parameters:
  • raytrace – A function that maps the lensing plane coordinates to the source plane coordinates.

  • x (Tensor) – The x-coordinate on the lensing plane.

  • y (Tensor) – The y-coordinate on the lensing plane.

  • z_s (Tensor) – The redshift of the source.

Returns:

The magnification at the given point on the lensing plane.

Module contents