lenstronomy.GalKin package¶
Submodules¶
lenstronomy.GalKin.analytic_kinematics module¶
-
class
lenstronomy.GalKin.analytic_kinematics.
AnalyticKinematics
(kwargs_cosmo, interpol_grid_num=100, log_integration=False, max_integrate=100, min_integrate=0.001)[source]¶ Bases:
lenstronomy.GalKin.anisotropy.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.
all units are meant to be in angular arc seconds. The physical units are fold in through the angular diameter distances
-
static
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
-
grav_potential
(r, kwargs_mass)[source]¶ Gravitational potential in SI units
Parameters: - r – radius (arc seconds)
- kwargs_mass –
Returns: gravitational potential
-
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]
-
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
lenstronomy.GalKin.anisotropy module¶
-
class
lenstronomy.GalKin.anisotropy.
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
-
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)
-
-
class
lenstronomy.GalKin.anisotropy.
Const
[source]¶ Bases:
object
constant anisotropy model class See Mamon & Lokas 2005 for details
-
static
K
(r, R, beta)[source]¶ equation A16 im Mamon & Lokas for constant anisotropy
Parameters: - r – 3d radius
- R – projected 2d radius
- beta – anisotropy
Returns: K(r, R, beta)
-
static
-
class
lenstronomy.GalKin.anisotropy.
Isotropic
[source]¶ Bases:
object
class for isotropic (beta=0) stellar orbits See Mamon & Lokas 2005 for details
-
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
-
class
lenstronomy.GalKin.anisotropy.
Radial
[source]¶ Bases:
object
class for radial (beta=1) stellar orbits See Mamon & Lokas 2005 for details
-
K
(r, R)[source]¶ equation A16 im Mamon & Lokas for constant anisotropy
Parameters: - r – 3d radius
- R – projected 2d radius
Returns: K(r, R)
-
-
class
lenstronomy.GalKin.anisotropy.
OsipkovMerritt
[source]¶ Bases:
object
class for Osipkov&Merrit stellar orbits See Mamon & Lokas 2005 for details
-
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)
-
-
class
lenstronomy.GalKin.anisotropy.
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)
-
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)
-
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)
-
-
class
lenstronomy.GalKin.anisotropy.
Colin
[source]¶ Bases:
object
class for stellar orbits anisotropy parameter based on Colin et al. (2000) See Mamon & Lokas 2005 for details
lenstronomy.GalKin.aperture module¶
-
class
lenstronomy.GalKin.aperture.
Aperture
(aperture_type, **kwargs_aperture)[source]¶ Bases:
object
defines mask(s) of spectra, can handle IFU and single slit/box type data.
-
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
-
num_segments
¶
-
lenstronomy.GalKin.aperture_types module¶
-
class
lenstronomy.GalKin.aperture_types.
Slit
(length, width, center_ra=0, center_dec=0, angle=0)[source]¶ Bases:
object
Slit aperture description
-
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
-
num_segments
¶ number of segments with separate measurements of the velocity dispersion :return: int
-
-
lenstronomy.GalKin.aperture_types.
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
lenstronomy.GalKin.aperture_types.
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
-
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
-
num_segments
¶ number of segments with separate measurements of the velocity dispersion :return: int
-
-
lenstronomy.GalKin.aperture_types.
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
lenstronomy.GalKin.aperture_types.
Shell
(r_in, r_out, center_ra=0, center_dec=0)[source]¶ Bases:
object
Shell aperture
-
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
-
num_segments
¶ number of segments with separate measurements of the velocity dispersion :return: int
-
-
lenstronomy.GalKin.aperture_types.
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
lenstronomy.GalKin.aperture_types.
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
-
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
-
num_segments
¶ number of segments with separate measurements of the velocity dispersion :return: int
-
-
lenstronomy.GalKin.aperture_types.
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
lenstronomy.GalKin.cosmo.
Cosmo
(d_d, d_s, d_ds)[source]¶ Bases:
object
cosmological quantities
-
arcsec2phys_lens
(theta)[source]¶ converts are seconds to physical units on the deflector :param theta: angle observed on the sky in units of arc seconds :return: pyhsical distance of the angle in units of Mpc
-
epsilon_crit
¶ returns the critical projected mass density in units of M_sun/Mpc^2 (physical units)
-
lenstronomy.GalKin.galkin module¶
-
class
lenstronomy.GalKin.galkin.
Galkin
(kwargs_model, kwargs_aperture, kwargs_psf, kwargs_cosmo, kwargs_numerics=None, analytic_kinematics=False)[source]¶ Bases:
lenstronomy.GalKin.galkin_model.GalkinModel
,lenstronomy.GalKin.observation.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.
-
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
lenstronomy.GalKin.galkin_model module¶
-
class
lenstronomy.GalKin.galkin_model.
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.
-
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
- theta_E – Einstein radius in arc seconds
- gamma – power-law slope
- a_ani – scaled transition radius of the OM anisotropy distribution
- r_eff – half-light radius in arc seconds
Returns: equation (10) >= 0 for physical interpretation
lenstronomy.GalKin.light_profile module¶
-
class
lenstronomy.GalKin.light_profile.
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)
-
delete_cache
()[source]¶ deletes cached interpolation function of the CDF for a specific light profile
Returns: None
-
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_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_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
-
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
-
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
lenstronomy.GalKin.numeric_kinematics module¶
-
class
lenstronomy.GalKin.numeric_kinematics.
NumericKinematics
(kwargs_model, kwargs_cosmo, interpol_grid_num=1000, log_integration=True, max_integrate=1000, min_integrate=0.0001, max_light_draw=None, lum_weight_int_method=True)[source]¶ Bases:
lenstronomy.GalKin.anisotropy.Anisotropy
-
delete_cache
()[source]¶ delete interpolation function for a specific mass and light profile as well as for a specific anisotropy model
Returns:
-
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
-
grav_potential
(r, kwargs_mass)[source]¶ Gravitational potential in SI units
Parameters: - r – radius (arc seconds)
- kwargs_mass –
Returns: gravitational potential
-
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
-
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
-
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
-
lenstronomy.GalKin.observation module¶
-
class
lenstronomy.GalKin.observation.
GalkinObservation
(kwargs_aperture, kwargs_psf)[source]¶ Bases:
lenstronomy.GalKin.psf.PSF
,lenstronomy.GalKin.aperture.Aperture
this class sets the base for the observational properties (aperture and seeing condition)
lenstronomy.GalKin.psf module¶
-
class
lenstronomy.GalKin.psf.
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
lenstronomy.GalKin.velocity_util module¶
-
lenstronomy.GalKin.velocity_util.
hyp_2F1
(a, b, c, z)[source]¶ http://docs.sympy.org/0.7.1/modules/mpmath/functions/hypergeometric.html
-
lenstronomy.GalKin.velocity_util.
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
-
lenstronomy.GalKin.velocity_util.
moffat_r
(r, alpha, beta)[source]¶ Moffat profile
Parameters: - r – radial coordinate
- alpha – Moffat parameter
- beta – exponent
Returns: Moffat profile
-
lenstronomy.GalKin.velocity_util.
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
-
lenstronomy.GalKin.velocity_util.
draw_moffat_r
(FWHM, beta)[source]¶ Parameters: - FWHM – full width at half maximum
- beta – Moffat beta parameter
Returns: draw from radial Moffat distribution
-
lenstronomy.GalKin.velocity_util.
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
-
lenstronomy.GalKin.velocity_util.
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: