lenstronomy.LensModel package

Subpackages

Submodules

lenstronomy.LensModel.convergence_integrals module

potential_from_kappa_grid(kappa, grid_spacing)[source]

Lensing potential \(\psi ({\vec {\theta }})\) on the convergence grid \(\kappa\).

\[\psi ({\vec {\theta }})={\frac {1}{\pi }}\int d^{2}\theta ^{\prime } \kappa ({\vec {\theta }}^{\prime })\ln |{\vec {\theta }}-{\vec {\theta }}^{\prime }|\]

The computation is performed as a convolution of the Green’s function with the convergence map using FFT.

Parameters:
  • kappa – 2d grid of convergence values

  • grid_spacing – scale of an individual pixel (per axis) of grid

Returns:

lensing potential in a 2d grid at positions x_grid, y_grid

potential_from_kappa_grid_adaptive(kappa_high_res, grid_spacing, low_res_factor, high_res_kernel_size)[source]

Lensing potential on the convergence grid the computation is performed as a convolution of the Green’s function with the convergence map using FFT.

Parameters:
  • kappa_high_res – 2d grid of convergence values

  • grid_spacing – scale of an individual pixel (per axis) of grid

  • low_res_factor – lower resolution factor of larger scale kernel.

  • high_res_kernel_size – int, size of high resolution kernel in units of degraded pixels

Returns:

lensing potential in a 2d grid at positions x_grid, y_grid

deflection_from_kappa_grid(kappa, grid_spacing)[source]

Deflection angle \(\vec {\alpha }}\) from a convergence grid \(\kappa\).

\[{\vec {\alpha }}({\vec {\theta }})={\frac {1}{\pi }} \int d^{2}\theta ^{\prime }{\frac {({\vec {\theta }}-{\vec {\theta }}^{\prime }) \kappa ({\vec {\theta }}^{\prime })}{|{\vec {\theta }}-{\vec {\theta }}^{\prime }|^{2}}}\]

The computation is performed as a convolution of the Green’s function with the convergence map using FFT.

Parameters:
  • kappa – convergence values for each pixel (2-d array)

  • grid_spacing – scale of an individual pixel (per axis) of grid

Returns:

numerical deflection angles in x- and y- direction over the convergence grid points

deflection_from_kappa_grid_adaptive(kappa_high_res, grid_spacing, low_res_factor, high_res_kernel_size)[source]

Deflection angles on the convergence grid with adaptive FFT the computation is performed as a convolution of the Green’s function with the convergence map using FFT The grid is returned in the lower resolution grid.

Parameters:
  • kappa_high_res – convergence values for each pixel (2-d array)

  • grid_spacing – pixel size of high resolution grid

  • low_res_factor – lower resolution factor of larger scale kernel.

  • high_res_kernel_size – int, size of high resolution kernel in units of degraded pixels

Returns:

numerical deflection angles in x- and y- direction

potential_kernel(num_pix, delta_pix)[source]

Numerical gridded integration kernel for convergence to lensing kernel with given pixel size.

Parameters:
  • num_pix – integer; number of pixels of kernel per axis

  • delta_pix – pixel size (per dimension in units of angle)

Returns:

kernel for lensing potential

deflection_kernel(num_pix, delta_pix)[source]

Numerical gridded integration kernel for convergence to deflection angle with given pixel size.

Parameters:
  • num_pix – integer; number of pixels of kernel per axis, should be odd number to have a defined center

  • delta_pix – pixel size (per dimension in units of angle)

Returns:

kernel for x-direction and kernel of y-direction deflection angles

lenstronomy.LensModel.lens_model module

class LensModel(lens_model_list, z_lens=None, z_source=None, lens_redshift_list=None, cosmo=None, multi_plane=False, numerical_alpha_class=None, observed_convention_index=None, z_source_convention=None, cosmo_interp=False, z_interp_stop=None, num_z_interp=100, kwargs_interp=None, kwargs_synthesis=None, distance_ratio_sampling=False)[source]

Bases: object

Class to handle an arbitrary list of lens models.

This is the main lenstronomy LensModel API for all other modules.

__init__(lens_model_list, z_lens=None, z_source=None, lens_redshift_list=None, cosmo=None, multi_plane=False, numerical_alpha_class=None, observed_convention_index=None, z_source_convention=None, cosmo_interp=False, z_interp_stop=None, num_z_interp=100, kwargs_interp=None, kwargs_synthesis=None, distance_ratio_sampling=False)[source]
Parameters:
  • lens_model_list – list of strings with lens model names

  • z_lens – redshift of the deflector (only considered when operating in single plane mode). Is only needed for specific functions that require a cosmology.

  • z_source – redshift of the source: Needed in multi_plane option only, not required for the core functionalities in the single plane mode.

  • lens_redshift_list – list of deflector redshift (corresponding to the lens model list), only applicable in multi_plane mode.

  • cosmo – instance of the astropy cosmology class. If not specified, uses the default cosmology.

  • multi_plane – bool, if True, uses multi-plane mode. Default is False.

  • numerical_alpha_class – an instance of a custom class for use in TabulatedDeflections() lens model (see documentation in Profiles/numerical_deflections)

  • kwargs_interp – interpolation keyword arguments specifying the numerics. See description in the Interpolate() class. Only applicable for ‘INTERPOL’ and ‘INTERPOL_SCALED’ models.

  • observed_convention_index – a list of indices, corresponding to the lens_model_list element with same index, where the ‘center_x’ and ‘center_y’ kwargs correspond to observed (lensed) positions, not physical positions. The code will compute the physical locations when performing computations

  • z_source_convention – float, redshift of a source to define the reduced deflection angles of the lens models. If None, ‘z_source’ is used.

  • cosmo_interp – boolean (only employed in multi-plane mode), interpolates astropy.cosmology distances for faster calls when accessing several lensing planes

  • z_interp_stop – (only in multi-plane with cosmo_interp=True); maximum redshift for distance interpolation This number should be higher or equal the maximum of the source redshift and/or the z_source_convention

  • num_z_interp – (only in multi-plane with cosmo_interp=True); number of redshift bins for interpolating

distances :param distance_ratio_sampling: bool, if True, will use sampled distance ratios to update T_ij value in multi-lens plane computation.

ray_shooting(x, y, kwargs, k=None)[source]

Maps image to source position (inverse deflection)

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

Returns:

source plane positions corresponding to (x, y) in the image plane

fermat_potential(x_image, y_image, kwargs_lens, x_source=None, y_source=None)[source]

Fermat potential (negative sign means earlier arrival time) for Multi-plane lensing, it computes the effective Fermat potential (derived from the arrival time and subtracted off the time-delay distance for the given cosmology). The units are given in arcsecond square.

Parameters:
  • x_image – image position

  • y_image – image position

  • x_source – source position

  • y_source – source position

  • kwargs_lens – list of keyword arguments of lens model parameters matching the lens model classes

Returns:

fermat potential in arcsec**2 without geometry term (second part of Eqn 1 in Suyu et al. 2013) as a list

arrival_time(x_image, y_image, kwargs_lens, kappa_ext=0, x_source=None, y_source=None)[source]

Arrival time of images relative to a straight line without lensing. Negative values correspond to images arriving earlier, and positive signs correspond to images arriving later.

Parameters:
  • x_image – image position

  • y_image – image position

  • kwargs_lens – lens model parameter keyword argument list

  • kappa_ext – external convergence contribution not accounted in the lens model that leads to the same observables in position and relative fluxes but rescales the time delays

  • x_source – source position (optional), otherwise computed with ray-tracing

  • y_source – source position (optional), otherwise computed with ray-tracing

Returns:

arrival time of image positions in units of days

potential(x, y, kwargs, k=None)[source]

Lensing potential.

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

Returns:

lensing potential in units of arcsec^2

alpha(x, y, kwargs, k=None, diff=None)[source]

Deflection angles.

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

  • diff – None or float. If set, computes the deflection as a finite numerical differential of the lensing potential. This differential is only applicable in the single lensing plane where the form of the lensing potential is analytically known

Returns:

deflection angles in units of arcsec

hessian(x, y, kwargs, k=None, diff=None, diff_method='square')[source]

Hessian matrix.

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

  • diff – float, scale over which the finite numerical differential is computed. If None, then using the exact (if available) differentials.

  • diff_method – string, ‘square’ or ‘cross’, indicating whether finite differentials are computed from a cross or a square of points around (x, y)

Returns:

f_xx, f_xy, f_yx, f_yy components

kappa(x, y, kwargs, k=None, diff=None, diff_method='square')[source]

Lensing convergence k = 1/2 laplacian(phi)

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

  • diff – float, scale over which the finite numerical differential is computed. If None, then using the exact (if available) differentials.

  • diff_method – string, ‘square’ or ‘cross’, indicating whether finite differentials are computed from a cross or a square of points around (x, y)

Returns:

lensing convergence

curl(x, y, kwargs, k=None, diff=None, diff_method='square')[source]

curl computation F_xy - F_yx

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

  • diff – float, scale over which the finite numerical differential is computed. If None, then using the exact (if available) differentials.

  • diff_method – string, ‘square’ or ‘cross’, indicating whether finite differentials are computed from a cross or a square of points around (x, y)

Returns:

curl at position (x, y)

gamma(x, y, kwargs, k=None, diff=None, diff_method='square')[source]

shear computation g1 = 1/2(d^2phi/dx^2 - d^2phi/dy^2) g2 = d^2phi/dxdy

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

  • diff – float, scale over which the finite numerical differential is computed. If None, then using the exact (if available) differentials.

  • diff_method – string, ‘square’ or ‘cross’, indicating whether finite differentials are computed from a cross or a square of points around (x, y)

Returns:

gamma1, gamma2

magnification(x, y, kwargs, k=None, diff=None, diff_method='square')[source]

magnification.

mag = 1/det(A) A = 1 - d^2phi/d_ij

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

  • diff – float, scale over which the finite numerical differential is computed. If None, then using the exact (if available) differentials.

  • diff_method – string, ‘square’ or ‘cross’, indicating whether finite differentials are computed from a cross or a square of points around (x, y)

Returns:

magnification

flexion(x, y, kwargs, k=None, diff=1e-06, hessian_diff=True)[source]

Third derivatives (flexion)

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – int or None, if set, only evaluates the differential from one model component

  • diff – numerical differential length of Flexion

  • hessian_diff – boolean, if true also computes the numerical differential length of Hessian (optional)

Returns:

f_xxx, f_xxy, f_xyy, f_yyy

set_static(kwargs)[source]

Set this instance to a static lens model. This can improve the speed in evaluating lensing quantities at different positions but must not be used with different lens model parameters!

Parameters:

kwargs – lens model keyword argument list

Returns:

kwargs_updated (in case of image position convention in multiplane lensing this is changed)

set_dynamic()[source]

Deletes cache for static setting and makes sure the observed convention in the position of lensing profiles in the multi-plane setting is enabled. Dynamic is the default setting of this class enabling an accurate computation of lensing quantities with different parameters in the lensing profiles.

Returns:

None

hessian_z1z2(z1, z2, theta_x, theta_y, kwargs_lens, diff=1e-08)[source]

Computes Hessian matrix when Observed at z1 with rays going to z2 with z1 < z2 for multi_plane.

Parameters:
  • z1 – Observer redshift

  • z2 – source redshift

  • theta_x – angular position and direction of the ray

  • theta_y – angular position and direction of the ray

  • kwargs_lens – list of keyword arguments of lens model parameters matching the lens model classes

  • diff – numerical differential step (float)

Returns:

f_xx, f_xy, f_yx, f_yy

lenstronomy.LensModel.lens_model_extensions module

class LensModelExtensions(lensModel)[source]

Bases: object

Class with extension routines not part of the LensModel core routines.

__init__(lensModel)[source]
Parameters:

lensModel – instance of the LensModel() class, or with same functionalities. In particular, the following definitions are required to execute all functionalities presented in this class: def ray_shooting() def magnification() def kappa() def alpha() def hessian()

magnification_finite_adaptive(x_image, y_image, kwargs_lens, source_model, kwargs_source, grid_resolution, grid_radius_arcsec, axis_ratio=0.5, tol=0.001, step_size=0.05, use_largest_eigenvalue=True, fixed_aperture_size=False)[source]

This method computes image magnifications with a finite-size background source.

This new updates allows for more flexibility in the source light model by requiring the user to specify the source light mode, grid size and grid resolution before calling the function.

The functions auto_raytrracing_grid_size and auto_raytracing_grid_resolution give good estimates for appropriate parameter choices for grid_radius_arcsec and grid_resolution. It can be much faster that magnification_finite for lens models with many deflectors and a compact source. This is because most pixels in a rectangular window around a lensed image of a compact source do not map onto the source, and therefore don’t contribute to the integrated flux in the image plane.

Rather than ray tracing through a rectangular grid, this routine accelerates the computation of image magnifications with finite-size sources by ray tracing through an elliptical region oriented such that tracks the surface brightness of the lensed image. The aperture size is initially quite small, and increases in size until the flux inside of it (and hence the magnification) converges. The orientation of the elliptical aperture is computed from the magnification tensor evaluated at the image coordinate.

If for whatever reason you prefer a circular aperture to the elliptical approximation using the hessian eigenvectors, you can just set axis_ratio = 1.

To use the eigenvalues of the hessian matrix to estimate the optimum axis ratio, set axis_ratio = 0.

The default settings for the grid resolution and ray tracing window size work well for sources with fwhm between 0.5 - 100 pc.

Parameters:
  • x_image – a list or array of x coordinates [units arcsec]

  • y_image – a list or array of y coordinates [units arcsec]

  • kwargs_lens – keyword arguments for the lens model

  • source_model – instance of LightModel for the source

  • grid_resolution – the grid resolution in units arcsec/pixel; if not specified, an appropriate value will be estimated from the source size

  • grid_radius_arcsec – the size of the ray tracing region in arcsec; if not specified, an appropriate value will be estimated from the source size

  • axis_ratio – the axis ratio of the ellipse used for ray tracing; if axis_ratio = 0, then the eigenvalues the hessian matrix will be used to estimate an appropriate axis ratio. Be warned: if the image is highly magnified it will tend to curve out of the resulting ellipse

  • tol – tolerance for convergence in the magnification

  • step_size – sets the increment for the successively larger ray tracing windows

  • use_largest_eigenvalue – bool; if True, then the major axis of the ray tracing ellipse region will be aligned with the eigenvector corresponding to the largest eigenvalue of the hessian matrix

  • fixed_aperture_size – bool, if True the flux is computed inside a fixed aperture size with radius grid_radius_arcsec

Kwargs_source:

keyword arguments for the light profile of the source (list of dictionary)

Returns:

an array of image magnifications

magnification_finite(x_pos, y_pos, kwargs_lens, source_sigma=0.003, window_size=0.1, grid_number=100, polar_grid=False, aspect_ratio=0.5)[source]

Returns the magnification of an extended source with Gaussian light profile.

Parameters:
  • x_pos – x-axis positons of point sources

  • y_pos – y-axis position of point sources

  • kwargs_lens – lens model kwargs

  • source_sigma – Gaussian sigma in arc sec in source

  • window_size – size of window to compute the finite flux

  • grid_number – number of grid cells per axis in the window to numerically compute the flux

Returns:

numerically computed brightness of the sources

zoom_source(x_pos, y_pos, kwargs_lens, source_sigma=0.003, window_size=0.1, grid_number=100, shape='GAUSSIAN')[source]

Computes the surface brightness on an image with a zoomed window.

Parameters:
  • x_pos – angular coordinate of center of image

  • y_pos – angular coordinate of center of image

  • kwargs_lens – lens model parameter list

  • source_sigma – source size (in angular units)

  • window_size – window size in angular units

  • grid_number – number of grid points per axis

  • shape – string, shape of source, supports ‘GAUSSIAN’ and ‘TORUS

Returns:

2d numpy array

critical_curve_tiling(kwargs_lens, compute_window=5, start_scale=0.5, max_order=10, center_x=0, center_y=0)[source]
Parameters:
  • kwargs_lens – lens model keyword argument list

  • compute_window – total window in the image plane where to search for critical curves

  • start_scale – float, angular scale on which to start the tiling from (if there are two distinct curves in a region, it might only find one.

  • max_order – int, maximum order in the tiling to compute critical curve triangles

  • center_x – float, center of the window to compute critical curves and caustics

  • center_y – float, center of the window to compute critical curves and caustics

Returns:

list of positions representing coordinates of the critical curve (in RA and DEC)

caustic_area(kwargs_lens, kwargs_caustic_num, index_vertices=0)[source]

Computes the area inside a connected caustic curve.

Parameters:
  • kwargs_lens – lens model keyword argument list

  • kwargs_caustic_num – keyword arguments for the numerical calculation of the caustics, as input of self.critical_curve_caustics()

  • index_vertices – integer, index of connected vortex from the output of self.critical_curve_caustics() of disconnected curves.

Returns:

area within the caustic curve selected

critical_curve_caustics(kwargs_lens, compute_window=5, grid_scale=0.01, center_x=0, center_y=0)[source]
Parameters:
  • kwargs_lens – lens model kwargs

  • compute_window – window size in arcsec where the critical curve is computed

  • grid_scale – numerical grid spacing of the computation of the critical curves

  • center_x – float, center of the window to compute critical curves and caustics

  • center_y – float, center of the window to compute critical curves and caustics

Returns:

lists of ra and dec arrays corresponding to different disconnected critical curves and their caustic counterparts

hessian_eigenvectors(x, y, kwargs_lens, diff=None)[source]

Computes magnification eigenvectors at position (x, y)

Parameters:
  • x – x-position

  • y – y-position

  • kwargs_lens – lens model keyword arguments

Returns:

radial stretch, tangential stretch

radial_tangential_stretch(x, y, kwargs_lens, diff=None, ra_0=0, dec_0=0, coordinate_frame_definitions=False)[source]

Computes the radial and tangential stretches at a given position.

Parameters:
  • x – x-position

  • y – y-position

  • kwargs_lens – lens model keyword arguments

  • diff – float or None, finite average differential scale

Returns:

radial stretch, tangential stretch

radial_tangential_differentials(x, y, kwargs_lens, center_x=0, center_y=0, smoothing_3rd=0.001, smoothing_2nd=None)[source]

Computes the differentials in stretches and directions.

Parameters:
  • x – x-position

  • y – y-position

  • kwargs_lens – lens model keyword arguments

  • center_x – x-coord of center towards which the rotation direction is defined

  • center_y – x-coord of center towards which the rotation direction is defined

  • smoothing_3rd – finite differential length of third order in units of angle

  • smoothing_2nd – float or None, finite average differential scale of Hessian

Returns:

curved_arc_estimate(x, y, kwargs_lens, smoothing=None, smoothing_3rd=0.001, tan_diff=False)[source]

Performs the estimation of the curved arc description at a particular position of an arbitrary lens profile.

Parameters:
  • x – float, x-position where the estimate is provided

  • y – float, y-position where the estimate is provided

  • kwargs_lens – lens model keyword arguments

  • smoothing – (optional) finite differential of second derivative (radial and tangential stretches)

  • smoothing_3rd – differential scale for third derivative to estimate the tangential curvature

  • tan_diff – boolean, if True, also returns the relative tangential stretch differential in tangential direction

Returns:

keyword argument list corresponding to a CURVED_ARC profile at (x, y) given the initial lens model

tangential_average(x, y, kwargs_lens, dr, smoothing=None, num_average=9)[source]

Computes average tangential stretch around position (x, y) within dr in radial direction.

Parameters:
  • x – x-position (float)

  • y – y-position (float)

  • kwargs_lens – lens model keyword argument list

  • dr – averaging scale in radial direction

  • smoothing – smoothing scale of derivative

  • num_average – integer, number of points averaged over within dr in the radial direction

Returns:

curved_arc_finite_area(x, y, kwargs_lens, dr)[source]

Computes an estimated curved arc over a finite extent mimicking the appearance of a finite source with radius dr.

Parameters:
  • x – x-position (float)

  • y – y-position (float)

  • kwargs_lens – lens model keyword argument list

  • dr – radius of finite source

Returns:

keyword arguments of curved arc

lenstronomy.LensModel.lens_param module

class LensParam(lens_model_list, kwargs_fixed, kwargs_lower=None, kwargs_upper=None, kwargs_logsampling=None, num_images=0, solver_type='NONE', num_shapelet_lens=0)[source]

Bases: object

Class to handle the lens model parameter.

__init__(lens_model_list, kwargs_fixed, kwargs_lower=None, kwargs_upper=None, kwargs_logsampling=None, num_images=0, solver_type='NONE', num_shapelet_lens=0)[source]
Parameters:
  • lens_model_list – list of strings of lens model names

  • kwargs_fixed – list of keyword arguments for model parameters to be held fixed

  • kwargs_lower – list of keyword arguments of the lower bounds of the model parameters

  • kwargs_upper – list of keyword arguments of the upper bounds of the model parameters

  • kwargs_logsampling – list of keyword arguments of parameters to be sampled in log10 space

  • num_images – number of images to be constrained by a non-linear solver (only relevant when shapelet potential functions are used)

  • solver_type – string, type of non-linear solver (only relevant in this class when ‘SHAPELETS’ is the solver type)

  • num_shapelet_lens – integer, number of shapelets in the lensing potential (only relevant when ‘SHAPELET’ lens model is used)

get_params(args, i)[source]
Parameters:
  • args – tuple of individual floats of sampling argument

  • i – integer, index at the beginning of the tuple for read out to keyword argument convention

Returns:

kwargs_list, index at the end of read out of this model component

set_params(kwargs_list)[source]
Parameters:

kwargs_list – keyword argument list of lens model components

Returns:

tuple of arguments (floats) that are being sampled

num_param()[source]
Returns:

integer, number of free parameters being sampled from the lens model components

lenstronomy.LensModel.profile_integrals module

class ProfileIntegrals(profile_class)[source]

Bases: object

class to perform integrals of spherical profiles to compute: - projected densities - enclosed densities - projected enclosed densities

__init__(profile_class)[source]
Parameters:

profile_class – list of lens models

mass_enclosed_3d(r, kwargs_profile, lens_param=False)[source]

Computes the mass enclosed within a sphere of radius r.

Parameters:
  • r – radius (arcsec)

  • kwargs_profile – keyword argument list with lens model parameters

  • lens_param – boolean, if True uses the lens model parameterization in computing the 3d density convention and the return is the convergence

Returns:

3d mass enclosed of r

density_2d(r, kwargs_profile, lens_param=False)[source]

Computes the projected density along the line-of-sight.

Parameters:
  • r – radius (arcsec)

  • kwargs_profile – keyword argument list with lens model parameters

  • lens_param – boolean, if True uses the lens model parameterization in computing the 3d density convention and the return is the convergence

Returns:

2d projected density at projected radius r

mass_enclosed_2d(r, kwargs_profile)[source]

computes the mass enclosed the projected line-of-sight :param r: radius (arcsec) :param kwargs_profile: keyword argument list with lens model parameters :return: projected mass enclosed radius r

lenstronomy.LensModel.profile_list_base module

class ProfileListBase(lens_model_list, numerical_alpha_class=None, lens_redshift_list=None, z_source_convention=None, kwargs_interp=None, kwargs_synthesis=None)[source]

Bases: object

Class that manages the list of lens model class instances.

This class is applicable for single plane and multi plane lensing

__init__(lens_model_list, numerical_alpha_class=None, lens_redshift_list=None, z_source_convention=None, kwargs_interp=None, kwargs_synthesis=None)[source]
Parameters:
  • lens_model_list – list of strings with lens model names

  • numerical_alpha_class – an instance of a custom class for use in NumericalAlpha() lens model deflection angles as a lens model. See the documentation in Profiles.numerical_deflections

  • kwargs_interp – interpolation keyword arguments specifying the numerics. See description in the Interpolate() class. Only applicable for ‘INTERPOL’ and ‘INTERPOL_SCALED’ models.

  • kwargs_synthesis – keyword arguments for the ‘SYNTHESIS’ lens model, if applicable

set_static(kwargs_list)[source]
Parameters:

kwargs_list – list of keyword arguments for each profile

Returns:

kwargs_list

set_dynamic()[source]

Frees cache set by static model (if exists) and re-computes all lensing quantities each time a definition is called assuming different parameters are executed. This is the default mode if not specified as set_static()

Returns:

None

lenstronomy.LensModel.single_plane module

class SinglePlane(lens_model_list, numerical_alpha_class=None, lens_redshift_list=None, z_source_convention=None, kwargs_interp=None, kwargs_synthesis=None)[source]

Bases: ProfileListBase

Class to handle an arbitrary list of lens models in a single lensing plane.

ray_shooting(x, y, kwargs, k=None)[source]

Maps image to source position (inverse deflection).

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

Returns:

source plane positions corresponding to (x, y) in the image plane

fermat_potential(x_image, y_image, kwargs_lens, x_source=None, y_source=None, k=None)[source]

Fermat potential (negative sign means earlier arrival time)

Parameters:
  • x_image – image position

  • y_image – image position

  • x_source – source position

  • y_source – source position

  • kwargs_lens – list of keyword arguments of lens model parameters matching the lens model classes

  • k

Returns:

fermat potential in arcsec**2 without geometry term (second part of Eqn 1 in Suyu et al. 2013) as a list

potential(x, y, kwargs, k=None)[source]

Lensing potential.

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

Returns:

lensing potential in units of arcsec^2

alpha(x, y, kwargs, k=None)[source]

Deflection angles.

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

Returns:

deflectionangles in units of arcsec

hessian(x, y, kwargs, k=None)[source]

Hessian matrix.

Parameters:
  • x (numpy array) – x-position (preferentially arcsec)

  • y (numpy array) – y-position (preferentially arcsec)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • k – only evaluate the k-th lens model

Returns:

f_xx, f_xy, f_yx, f_yy components

mass_3d(r, kwargs, bool_list=None)[source]

Computes the mass within a 3d sphere of radius r.

if you want to have physical units of kg, you need to multiply by this factor: const.arcsec ** 2 * self._cosmo.dd * self._cosmo.ds / self._cosmo.dds * const.Mpc * const.c ** 2 / (4 * np.pi * const.G) grav_pot = -const.G * mass_dim / (r * const.arcsec * self._cosmo.dd * const.Mpc)

Parameters:
  • r – radius (in angular units)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • bool_list – list of bools that are part of the output

Returns:

mass (in angular units, modulo epsilon_crit)

mass_2d(r, kwargs, bool_list=None)[source]

Computes the mass enclosed a projected (2d) radius r.

The mass definition is such that:

\[\alpha = mass_2d / r / \pi\]

with alpha is the deflection angle

Parameters:
  • r – radius (in angular units)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • bool_list – list of bools that are part of the output

Returns:

projected mass (in angular units, modulo epsilon_crit)

density(r, kwargs, bool_list=None)[source]

3d mass density at radius r The integral in the LOS projection of this quantity results in the convergence quantity.

Parameters:
  • r – radius (in angular units)

  • kwargs – list of keyword arguments of lens model parameters matching the lens model classes

  • bool_list – list of bools that are part of the output

Returns:

mass density at radius r (in angular units, modulo epsilon_crit)

Module contents