lenstronomy.Analysis package

Submodules

lenstronomy.Analysis.kinematics_api module

class lenstronomy.Analysis.kinematics_api.KinematicsAPI(z_lens, z_source, kwargs_model, kwargs_aperture, kwargs_seeing, anisotropy_model, cosmo=None, lens_model_kinematics_bool=None, light_model_kinematics_bool=None, multi_observations=False, kwargs_numerics_galkin=None, analytic_kinematics=False, Hernquist_approx=False, MGE_light=False, MGE_mass=False, kwargs_mge_light=None, kwargs_mge_mass=None, sampling_number=1000, num_kin_sampling=1000, num_psf_sampling=100)[source]

Bases: object

this class contains routines to compute time delays, magnification ratios, line of sight velocity dispersions etc for a given lens model

galkin_settings(kwargs_lens, kwargs_lens_light, r_eff=None, theta_E=None, gamma=None)[source]
Parameters:
  • kwargs_lens – lens model keyword argument list
  • kwargs_lens_light – deflector light keyword argument list
  • r_eff – half-light radius (optional)
  • theta_E – Einstein radius (optional)
  • gamma – local power-law slope at the Einstein radius (optional)
Returns:

Galkin() instance and mass and light profiles configured for the Galkin module

kinematic_lens_profiles(kwargs_lens, MGE_fit=False, model_kinematics_bool=None, theta_E=None, gamma=None, kwargs_mge=None, analytic_kinematics=False)[source]

translates the lenstronomy lens and mass profiles into a (sub) set of profiles that are compatible with the GalKin module to compute the kinematics thereof. The requirement is that the profiles are centered at (0, 0) and that for all profile types there exists a 3d de-projected analytical representation.

Parameters:
  • kwargs_lens – lens model parameters
  • MGE_fit – bool, if true performs the MGE for the mass distribution
  • model_kinematics_bool – bool list of length of the lens model. Only takes a subset of all the models as part of the kinematics computation (can be used to ignore substructure, shear etc that do not describe the main deflector potential
  • theta_E – (optional float) estimate of the Einstein radius. If present, does not numerically compute this quantity in this routine numerically
  • gamma – local power-law slope at the Einstein radius (optional)
  • kwargs_mge – keyword arguments that go into the MGE decomposition routine
  • analytic_kinematics – bool, if True, solves the Jeans equation analytically for the power-law mass profile with Hernquist light profile
Returns:

mass_profile_list, keyword argument list

kinematic_light_profile(kwargs_lens_light, r_eff=None, MGE_fit=False, model_kinematics_bool=None, Hernquist_approx=False, kwargs_mge=None, analytic_kinematics=False)[source]

setting up of the light profile to compute the kinematics in the GalKin module. The requirement is that the profiles are centered at (0, 0) and that for all profile types there exists a 3d de-projected analytical representation.

Parameters:
  • kwargs_lens_light – deflector light model keyword argument list
  • r_eff – (optional float, else=None) Pre-calculated projected half-light radius of the deflector profile. If not provided, numerical calculation is done in this routine if required.
  • MGE_fit – boolean, if True performs a Multi-Gaussian expansion of the radial light profile and returns this solution.
  • model_kinematics_bool – list of booleans to indicate a subset of light profiles to be part of the physical deflector light.
  • Hernquist_approx – boolean, if True replaces the actual light profile(s) with a Hernquist model with matched half-light radius.
  • kwargs_mge – keyword arguments that go into the MGE decomposition routine
  • analytic_kinematics – bool, if True, solves the Jeans equation analytically for the power-law mass profile with Hernquist light profile and adjust the settings accordingly
Returns:

deflector type list, keyword arguments list

kinematics_modeling_settings(anisotropy_model, kwargs_numerics_galkin, analytic_kinematics=False, Hernquist_approx=False, MGE_light=False, MGE_mass=False, kwargs_mge_light=None, kwargs_mge_mass=None, sampling_number=1000, num_kin_sampling=1000, num_psf_sampling=100)[source]
Parameters:
  • anisotropy_model – type of stellar anisotropy model. See details in MamonLokasAnisotropy() class of lenstronomy.GalKin.anisotropy
  • analytic_kinematics – boolean, if True, used the analytic JAM modeling for a power-law profile on top of a Hernquist light profile ATTENTION: This may not be accurate for your specific problem!
  • Hernquist_approx – bool, if True, uses a Hernquist light profile matched to the half light radius of the deflector light profile to compute the kinematics
  • MGE_light – bool, if true performs the MGE for the light distribution
  • MGE_mass – bool, if true performs the MGE for the mass distribution
  • kwargs_numerics_galkin – numerical settings for the integrated line-of-sight velocity dispersion
  • kwargs_mge_mass – keyword arguments that go into the MGE decomposition routine
  • kwargs_mge_light – keyword arguments that go into the MGE decomposition routine
  • sampling_number – number of spectral rendering on a single slit
  • num_kin_sampling – number of kinematic renderings on a total IFU
  • num_psf_sampling – number of PSF displacements for each kinematic rendering on the IFU
Returns:

static transform_kappa_ext(sigma_v, kappa_ext=0)[source]
Parameters:
  • sigma_v – velocity dispersion estimate of the lensing deflector without considering external convergence
  • kappa_ext – external convergence to be used in the mass-sheet degeneracy
Returns:

transformed velocity dispersion

velocity_dispersion(kwargs_lens, kwargs_lens_light, kwargs_anisotropy, r_eff=None, theta_E=None, gamma=None, kappa_ext=0)[source]

API for both, analytic and numerical JAM to compute the velocity dispersion [km/s] This routine uses the galkin_setting() routine for the Galkin configurations (see there what options and input is relevant.

Parameters:
  • kwargs_lens – lens model keyword arguments
  • kwargs_lens_light – lens light model keyword arguments
  • kwargs_anisotropy – stellar anisotropy keyword arguments
  • r_eff – projected half-light radius of the stellar light associated with the deflector galaxy, optional, if set to None will be computed in this function with default settings that may not be accurate.
  • theta_E – Einstein radius (optional)
  • gamma – power-law slope (optional)
  • kappa_ext – external convergence (optional)
Returns:

velocity dispersion [km/s]

velocity_dispersion_analytical(theta_E, gamma, r_eff, r_ani, kappa_ext=0)[source]

computes the LOS velocity dispersion of the lens within a slit of size R_slit x dR_slit and seeing psf_fwhm. The assumptions are a Hernquist light profile and the spherical power-law lens model at the first position and an Osipkov and Merritt (‘OM’) stellar anisotropy distribution.

Further information can be found in the AnalyticKinematics() class.

Parameters:
  • theta_E – Einstein radius
  • gamma – power-low slope of the mass profile (=2 corresponds to isothermal)
  • r_ani – anisotropy radius in units of angles
  • r_eff – projected half-light radius
  • kappa_ext – external convergence not accounted in the lens models
Returns:

velocity dispersion in units [km/s]

velocity_dispersion_map(kwargs_lens, kwargs_lens_light, kwargs_anisotropy, r_eff=None, theta_E=None, gamma=None, kappa_ext=0)[source]

API for both, analytic and numerical JAM to compute the velocity dispersion map with IFU data [km/s]

Parameters:
  • kwargs_lens – lens model keyword arguments
  • kwargs_lens_light – lens light model keyword arguments
  • kwargs_anisotropy – stellar anisotropy keyword arguments
  • r_eff – projected half-light radius of the stellar light associated with the deflector galaxy, optional, if set to None will be computed in this function with default settings that may not be accurate.
  • theta_E – circularized Einstein radius, optional, if not provided will either be computed in this function with default settings or not required
  • gamma – power-law slope at the Einstein radius, optional
  • kappa_ext – external convergence
Returns:

velocity dispersion [km/s]

lenstronomy.Analysis.lens_profile module

class lenstronomy.Analysis.lens_profile.LensProfileAnalysis(lens_model)[source]

Bases: object

class with analysis routines to compute derived properties of the lens model

convergence_peak(kwargs_lens, model_bool_list=None, grid_num=200, grid_spacing=0.01, center_x_init=0, center_y_init=0)[source]

computes the maximal convergence position on a grid and returns its coordinate

Parameters:
  • kwargs_lens – lens model keyword argument list
  • model_bool_list – bool list (optional) to include certain models or not
Returns:

center_x, center_y

effective_einstein_radius(kwargs_lens, center_x=None, center_y=None, model_bool_list=None, grid_num=200, grid_spacing=0.05, get_precision=False, verbose=True)[source]

computes the radius with mean convergence=1

Parameters:
  • kwargs_lens – list of lens model keyword arguments
  • center_x – position of the center (if not set, is attempting to find it from the parameters kwargs_lens)
  • center_y – position of the center (if not set, is attempting to find it from the parameters kwargs_lens)
  • model_bool_list – list of booleans indicating the addition (=True) of a model component in computing the Einstein radius
  • grid_num – integer, number of grid points to numerically evaluate the convergence and estimate the Einstein radius
  • grid_spacing – spacing in angular units of the grid
  • get_precision – If True, return the precision of estimated Einstein radius
  • verbose – boolean, if True prints warning if indication of insufficient result
Returns:

estimate of the Einstein radius

local_lensing_effect(kwargs_lens, ra_pos=0, dec_pos=0, model_list_bool=None)[source]

computes deflection, shear and convergence at (ra_pos,dec_pos) for those part of the lens model not included in the main deflector.

Parameters:
  • kwargs_lens – lens model keyword argument list
  • ra_pos – RA position where to compute the external effect
  • dec_pos – DEC position where to compute the external effect
  • model_list_bool – boolean list indicating which models effect to be added to the estimate
Returns:

alpha_x, alpha_y, kappa_ext, shear1, shear2

mass_fraction_within_radius(kwargs_lens, center_x, center_y, theta_E, numPix=100)[source]

computes the mean convergence of all the different lens model components within a spherical aperture

Parameters:
  • kwargs_lens – lens model keyword argument list
  • center_x – center of the aperture
  • center_y – center of the aperture
  • theta_E – radius of aperture
Returns:

list of average convergences for all the model components

mst_invariant_differential(kwargs_lens, radius, center_x=None, center_y=None, model_list_bool=None, num_points=10)[source]

Average of the radial stretch differential in radial direction, divided by the radial stretch factor.

\[\xi = \frac{\partial \lambda_{\rm rad}}{\partial r} \frac{1}{\lambda_{\rm rad}}\]

This quantity is invariant under the MST. The specific definition is provided by Birrer 2021. Equivalent (proportional) definitions are provided by e.g. Kochanek 2020, Sonnenfeld 2018.

Parameters:
  • kwargs_lens – lens model keyword argument list
  • radius – radius from the center where to compute the MST invariant differential
  • center_x – center position
  • center_y – center position
  • model_list_bool – indicate which part of the model to consider
  • num_points – number of estimates around the radius
Returns:

xi

multi_gaussian_lens(kwargs_lens, center_x=None, center_y=None, model_bool_list=None, n_comp=20)[source]

multi-gaussian lens model in convergence space

Parameters:
  • kwargs_lens
  • n_comp
Returns:

profile_slope(kwargs_lens, radius, center_x=None, center_y=None, model_list_bool=None, num_points=10)[source]

computes the logarithmic power-law slope of a profile. ATTENTION: this is not an observable!

Parameters:
  • kwargs_lens – lens model keyword argument list
  • radius – radius from the center where to compute the logarithmic slope (angular units
  • model_list_bool – bool list, indicate which part of the model to consider
  • num_points – number of estimates around the Einstein radius
Returns:

logarithmic power-law slope

radial_lens_profile(r_list, kwargs_lens, center_x=None, center_y=None, model_bool_list=None)[source]
Parameters:
  • r_list – list of radii to compute the spherically averaged lens light profile
  • center_x – center of the profile
  • center_y – center of the profile
  • kwargs_lens – lens parameter keyword argument list
  • model_bool_list – bool list or None, indicating which profiles to sum over
Returns:

flux amplitudes at r_list radii spherically averaged

lenstronomy.Analysis.light2mass module

lenstronomy.Analysis.light2mass.light2mass_interpol(lens_light_model_list, kwargs_lens_light, numPix=100, deltaPix=0.05, subgrid_res=5, center_x=0, center_y=0)[source]

takes a lens light model and turns it numerically in a lens model (with all lensmodel quantities computed on a grid). Then provides an interpolated grid for the quantities.

Parameters:
  • kwargs_lens_light – lens light keyword argument list
  • numPix – number of pixels per axis for the return interpolation
  • deltaPix – interpolation/pixel size
  • center_x – center of the grid
  • center_y – center of the grid
  • subgrid_res – subgrid for the numerical integrals
Returns:

keyword arguments for ‘INTERPOL’ lens model

lenstronomy.Analysis.light_profile module

class lenstronomy.Analysis.light_profile.LightProfileAnalysis(light_model)[source]

Bases: object

class with analysis routines to compute derived properties of the lens model

ellipticity(kwargs_light, grid_spacing, grid_num, center_x=None, center_y=None, model_bool_list=None)[source]

make sure that the window covers all the light, otherwise the moments may give a too low answers.

Parameters:
  • kwargs_light – keyword argument list of profiles
  • center_x – center of profile, if None takes it from the first profile in kwargs_light
  • center_y – center of profile, if None takes it from the first profile in kwargs_light
  • model_bool_list – list of booleans to select subsets of the profile
  • grid_spacing – grid spacing over which the moments are computed
  • grid_num – grid size over which the moments are computed
Returns:

eccentricities e1, e2

flux_components(kwargs_light, grid_num=400, grid_spacing=0.01)[source]

computes the total flux in each component of the model

Parameters:
  • kwargs_light
  • grid_num
  • grid_spacing
Returns:

half_light_radius(kwargs_light, grid_spacing, grid_num, center_x=None, center_y=None, model_bool_list=None)[source]

computes numerically the half-light-radius of the deflector light and the total photon flux

Parameters:
  • kwargs_light – keyword argument list of profiles
  • center_x – center of profile, if None takes it from the first profile in kwargs_light
  • center_y – center of profile, if None takes it from the first profile in kwargs_light
  • model_bool_list – list of booleans to select subsets of the profile
  • grid_spacing – grid spacing over which the moments are computed
  • grid_num – grid size over which the moments are computed
Returns:

half-light radius

multi_gaussian_decomposition(kwargs_light, grid_spacing=0.01, grid_num=100, model_bool_list=None, n_comp=20, center_x=None, center_y=None)[source]

multi-gaussian decomposition of the lens light profile (in 1-dimension)

Parameters:
  • kwargs_light – keyword argument list of profiles
  • center_x – center of profile, if None takes it from the first profile in kwargs_light
  • center_y – center of profile, if None takes it from the first profile in kwargs_light
  • model_bool_list – list of booleans to select subsets of the profile
  • grid_spacing – grid spacing over which the moments are computed
  • grid_num – grid size over which the moments are computed
  • n_comp – maximum number of Gaussian’s in the MGE
Returns:

amplitudes, sigmas, center_x, center_y

multi_gaussian_decomposition_ellipse(kwargs_light, model_bool_list=None, center_x=None, center_y=None, grid_num=100, grid_spacing=0.05, n_comp=20)[source]

MGE with ellipticity estimate. Attention: numerical grid settings for ellipticity estimate and radial MGE may not necessarily be the same!

Parameters:
  • kwargs_light – keyword argument list of profiles
  • center_x – center of profile, if None takes it from the first profile in kwargs_light
  • center_y – center of profile, if None takes it from the first profile in kwargs_light
  • model_bool_list – list of booleans to select subsets of the profile
  • grid_spacing – grid spacing over which the moments are computed
  • grid_num – grid size over which the moments are computed
  • n_comp – maximum number of Gaussians in the MGE
Returns:

keyword arguments of the elliptical multi Gaussian profile in lenstronomy conventions

radial_light_profile(r_list, kwargs_light, center_x=None, center_y=None, model_bool_list=None)[source]
Parameters:
  • r_list – list of radii to compute the spherically averaged lens light profile
  • center_x – center of the profile
  • center_y – center of the profile
  • kwargs_light – lens light parameter keyword argument list
  • model_bool_list – bool list or None, indicating which profiles to sum over
Returns:

flux amplitudes at r_list radii spherically averaged

lenstronomy.Analysis.td_cosmography module

class lenstronomy.Analysis.td_cosmography.TDCosmography(z_lens, z_source, kwargs_model, cosmo_fiducial=None, lens_model_kinematics_bool=None, light_model_kinematics_bool=None, kwargs_seeing={}, kwargs_aperture={}, anisotropy_model=None, multi_observations=False)[source]

Bases: lenstronomy.Analysis.kinematics_api.KinematicsAPI

class equipped to perform a cosmographic analysis from a lens model with added measurements of time delays and kinematics. This class does not require any cosmological knowledge and can return angular diameter distance estimates self-consistently integrating the kinematics routines and time delay estimates in the lens modeling. This description follows Birrer et al. 2016, 2019.

ddt_dd_from_time_delay_and_kinematics(d_fermat_model, dt_measured, sigma_v_measured, J, kappa_s=0, kappa_ds=0, kappa_d=0)[source]
Parameters:
  • d_fermat_model – relative Fermat potential in units arcsec^2
  • dt_measured – measured relative time delay [days]
  • sigma_v_measured – 1-sigma Gaussian uncertainty in the measured velocity dispersion
  • J – modeled dimensionless kinematic estimate
  • kappa_s – LOS convergence from observer to source
  • kappa_ds – LOS convergence from deflector to source
  • kappa_d – LOS convergence from observer to deflector
Returns:

D_dt, D_d

static ddt_from_time_delay(d_fermat_model, dt_measured, kappa_s=0, kappa_ds=0, kappa_d=0)[source]

Time-delay distance in units of Mpc from the modeled Fermat potential and measured time delay from an image pair.

Parameters:
  • d_fermat_model – relative Fermat potential between two images from the same source in units arcsec^2
  • dt_measured – measured time delay between the same image pair in units of days
Returns:

D_dt, time-delay distance

static ds_dds_from_kinematics(sigma_v, J, kappa_s=0, kappa_ds=0)[source]

computes the estimate of the ratio of angular diameter distances Ds/Dds from the kinematic estimate of the lens and the measured dispersion.

Parameters:
  • sigma_v – velocity dispersion [km/s]
  • J – dimensionless kinematic constraint (see Birrer et al. 2016, 2019)
Returns:

Ds/Dds

fermat_potential(kwargs_lens, kwargs_ps, original_ps_position=False)[source]
Parameters:
  • kwargs_lens – lens model keyword argument list
  • kwargs_ps – point source keyword argument list
  • original_ps_position – boolean (only applies when first point source model is of type ‘LENSED_POSITION’), uses the image positions in the model parameters and does not re-compute images (which might be differently ordered) in case of the lens equation solver
Returns:

Fermat potential of all the image positions in the first point source list entry

time_delays(kwargs_lens, kwargs_ps, kappa_ext=0, original_ps_position=False)[source]

predicts the time delays of the image positions given the fiducial cosmology

Parameters:
  • kwargs_lens – lens model parameters
  • kwargs_ps – point source parameters
  • kappa_ext – external convergence (optional)
  • original_ps_position – boolean (only applies when first point source model is of type ‘LENSED_POSITION’), uses the image positions in the model parameters and does not re-compute images (which might be differently ordered) in case of the lens equation solver
Returns:

time delays at image positions for the fixed cosmology

velocity_dispersion_dimension_less(kwargs_lens, kwargs_lens_light, kwargs_anisotropy, r_eff=None, theta_E=None, gamma=None)[source]

sigma**2 = Dd/Dds * c**2 * J(kwargs_lens, kwargs_light, anisotropy) (Equation 4.11 in Birrer et al. 2016 or Equation 6 in Birrer et al. 2019) J() is a dimensionless and cosmological independent quantity only depending on angular units. This function returns J given the lens and light parameters and the anisotropy choice without an external mass sheet correction.

Parameters:
  • kwargs_lens – lens model keyword arguments
  • kwargs_lens_light – lens light model keyword arguments
  • kwargs_anisotropy – stellar anisotropy keyword arguments
  • r_eff – projected half-light radius of the stellar light associated with the deflector galaxy, optional, if set to None will be computed in this function with default settings that may not be accurate.
  • theta_E – pre-computed Einstein radius (optional)
  • gamma – pre-computed power-law slope of mass profile
Returns:

dimensionless velocity dispersion (see e.g. Birrer et al. 2016, 2019)

velocity_dispersion_map_dimension_less(kwargs_lens, kwargs_lens_light, kwargs_anisotropy, r_eff=None, theta_E=None, gamma=None)[source]

sigma**2 = Dd/Dds * c**2 * J(kwargs_lens, kwargs_light, anisotropy) (Equation 4.11 in Birrer et al. 2016 or Equation 6 in Birrer et al. 2019) J() is a dimensionless and cosmological independent quantity only depending on angular units. This function returns J given the lens and light parameters and the anisotropy choice without an external mass sheet correction. This routine computes the IFU map of the kinematic quantities.

Parameters:
  • kwargs_lens – lens model keyword arguments
  • kwargs_lens_light – lens light model keyword arguments
  • kwargs_anisotropy – stellar anisotropy keyword arguments
  • r_eff – projected half-light radius of the stellar light associated with the deflector galaxy, optional, if set to None will be computed in this function with default settings that may not be accurate.
Returns:

dimensionless velocity dispersion (see e.g. Birrer et al. 2016, 2019)

Module contents