lenstronomy.GalKin package

Submodules

lenstronomy.GalKin.analytic_kinematics module

class AnalyticKinematics(kwargs_cosmo, interpol_grid_num=200, log_integration=True, max_integrate=100, min_integrate=0.0001)[source]

Bases: Anisotropy

Class to compute eqn 20 in Suyu+2010 with a Monte-Carlo from rendering from the light profile distribution and displacing them with a Gaussian seeing convolution.

This class assumes spherical symmetry in light and mass distribution and
  • a Hernquist light profile (parameterised by the half-light radius)

  • a power-law mass profile (parameterized by the Einstein radius and logarithmic slop)

The analytic equations for the kinematics in this approximation are presented e.g. in Suyu et al. 2010 and the spectral rendering approach to compute the seeing convolved slit measurement is presented in Birrer et al. 2016. The stellar anisotropy is parameterised based on Osipkov 1979; Merritt 1985.

WARNING!!! Only supports Osipkov-Merritt anisotropy for now!

All units are meant to be in angular arc seconds. The physical units are fold in through the angular diameter distances

__init__(kwargs_cosmo, interpol_grid_num=200, log_integration=True, max_integrate=100, min_integrate=0.0001)[source]
Parameters:
  • kwargs_cosmo – keyword argument with angular diameter distances entering the Galkin.cosmo class

  • interpol_grid_num – number of interpolations in radius to compute radial velocity dispersion

  • log_integration – perform numerical integration in logarithmic space, setting False may lead to less accurate results

  • max_integrate – maximum radius of integration (in projected arc seconds)

  • min_integrate – minimum drawing/calculation of velocity dispersion (in projected arc seconds)

property max_integrate

Get the maximum range of integration.

property min_integrate

Get the minimum range of integration.

classmethod draw_light(kwargs_light)[source]

Draws a random light tracer particle from the Hernquist light profile.

Parameters:

kwargs_light – keyword argument (list) of the light model

Returns:

3d radius (if possible), 2d projected radius, x-projected coordinate, y-projected coordinate

sigma_s2(r, R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Returns unweighted los velocity dispersion for a specified projected radius, with weight 1.

Parameters:
  • r – 3d radius (not needed for this calculation)

  • R – 2d projected radius (in angular units of arcsec)

  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

Returns:

line-of-sight projected velocity dispersion at projected radius R from 3d radius r

sigma_r2(r, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Equation (19) in Suyu+ 2010.

Parameters:
  • r – 3d radius

  • kwargs_mass – mass profile keyword arguments

  • kwargs_light – light profile keyword arguments

  • kwargs_anisotropy – anisotropy keyword arguments

Returns:

velocity dispersion in [m/s]

I_R_sigma2_and_IR(R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Return I(R)*sigma^2 equation 20 in Suyu 2010 as interpolation in log space, and I(R)

Parameters:
  • R – projected radius

  • kwargs_mass – mass profile keyword arguments

  • kwargs_light – light model keyword arguments

  • kwargs_anisotropy – stellar anisotropy keyword arguments

Returns:

grav_potential(r, kwargs_mass)[source]

Gravitational potential in SI units.

Parameters:
  • r – radius (arc seconds)

  • kwargs_mass

Returns:

gravitational potential

delete_cache()[source]

Deletes temporary cache tight to a specific model.

Returns:

lenstronomy.GalKin.anisotropy module

class Anisotropy(anisotropy_type)[source]

Bases: object

class that handles the kinematic anisotropy sources: Mamon & Lokas 2005 https://arxiv.org/pdf/astro-ph/0405491.pdf

Agnello et al. 2014 https://arxiv.org/pdf/1401.4462.pdf

__init__(anisotropy_type)[source]
Parameters:

anisotropy_type – string, anisotropy model type

beta_r(r, **kwargs)[source]

Returns the anisotropy parameter at a given radius.

Parameters:
  • r – 3d radius

  • kwargs – parameters of the specified anisotropy model

Returns:

beta(r)

K(r, R, **kwargs)[source]

Equation A16 im Mamon & Lokas for Osipkov&Merrit anisotropy.

Parameters:
  • r – 3d radius

  • R – projected 2d radius

  • kwargs – parameters of the specified anisotropy model

Returns:

K(r, R)

anisotropy_solution(r, **kwargs)[source]

The solution to d ln(f)/ d ln(r) = 2 beta(r)

Parameters:
  • r – 3d radius

  • kwargs – parameters of the specified anisotropy model

Returns:

f(r)

delete_anisotropy_cache()[source]

Deletes cached interpolations for a fixed anisotropy model.

Returns:

None

class Const[source]

Bases: object

Constant anisotropy model class See Mamon & Lokas 2005 for details.

__init__()[source]
static K(r, R, beta)[source]

Equation A16 im Mamon & Lokas for constant anisotropy.

Parameters:
  • r – 3d radius

  • R – projected 2d radius

  • beta – anisotropy, float >-0.5

Returns:

K(r, R, beta)

static beta_r(r, beta)[source]

Anisotropy as a function of radius.

Parameters:
  • r – 3d radius

  • beta – anisotropy

Returns:

beta

anisotropy_solution(r, **kwargs)[source]

The solution to d ln(f)/ d ln(r) = 2 beta(r)

Parameters:
  • r – 3d radius

  • kwargs – parameters of the specified anisotropy model

Returns:

f(r)

class Isotropic[source]

Bases: object

Class for isotropic (beta=0) stellar orbits See Mamon & Lokas 2005 for details.

__init__()[source]
static K(r, R)[source]

Equation A16 im Mamon & Lokas for constant anisotropy.

Parameters:
  • r – 3d radius

  • R – projected 2d radius

Returns:

K(r, R)

static beta_r(r)[source]

Anisotropy as a function of radius.

Parameters:

r – 3d radius

Returns:

beta

static anisotropy_solution(r, **kwargs)[source]

The solution to d ln(f)/ d ln(r) = 2 beta(r) See e.g. A3 in Mamon & Lokas.

Parameters:
  • r – 3d radius

  • kwargs – parameters of the specified anisotropy model

Returns:

f(r)

class Radial[source]

Bases: object

Class for radial (beta=1) stellar orbits See Mamon & Lokas 2005 for details.

__init__()[source]
static K(r, R)[source]

Equation A16 im Mamon & Lokas for constant anisotropy.

Parameters:
  • r – 3d radius

  • R – projected 2d radius

Returns:

K(r, R)

static beta_r(r)[source]

Anisotropy as a function of radius.

Parameters:

r – 3d radius

Returns:

beta

static anisotropy_solution(r)[source]

The solution to d ln(f)/ d ln(r) = 2 beta(r) See e.g. A4 in Mamon & Lokas.

Parameters:

r – 3d radius

Returns:

f(r)

class OsipkovMerritt[source]

Bases: object

Class for Osipkov&Merrit stellar orbits See Mamon & Lokas 2005 for details.

__init__()[source]
static K(r, R, r_ani)[source]

Equation A16 im Mamon & Lokas 2005 for Osipkov&Merrit anisotropy.

Parameters:
  • r – 3d radius

  • R – projected 2d radius

  • r_ani – anisotropy radius

Returns:

K(r, R)

static beta_r(r, r_ani)[source]

Anisotropy as a function of radius.

Parameters:
  • r – 3d radius

  • r_ani – anisotropy radius

Returns:

beta

static anisotropy_solution(r, r_ani)[source]

The solution to d ln(f)/ d ln(r) = 2 beta(r) See e.g. A5 in Mamon & Lokas.

Parameters:
  • r – 3d radius

  • r_ani – anisotropy radius

Returns:

f(r)

class GeneralizedOM[source]

Bases: object

Generalized Osipkov&Merrit profile.

see Agnello et al. 2014 https://arxiv.org/pdf/1401.4462.pdf b(r) = beta_inf * r^2 / (r^2 + r_ani^2)

__init__()[source]
static beta_r(r, r_ani, beta_inf)[source]

Anisotropy as a function of radius.

Parameters:
  • r – 3d radius

  • r_ani – anisotropy radius

  • beta_inf – anisotropy at infinity

Returns:

beta

K(r, R, r_ani, beta_inf)[source]

equation19 in Agnello et al. 2014 for k_beta(R, r) such that K(R, r) = (sqrt(r^2 - R^2) + k_beta(R, r)) / r

Parameters:
  • r – 3d radius

  • R – projected 2d radius

  • r_ani – anisotropy radius

  • beta_inf – anisotropy at infinity

Returns:

K(r, R)

static anisotropy_solution(r, r_ani, beta_inf)[source]

The solution to d ln(f)/ d ln(r) = 2 beta(r) See e.g. A5 in Mamon & Lokas with a scaling (nominator of Agnello et al. 2014 Equation (12)

Parameters:
  • r – 3d radius

  • r_ani – anisotropy radius

  • beta_inf – anisotropy at infinity

Returns:

f(r)

delete_cache()[source]

Deletes the interpolation function of the hypergeometic function for a specific beta_inf.

Returns:

deleted self variables

class Colin[source]

Bases: object

Class for stellar orbits anisotropy parameter based on Colin et al.

  1. See Mamon & Lokas 2005 for details

__init__()[source]
static K(r, R, r_ani)[source]

Equation A16 im Mamon & Lokas for Osipkov&Merrit anisotropy.

Parameters:
  • r – 3d radius

  • R – projected 2d radius

  • r_ani – anisotropy radius

Returns:

K(r, R)

static beta_r(r, r_ani)[source]

Anisotropy as a function of radius.

Parameters:
  • r – 3d radius

  • r_ani – anisotropy radius

Returns:

beta

lenstronomy.GalKin.aperture module

class Aperture(aperture_type, **kwargs_aperture)[source]

Bases: object

Defines mask(s) of spectra, can handle IFU and single slit/box type data.

__init__(aperture_type, **kwargs_aperture)[source]
Parameters:
  • aperture_type – string

  • kwargs_aperture – keyword arguments reflecting the aperture type chosen. We refer to the specific class instances for documentation.

aperture_select(ra, dec)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

Returns:

bool, True if photon/ray is within the slit, False otherwise, int of the segment of the IFU

property num_segments

lenstronomy.GalKin.aperture_types module

class Slit(length, width, center_ra=0, center_dec=0, angle=0)[source]

Bases: object

Slit aperture description.

__init__(length, width, center_ra=0, center_dec=0, angle=0)[source]
Parameters:
  • length – length of slit

  • width – width of slit

  • center_ra – center of slit

  • center_dec – center of slit

  • angle – orientation angle of slit, angle=0 corresponds length in RA direction

aperture_select(ra, dec)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

Returns:

bool, True if photon/ray is within the slit, False otherwise

property num_segments

Number of segments with separate measurements of the velocity dispersion.

Returns:

int

slit_select(ra, dec, length, width, center_ra=0, center_dec=0, angle=0)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

  • length – length of slit

  • width – width of slit

  • center_ra – center of slit

  • center_dec – center of slit

  • angle – orientation angle of slit, angle=0 corresponds length in RA direction

Returns:

bool, True if photon/ray is within the slit, False otherwise

class Frame(width_outer, width_inner, center_ra=0, center_dec=0, angle=0)[source]

Bases: object

Rectangular box with a hole in the middle (also rectangular), effectively a frame.

__init__(width_outer, width_inner, center_ra=0, center_dec=0, angle=0)[source]
Parameters:
  • width_outer – width of box to the outer parts

  • width_inner – width of inner removed box

  • center_ra – center of slit

  • center_dec – center of slit

  • angle – orientation angle of slit, angle=0 corresponds length in RA direction

aperture_select(ra, dec)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

Returns:

bool, True if photon/ray is within the slit, False otherwise

property num_segments

Number of segments with separate measurements of the velocity dispersion.

Returns:

int

frame_select(ra, dec, width_outer, width_inner, center_ra=0, center_dec=0, angle=0)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

  • width_outer – width of box to the outer parts

  • width_inner – width of inner removed box

  • center_ra – center of slit

  • center_dec – center of slit

  • angle – orientation angle of slit, angle=0 corresponds length in RA direction

Returns:

bool, True if photon/ray is within the box with a hole, False otherwise

class Shell(r_in, r_out, center_ra=0, center_dec=0)[source]

Bases: object

Shell aperture.

__init__(r_in, r_out, center_ra=0, center_dec=0)[source]
Parameters:
  • r_in – innermost radius to be selected

  • r_out – outermost radius to be selected

  • center_ra – center of the sphere

  • center_dec – center of the sphere

aperture_select(ra, dec)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

Returns:

bool, True if photon/ray is within the slit, False otherwise

property num_segments

Number of segments with separate measurements of the velocity dispersion.

Returns:

int

shell_select(ra, dec, r_in, r_out, center_ra=0, center_dec=0)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

  • r_in – innermost radius to be selected

  • r_out – outermost radius to be selected

  • center_ra – center of the sphere

  • center_dec – center of the sphere

Returns:

boolean, True if within the radial range, False otherwise

class IFUShells(r_bins, center_ra=0, center_dec=0)[source]

Bases: object

Class for an Integral Field Unit spectrograph with azimuthal shells where the kinematics are measured.

__init__(r_bins, center_ra=0, center_dec=0)[source]
Parameters:
  • r_bins – array of radial bins to average the dispersion spectra in ascending order. It starts with the innermost edge to the outermost edge.

  • center_ra – center of the sphere

  • center_dec – center of the sphere

aperture_select(ra, dec)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

Returns:

bool, True if photon/ray is within the slit, False otherwise, index of shell

property num_segments

Number of segments with separate measurements of the velocity dispersion :return: int.

class IFUGrid(x_grid, y_grid)[source]

Bases: object

Class for an Integral Field Unit spectrograph with rectangular grid where the kinematics are measured.

__init__(x_grid, y_grid)[source]
Parameters:
  • x_grid – x coordinates of the grid

  • y_grid – y coordinates of the grid

aperture_select(ra, dec)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

Returns:

bool, True if photon/ray is within the slit, False otherwise, index of shell

property num_segments

Number of segments with separate measurements of the velocity dispersion.

Returns:

int

property x_grid

X coordinates of the grid.

property y_grid

Y coordinates of the grid.

grid_ifu_select(ra, dec, x_grid, y_grid)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

  • x_grid – array of x_grid bins

  • y_grid – array of y_grid bins

Returns:

boolean, True if within the grid range, False otherwise

shell_ifu_select(ra, dec, r_bin, center_ra=0, center_dec=0)[source]
Parameters:
  • ra – angular coordinate of photon/ray

  • dec – angular coordinate of photon/ray

  • r_bin – array of radial bins to average the dispersion spectra in ascending order. It starts with the inner-most edge to the outermost edge.

  • center_ra – center of the sphere

  • center_dec – center of the sphere

Returns:

boolean, True if within the radial range, False otherwise

lenstronomy.GalKin.cosmo module

class Cosmo(d_d, d_s, d_ds)[source]

Bases: object

Cosmological quantities.

__init__(d_d, d_s, d_ds)[source]
Parameters:
  • d_d – angular diameter distance to the deflector

  • d_s – angular diameter distance to the source

  • d_ds – angular diameter distance between deflector and source

arcsec2phys_lens(theta)[source]

Converts are seconds to physical units on the deflector.

Parameters:

theta – angle observed on the sky in units of arc seconds

Returns:

physical distance of the angle in units of Mpc

property epsilon_crit

Returns the critical projected mass density in units of M_sun/Mpc^2 (physical units)

lenstronomy.GalKin.galkin module

class Galkin(kwargs_model, kwargs_aperture, kwargs_psf, kwargs_cosmo, kwargs_numerics=None, analytic_kinematics=False)[source]

Bases: GalkinModel, GalkinObservation

Major class to compute velocity dispersion measurements given light and mass models.

The class supports any mass and light distribution (and superposition thereof) that has a 3d correspondance in their 2d lens model distribution. For models that do not have this correspondance, you may want to apply a Multi-Gaussian Expansion (MGE) on their models and use the MGE to be de-projected to 3d.

The computation follows Mamon&Lokas 2005 and performs the spectral rendering of the seeing convolved apperture with the method introduced by Birrer et al. 2016.

The class supports various types of anisotropy models (see Anisotropy class) and aperture types (see Aperture class).

Solving the Jeans Equation requires a numerical integral over the 3d light and mass profile (see Mamon&Lokas 2005). This class (as well as the dedicated LightModel and MassModel classes) perform those integral numerically with an interpolated grid.

The seeing convolved integral over the aperture is computed by rendering spectra (light weighted LOS kinematics) from the light distribution.

The cosmology assumed to compute the physical mass and distances are set via the kwargs_cosmo keyword arguments.

d_d: Angular diameter distance to the deflector (in Mpc) d_s: Angular diameter distance to the source (in Mpc) d_ds: Angular diameter distance from the deflector to the source (in Mpc)

The numerical options can be chosen through the kwargs_numerics keywords

interpol_grid_num: number of interpolation points in the light and mass profile (radially). This number should be chosen high enough to accurately describe the true light profile underneath. log_integration: bool, if True, performs the interpolation and numerical integration in log-scale.

max_integrate: maximum 3d radius to where the numerical integration of the Jeans Equation solver is made. This value should be large enough to contain most of the light and to lead to a converged result. min_integrate: minimal integration value. This value should be very close to zero but some mass and light profiles are diverging and a numerically stabel value should be chosen.

These numerical options should be chosen to allow for a converged result (within your tolerance) but not too conservative to impact too much the computational cost. Reasonable values might depend on the specific problem.

__init__(kwargs_model, kwargs_aperture, kwargs_psf, kwargs_cosmo, kwargs_numerics=None, analytic_kinematics=False)[source]
Parameters:
  • kwargs_model – keyword arguments describing the model components

  • kwargs_aperture – keyword arguments describing the spectroscopic aperture, see Aperture() class

  • kwargs_psf – keyword argument specifying the PSF of the observation

  • kwargs_cosmo – keyword arguments that define the cosmology in terms of the angular diameter distances involved

  • kwargs_numerics – numerics keyword arguments

  • analytic_kinematics – bool, if True uses the analytic kinematic model

dispersion(kwargs_mass, kwargs_light, kwargs_anisotropy, sampling_number=1000)[source]

Computes the averaged LOS velocity dispersion in the slit (convolved)

Parameters:
  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

  • sampling_number – int, number of spectral sampling of the light distribution

Returns:

integrated LOS velocity dispersion in units [km/s]

dispersion_map(kwargs_mass, kwargs_light, kwargs_anisotropy, num_kin_sampling=1000, num_psf_sampling=100)[source]

Computes the velocity dispersion in each Integral Field Unit.

Parameters:
  • kwargs_mass – keyword arguments of the mass model

  • kwargs_light – keyword argument of the light model

  • kwargs_anisotropy – anisotropy keyword arguments

  • num_kin_sampling – int, number of draws from a kinematic prediction of a LOS

  • num_psf_sampling – int, number of displacements/render from a spectra to be displaced on the IFU

Returns:

ordered array of velocity dispersions [km/s] for each unit

dispersion_map_grid_convolved(kwargs_mass, kwargs_light, kwargs_anisotropy, supersampling_factor=1, voronoi_bins=None)[source]

Computes the velocity dispersion in each Integral Field Unit.

Parameters:
  • kwargs_mass – keyword arguments of the mass model

  • kwargs_light – keyword argument of the light model

  • kwargs_anisotropy – anisotropy keyword arguments

  • supersampling_factor – sampling factor for the grid to do the 2D convolution on

  • voronoi_bins – mapping of the voronoi bins, bin indices should start from 0, -1 values for pixels not binned

Returns:

ordered array of velocity dispersions [km/s] for each unit

lenstronomy.GalKin.galkin_model module

class GalkinModel(kwargs_model, kwargs_cosmo, kwargs_numerics=None, analytic_kinematics=False)[source]

Bases: object

This class handles all the kinematic modeling aspects of Galkin Excluded are observational conditions (seeing, aperture etc) Major class to compute velocity dispersion measurements given light and mass models.

The class supports any mass and light distribution (and superposition thereof) that has a 3d correspondance in their 2d lens model distribution. For models that do not have this correspondence, you may want to apply a Multi-Gaussian Expansion (MGE) on their models and use the MGE to be de-projected to 3d.

The computation follows Mamon&Lokas 2005.

The class supports various types of anisotropy models (see Anisotropy class). Solving the Jeans Equation requires a numerical integral over the 3d light and mass profile (see Mamon&Lokas 2005). This class (as well as the dedicated LightModel and MassModel classes) perform those integral numerically with an interpolated grid.

The cosmology assumed to compute the physical mass and distances are set via the kwargs_cosmo keyword arguments.

d_d: Angular diameter distance to the deflector (in Mpc) d_s: Angular diameter distance to the source (in Mpc) d_ds: Angular diameter distance from the deflector to the source (in Mpc)

The numerical options can be chosen through the kwargs_numerics keywords

interpol_grid_num: number of interpolation points in the light and mass profile (radially). This number should be chosen high enough to accurately describe the true light profile underneath. log_integration: bool, if True, performs the interpolation and numerical integration in log-scale.

max_integrate: maximum 3d radius to where the numerical integration of the Jeans Equation solver is made. This value should be large enough to contain most of the light and to lead to a converged result. min_integrate: minimal integration value. This value should be very close to zero but some mass and light profiles are diverging and a numerically stable value should be chosen.

These numerical options should be chosen to allow for a converged result (within your tolerance) but not too conservative to impact too much the computational cost. Reasonable values might depend on the specific problem.

__init__(kwargs_model, kwargs_cosmo, kwargs_numerics=None, analytic_kinematics=False)[source]
Parameters:
  • kwargs_model – keyword arguments describing the model components

  • kwargs_cosmo – keyword arguments that define the cosmology in terms of the angular diameter distances involved

  • kwargs_numerics – numerics keyword arguments

  • analytic_kinematics – bool, if True uses the analytic kinematic model

check_df(r, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

checks whether the phase space distribution function of a given anisotropy model is positive. Currently this is implemented by the relation provided by Ciotti and Morganti 2010 equation (10) https://arxiv.org/pdf/1006.2344.pdf

Parameters:
  • r – 3d radius to check slope-anisotropy constraint

  • kwargs_mass – keyword arguments for mass (lens) profile

  • kwargs_light – keyword arguments for light profile

  • kwargs_anisotropy – keyword arguments for stellar anisotropy distribution

Returns:

equation (10) >= 0 for physical interpretation

lenstronomy.GalKin.light_profile module

class LightProfile(profile_list, interpol_grid_num=2000, max_interpolate=1000, min_interpolate=0.001, max_draw=None)[source]

Bases: object

Class to deal with the light distribution for GalKin.

In particular, this class allows for:
  • (faster) interpolated calculation for a given profile (for a range that the Jeans equation is computed)

  • drawing 3d and 2d distributions from a given (spherical) profile (within bounds where the Jeans equation is expected to be accurate)

  • 2d projected profiles within the 3d integration range (truncated)

__init__(profile_list, interpol_grid_num=2000, max_interpolate=1000, min_interpolate=0.001, max_draw=None)[source]
Parameters:
  • profile_list – list of light profiles for LightModel module (must support light_3d() functionalities)

  • interpol_grid_num – int; number of interpolation steps (logarithmically between min and max value)

  • max_interpolate – float; maximum interpolation of 3d light profile

  • min_interpolate – float; minimum interpolate (and also drawing of light profile)

  • max_draw – float; (optional) if set, draws up to this radius, else uses max_interpolate value

light_3d(r, kwargs_list)[source]

Three-dimensional light profile.

Parameters:
  • r – 3d radius

  • kwargs_list – list of keyword arguments of light profiles (see LightModule)

Returns:

flux per 3d volume at radius r

light_3d_interp(r, kwargs_list, new_compute=False)[source]

Interpolated three-dimensional light profile within bounds [min_interpolate, max_interpolate] in logarithmic units with interpol_grid_num numbers of interpolation steps.

Parameters:
  • r – 3d radius

  • kwargs_list – list of keyword arguments of light profiles (see LightModule)

  • new_compute – boolean, if True, re-computes the interpolation (becomes valid with updated kwargs_list argument)

Returns:

flux per 3d volume at radius r

light_2d(R, kwargs_list)[source]

Projected light profile (integrated to infinity in the projected axis)

Parameters:
  • R – projected 2d radius

  • kwargs_list – list of keyword arguments of light profiles (see LightModule)

Returns:

projected surface brightness

light_2d_finite(R, kwargs_list)[source]

Projected light profile (integrated to FINITE 3d boundaries from the max_interpolate)

Parameters:
  • R – projected 2d radius (between min_interpolate and max_interpolate

  • kwargs_list – list of keyword arguments of light profiles (see LightModule)

Returns:

projected surface brightness

draw_light_2d_linear(kwargs_list, n=1, new_compute=False)[source]

Constructs the CDF and draws from it random realizations of projected radii R The interpolation of the CDF is done in linear projected radius space.

Parameters:
  • kwargs_list – list of keyword arguments of light profiles (see LightModule)

  • n – int; number of draws

  • new_compute – boolean, if True, re-computes the interpolation (becomes valid with updated kwargs_list argument)

Returns:

draw of projected radius for the given light profile distribution

draw_light_2d(kwargs_list, n=1, new_compute=False)[source]

Constructs the CDF and draws from it random realizations of projected radii R CDF is constructed in logarithmic projected radius spacing.

Parameters:
  • kwargs_list – light model keyword argument list

  • n – int, number of draws per functino call

  • new_compute – re-computes the interpolated CDF

Returns:

realization of projected radius following the distribution of the light model

draw_light_3d(kwargs_list, n=1, new_compute=False)[source]

Constructs the CDF and draws from it random realizations of 3D radii r.

Parameters:
  • kwargs_list – light model keyword argument list

  • n – int, number of draws per function call

  • new_compute – re-computes the interpolated CDF

Returns:

realization of projected radius following the distribution of the light model

delete_cache()[source]

Deletes cached interpolation function of the CDF for a specific light profile.

Returns:

None

lenstronomy.GalKin.numeric_kinematics module

class NumericKinematics(kwargs_model, kwargs_cosmo, interpol_grid_num=200, log_integration=True, max_integrate=100, min_integrate=0.0001, max_light_draw=None, lum_weight_int_method=True)[source]

Bases: Anisotropy

__init__(kwargs_model, kwargs_cosmo, interpol_grid_num=200, log_integration=True, max_integrate=100, min_integrate=0.0001, max_light_draw=None, lum_weight_int_method=True)[source]

What we need: - max projected R to have ACCURATE I_R_sigma values - make sure everything outside cancels out (or is not rendered)

Parameters:
  • interpol_grid_num – number of interpolation bins for integrand and interpolated functions

  • log_integration – bool, if True, performs the numerical integral in log space distance (adviced) (only applies for lum_weight_int_method=True). If set to False, may lead to less accurate results

  • max_integrate – maximum radius (in arc seconds) of the Jeans equation integral (assumes zero tracer particles outside this radius)

  • max_light_draw – float; (optional) if set, draws up to this radius, else uses max_interpolate value

  • lum_weight_int_method – bool, luminosity weighted dispersion integral to calculate LOS projected Jean’s solution. ATTENTION: currently less accurate than 3d solution

  • min_integrate

property lum_weight_int_method

Get the luminosity weighted integration method.

property max_integrate

Get the maximum range of integration.

property min_integrate

Get the minimum range of integration.

lum_weighted_vel_disp(R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Luminosity-weighted line-of-sight velocity dispersion within a radius R.

Parameters:
  • R – 2d projected radius (in angular units of arcsec)

  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

Returns:

average velocity dispersion [km/s]

sigma_s2(r, R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Returns unweighted los velocity dispersion for a specified 3d and projected radius (if lum_weight_int_method=True then the 3d radius is not required and the function directly performs the luminosity weighted integral in projection at R)

Parameters:
  • r – 3d radius (not needed for this calculation)

  • R – 2d projected radius (in angular units of arcsec)

  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

Returns:

weighted line-of-sight projected velocity dispersion at projected radius R with weights I

sigma_s2_project(R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Returns luminosity-weighted los velocity dispersion for a specified projected radius R and weight.

Parameters:
  • R – 2d projected radius (in angular units of arcsec)

  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

Returns:

line-of-sight projected velocity dispersion at projected radius R

sigma_s2_r(r, R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Returns unweighted los velocity dispersion for a specified 3d radius r at projected radius R.

Parameters:
  • r – 3d radius (not needed for this calculation)

  • R – 2d projected radius (in angular units of arcsec)

  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

Returns:

line-of-sight projected velocity dispersion at projected radius R from 3d radius r

sigma_r2(r, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

computes numerically the solution of the Jeans equation for a specific 3d radius E.g. Equation (A1) of Mamon & Lokas https://arxiv.org/pdf/astro-ph/0405491.pdf

\[l(r) \sigma_r(r) ^ 2 = 1/f(r) \int_r^{\infty} f(s) l(s) G M(s) / s^2 ds\]

where l(r) is the 3d light profile M(s) is the enclosed 3d mass f is the solution to d ln(f)/ d ln(r) = 2 beta(r)

Parameters:
  • r – 3d radius

  • kwargs_mass – mass model parameters (following lenstronomy lens model conventions)

  • kwargs_light – deflector light parameters (following lenstronomy light model conventions)

  • kwargs_anisotropy – anisotropy parameters, may vary according to anisotropy type chosen. We refer to the Anisotropy() class for details on the parameters.

Returns:

sigma_r**2

mass_3d(r, kwargs)[source]

Mass enclosed a 3d radius.

Parameters:
  • r – in arc seconds

  • kwargs – lens model parameters in arc seconds

Returns:

mass enclosed physical radius in kg

grav_potential(r, kwargs_mass)[source]

Gravitational potential in SI units.

Parameters:
  • r – radius (arc seconds)

  • kwargs_mass

Returns:

gravitational potential

draw_light(kwargs_light)[source]
Parameters:

kwargs_light – keyword argument (list) of the light model

Returns:

3d radius (if possible), 2d projected radius, x-projected coordinate, y-projected coordinate

delete_cache()[source]

Delete interpolation function for a specific mass and light profile as well as for a specific anisotropy model.

Returns:

I_R_sigma2_and_IR(R, kwargs_mass, kwargs_light, kwargs_anisotropy)[source]

Return I(R)*sigma^2 equation A15 in Mamon&Lokas 2005 as interpolation in log space, and I(R)

Parameters:
  • R – projected radius

  • kwargs_mass – mass profile keyword arguments

  • kwargs_light – light model keyword arguments

  • kwargs_anisotropy – stellar anisotropy keyword arguments

Returns:

a tuple containing (I(R)*sigma^2, IR)

lenstronomy.GalKin.observation module

class GalkinObservation(kwargs_aperture, kwargs_psf)[source]

Bases: PSF, Aperture

This class sets the base for the observational properties (aperture and seeing condition)

__init__(kwargs_aperture, kwargs_psf)[source]
Parameters:
  • psf_type – string, point spread function type, current support for ‘GAUSSIAN’ and ‘MOFFAT’

  • kwargs_psf – keyword argument describing the relevant parameters of the PSF.

lenstronomy.GalKin.psf module

class PSF(psf_type, **kwargs_psf)[source]

Bases: object

General class to handle the PSF in the GalKin module for rendering the displacement of photons/spectro.

__init__(psf_type, **kwargs_psf)[source]
Parameters:
  • psf_type – string, point spread function type, current support for ‘GAUSSIAN’ and ‘MOFFAT’

  • kwargs_psf – keyword argument describing the relevant parameters of the PSF.

displace_psf(x, y)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

Returns:

x’, y’ displaced by the two-dimensional PSF distribution function

convolution_kernel(delta_pix, num_pix=21)[source]

Normalized convolution kernel.

Parameters:
  • delta_pix – pixel scale of kernel

  • num_pix – number of pixels per axis of the kernel

Returns:

2d numpy array of kernel

convolution_kernel_grid(x, y)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

Returns:

psf value at x and y grid positions

class PSFGaussian(fwhm)[source]

Bases: object

Gaussian PSF.

__init__(fwhm)[source]
Parameters:

fwhm – full width at half maximum seeing condition

displace_psf(x, y)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

Returns:

x’, y’ displaced by the two-dimensional PSF distribution function

convolution_kernel(delta_pix, num_pix=21)[source]

Normalized convolution kernel.

Parameters:
  • delta_pix – pixel scale of kernel

  • num_pix – number of pixels per axis of the kernel

Returns:

2d numpy array of kernel

convolution_kernel_grid(x, y)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

Returns:

psf value at x and y grid positions

property fwhm

Retrieve FWHM of PSF if stored as a private variable.

class PSFMoffat(fwhm, moffat_beta)[source]

Bases: object

Moffat PSF.

__init__(fwhm, moffat_beta)[source]
Parameters:
  • fwhm – full width at half maximum seeing condition

  • moffat_beta – float, beta parameter of Moffat profile

displace_psf(x, y)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

Returns:

x’, y’ displaced by the two-dimensional PSF distribution function

convolution_kernel(delta_pix, num_pix=21)[source]

Normalized convolution kernel.

Parameters:
  • delta_pix – pixel scale of kernel

  • num_pix – number of pixels per axis of the kernel

Returns:

2d numpy array of kernel

convolution_kernel_grid(x, y)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

Returns:

psf value at x and y grid positions

lenstronomy.GalKin.velocity_util module

hyp_2F1(a, b, c, z)[source]

http://docs.sympy.org/0.7.1/modules/mpmath/functions/hypergeometric.html

displace_PSF_gaussian(x, y, FWHM)[source]
Parameters:
  • x – x-coord (arc sec)

  • y – y-coord (arc sec)

  • FWHM – psf size (arc sec)

Returns:

x’, y’ random displaced according to psf

moffat_r(r, alpha, beta)[source]

Moffat profile.

Parameters:
  • r – radial coordinate

  • alpha – Moffat parameter

  • beta – exponent

Returns:

Moffat profile

moffat_fwhm_alpha(FWHM, beta)[source]

Computes alpha parameter from FWHM and beta for a Moffat profile.

Parameters:
  • FWHM – full width at half maximum

  • beta – beta parameter of Moffat profile

Returns:

alpha parameter of Moffat profile

draw_moffat_r(FWHM, beta)[source]
Parameters:
  • FWHM – full width at half maximum

  • beta – Moffat beta parameter

Returns:

draw from radial Moffat distribution

displace_PSF_moffat(x, y, FWHM, beta)[source]
Parameters:
  • x – x-coordinate of light ray

  • y – y-coordinate of light ray

  • FWHM – full width at half maximum

  • beta – Moffat beta parameter

Returns:

displaced ray by PSF

draw_cdf_Y(beta)[source]

Draw c.d.f for Moffat function according to Berge et al. Ufig paper, equation B2 cdf(Y) = 1-Y**(1-beta)

Returns:

project2d_random(r)[source]

Draws a random projection from radius r in 2d and 1d.

Parameters:

r – 3d radius

Returns:

R, x, y

draw_xy(R)[source]
Parameters:

R – projected radius

Returns:

draw_hernquist(a)[source]
Parameters:

a – 0.551*r_eff

Returns:

realisation of radius of Hernquist luminosity weighting in 3d

Module contents