lenstronomy.Cosmo package¶
Submodules¶
lenstronomy.Cosmo.background module¶
-
class
Background
(cosmo=None, interp=False, **kwargs_interp)[source]¶ Bases:
object
Class to compute cosmological distances.
-
T_xy
(z_observer, z_source)[source]¶ Parameters: - z_observer – observer
- z_source – source
Returns: transverse comoving distance in units of Mpc
-
__init__
(cosmo=None, interp=False, **kwargs_interp)[source]¶ Parameters: - cosmo – instance of astropy.cosmology
- interp – boolean, if True, uses interpolated cosmology to evaluate specific redshifts
- kwargs_interp – keyword arguments of CosmoInterp specifying the interpolation interval and maximum redshift
Returns: Background class with instance of astropy.cosmology
-
static
a_z
(z)[source]¶ Returns scale factor (a_0 = 1) for given redshift.
Parameters: z – redshift Returns: scale factor
-
d_xy
(z_observer, z_source)[source]¶ Parameters: - z_observer – observer redshift
- z_source – source redshift
Returns: angular diameter distance in units of Mpc
-
ddt
(z_lens, z_source)[source]¶ Time-delay distance.
Parameters: - z_lens – redshift of lens
- z_source – redshift of source
Returns: time-delay distance in units of proper Mpc
-
rho_crit
¶ Critical density.
Returns: value in M_sol/Mpc^3
-
lenstronomy.Cosmo.cosmo_solver module¶
-
cosmo2angular_diameter_distances
(H_0, omega_m, z_lens, z_source)[source]¶ Parameters: - H_0 – Hubble constant [km/s/Mpc]
- omega_m – dimensionless matter density at z=0
- z_lens – deflector redshift
- z_source – source redshift
Returns: angular diameter distances Dd and Ds/Dds
-
ddt2h0
(ddt, z_lens, z_source, cosmo)[source]¶ Converts time-delay distance to H0 for a given expansion history.
Parameters: - ddt – time-delay distance in Mpc
- z_lens – deflector redshift
- z_source – source redshift
- cosmo – astropy.cosmology class instance
Returns: h0 value which matches the cosmology class effectively replacing the h0 value used in the creation of this class
-
class
SolverFlatLCDM
(z_d, z_s)[source]¶ Bases:
object
Class to solve multidimensional non-linear equations to determine the cosmological parameters H0 and omega_m given the angular diameter distance relations.
-
class
InvertCosmo
(z_d, z_s, H0_range=None, omega_m_range=None)[source]¶ Bases:
object
Class to do an interpolation and call the inverse of this interpolation to get H_0 and omega_m.
lenstronomy.Cosmo.kde_likelihood module¶
-
class
KDELikelihood
(D_d_sample, D_delta_t_sample, kde_type='scipy_gaussian', bandwidth=1)[source]¶ Bases:
object
Class that samples the cosmographic likelihood given a distribution of points in the 2-dimensional distribution of D_d and D_delta_t.
-
__init__
(D_d_sample, D_delta_t_sample, kde_type='scipy_gaussian', bandwidth=1)[source]¶ Parameters: - D_d_sample – 1-d numpy array of angular diameter distances to the lens plane
- D_delta_t_sample – 1-d numpy array of time-delay distances
- kde_type (string) – The kernel to use. Valid kernels are ‘scipy_gaussian’ or [‘gaussian’|’tophat’|’epanechnikov’|’exponential’|’linear’|’cosine’] Default is ‘gaussian’.
- bandwidth – width of kernel (in same units as the angular diameter quantities)
-
logLikelihood
(D_d, D_delta_t)[source]¶ Likelihood of the data (represented in the distribution of this class) given a model with predicted angular diameter distances.
Parameters: - D_d – model predicted angular diameter distance
- D_delta_t – model predicted time-delay distance
Returns: loglikelihood (log of KDE value)
-
lenstronomy.Cosmo.lcdm module¶
-
class
LCDM
(z_lens, z_source, flat=True)[source]¶ Bases:
object
Flat LCDM cosmology background with free Hubble parameter and Omega_m at fixed lens redshift configuration.
-
D_d
(H_0, Om0, Ode0=None)[source]¶ Angular diameter to deflector.
Parameters: - H_0 – Hubble parameter [km/s/Mpc]
- Om0 – normalized matter density at present time
Returns: float [Mpc]
-
D_ds
(H_0, Om0, Ode0=None)[source]¶ Angular diameter from deflector to source.
Parameters: - H_0 – Hubble parameter [km/s/Mpc]
- Om0 – normalized matter density at present time
Returns: float [Mpc]
-
D_dt
(H_0, Om0, Ode0=None)[source]¶ Time-delay distance.
Parameters: - H_0 – Hubble parameter [km/s/Mpc]
- Om0 – normalized matter density at present time
Returns: float [Mpc]
-
lenstronomy.Cosmo.lens_cosmo module¶
-
class
LensCosmo
(z_lens, z_source, cosmo=None)[source]¶ Bases:
object
Class to manage the physical units and distances present in a single plane lens with fixed input cosmology.
-
__init__
(z_lens, z_source, cosmo=None)[source]¶ Parameters: - z_lens – redshift of lens
- z_source – redshift of source
- cosmo – astropy.cosmology instance
-
arcsec2phys_lens
(arcsec)[source]¶ Convert angular to physical quantities for lens plane.
Parameters: arcsec – angular size at lens plane [arcsec] Returns: physical size at lens plane [Mpc]
-
arcsec2phys_source
(arcsec)[source]¶ Convert angular to physical quantities for source plane.
Parameters: arcsec – angular size at source plane [arcsec] Returns: physical size at source plane [Mpc]
-
dd
¶ Returns: angular diameter distance to the deflector [Mpc]
-
dds
¶ Returns: angular diameter distance from deflector to source [Mpc]
-
ddt
¶ Returns: time delay distance [Mpc]
-
ds
¶ Returns: angular diameter distance to the source [Mpc]
-
h
¶
-
hernquist_angular2phys
(sigma0, rs_angle)[source]¶ ‘sigma0’ is defined such that the deflection at projected RS leads to alpha = 2./3 * Rs * sigma0.
Parameters: - sigma0 – convergence normalization
- rs_angle – rs in angular units [arcseconds]
Returns: mass [M_sun], rs [Mpc]
-
hernquist_phys2angular
(mass, rs)[source]¶ Translates physical mass definitions of the Hernquist profile to the angular units used in the Hernquist lens profile of lenstronomy.
‘sigma0’ is defined such that the deflection at projected RS leads to alpha = 2./3 * Rs * sigma0
Parameters: - mass – A spherical overdensity mass in M_sun corresponding to the mass definition mdef at redshift z
- rs – rs in units of physical Mpc
Returns: sigma0, Rs_angle
-
kappa2proj_mass
(kappa)[source]¶ Convert convergence to projected mass M_sun/Mpc^2.
Parameters: kappa – lensing convergence Returns: projected mass [M_sun/Mpc^2]
-
mass_in_coin
(theta_E)[source]¶ Parameters: theta_E – Einstein radius [arcsec] Returns: mass in coin calculated in mean density of the universe
-
mass_in_theta_E
(theta_E)[source]¶ Mass within Einstein radius (area * epsilon crit) [M_sun]
Parameters: theta_E – Einstein radius [arcsec] Returns: mass within Einstein radius [M_sun]
-
nfwParam_physical
(M, c)[source]¶ Returns the NFW parameters in physical units.
Parameters: - M – physical mass in M_sun in definition m200
- c – concentration
Returns: rho0 [Msun/Mpc^3], Rs [Mpc], r200 [Mpc]
-
nfw_M_theta_r200
(M)[source]¶ Returns r200 radius in angular units of arc seconds on the sky.
Parameters: M – physical mass in M_sun Returns: angle (in arc seconds) of the r200 radius
-
nfw_angle2physical
(Rs_angle, alpha_Rs)[source]¶ Converts the angular parameters into the physical ones for an NFW profile.
Parameters: - alpha_Rs – observed bending angle at the scale radius in units of arcsec
- Rs_angle – scale radius in units of arcsec
Returns: rho0 [Msun/Mpc^3], Rs [Mpc], c, r200 [Mpc], M200 [Msun]
-
nfw_physical2angle
(M, c)[source]¶ Converts the physical mass and concentration parameter of an NFW profile into the lensing quantities.
Parameters: - M – mass enclosed 200 rho_crit in units of M_sun (physical units, meaning no little h)
- c – NFW concentration parameter (r200/r_s)
Returns: Rs_angle (angle at scale radius) (in units of arcsec), alpha_Rs (observed bending angle at the scale radius
-
phys2arcsec_lens
(phys)[source]¶ Convert physical Mpc into arc seconds.
Parameters: phys – physical distance [Mpc] Returns: angular diameter [arcsec]
-
sersic_k_eff2m_star
(k_eff, R_sersic, n_sersic)[source]¶ Translates convergence at half-light radius to total integrated physical stellar mass for a Sersic profile.
Parameters: - k_eff – lensing convergence at half-light radius
- R_sersic – half-light radius in arc seconds
- n_sersic – Sersic index
Returns: stellar mass in physical Msun
-
sersic_m_star2k_eff
(m_star, R_sersic, n_sersic)[source]¶ Translates a total stellar mass into ‘k_eff’, the convergence at ‘R_sersic’ (effective radius or half-light radius) for a Sersic profile.
Parameters: - m_star – total stellar mass in physical Msun
- R_sersic – half-light radius in arc seconds
- n_sersic – Sersic index
Returns: k_eff
-
sigma_crit
¶ Returns the critical projected lensing mass density in units of M_sun/Mpc^2.
Returns: critical projected lensing mass density
-
sigma_crit_angle
¶ Returns the critical surface density in units of M_sun/arcsec^2 (in physical solar mass units) when provided a physical mass per physical Mpc^2.
Returns: critical projected mass density
-
sis_sigma_v2theta_E
(v_sigma)[source]¶ Converts the velocity dispersion into an Einstein radius for a SIS profile.
Parameters: v_sigma – velocity dispersion (km/s) Returns: theta_E (arcsec)
-
sis_theta_E2sigma_v
(theta_E)[source]¶ Converts the lensing Einstein radius into a physical velocity dispersion.
Parameters: theta_E – Einstein radius (in arcsec) Returns: velocity dispersion in units (km/s)
-
time_delay2fermat_pot
(dt)[source]¶ Parameters: dt – time delay in units of days Returns: Fermat potential in units arcsec**2 for a given cosmology
-
time_delay_units
(fermat_pot, kappa_ext=0)[source]¶ Parameters: - fermat_pot – in units of arcsec^2 (e.g. Fermat potential)
- kappa_ext – unit-less external shear not accounted for in the Fermat potential
Returns: time delay in days
-
uldm_angular2phys
(kappa_0, theta_c)[source]¶ Converts the anguar parameters entering the LensModel Uldm() (Ultra Light Dark Matter) class in physical masses, i.e. the total soliton mass and the mass of the particle.
Parameters: - kappa_0 – central convergence of profile
- theta_c – core radius (in arcseconds)
Returns: m_eV_log10, M_sol_log10, the log10 of the masses, m in eV and M in M_sun
-
uldm_mphys2angular
(m_log10, M_log10)[source]¶ Converts physical ULDM mass in the ones, in angular units, that enter the LensModel Uldm() class.
Parameters: - m_log10 – exponent of ULDM mass in eV
- M_log10 – exponent of soliton mass in M_sun
Returns: kappa_0, theta_c, the central convergence and core radius (in arcseconds)
-
lenstronomy.Cosmo.nfw_param module¶
-
class
NFWParam
(cosmo=None)[source]¶ Bases:
object
Class which contains a halo model parameters dependent on cosmology for NFW profile All distances are given in physical units.
Mass definitions are relative to 200 crit including redshift evolution. The redshift evolution is cosmology dependent (dark energy). The H0 dependence is propagated into the input and return units.
-
static
M200
(rs, rho0, c)[source]¶ Calculation of the mass enclosed r_200 for NFW profile defined as.
\[M_{200} = 4 \pi \rho_0^{3} * \left(\log(1+c) - c / (1 + c) \right))\]Parameters: - rs (float) – scale radius
- rho0 (float) – density normalization (characteristic density) in units mass/[distance unit of rs]^3
- c (float [4,40]) – concentration
Returns: M(R_200) mass in units of rho0 * rs^3
-
M_r200
(r200, z)[source]¶ Parameters: - r200 – r200 in physical Mpc/h
- z – redshift
Returns: M200 in M_sun/h
-
static
c_M_z
(M, z)[source]¶ fitting function of http://moriond.in2p3.fr/J08/proceedings/duffy.pdf for the mass and redshift dependence of the concentration parameter
Parameters: - M (float or numpy array) – halo mass in M_sun/h
- z (float >0) – redshift
Returns: concentration parameter as float
-
c_rho0
(rho0, z)[source]¶ Computes the concentration given density normalization rho_0 in h^2/Mpc^3 (physical) (inverse of function rho0_c)
Parameters: - rho0 – density normalization in h^2/Mpc^3 (physical)
- z – redshift
Returns: concentration parameter c
-
nfw_Mz
(M, z)[source]¶ Returns all needed parameter (in physical units modulo h) to draw the profile of the main halo r200 in physical Mpc/h rho_s in h^2/Mpc^3 (physical) Rs in Mpc/h physical c unit less.
Parameters: - M – Mass in physical M_sun/h
- z – redshift
-
r200_M
(M, z)[source]¶ Computes the radius R_200 crit of a halo of mass M in physical mass M/h.
Parameters: - M (float or numpy array) – halo mass in M_sun/h
- z (float) – redshift
Returns: radius R_200 in physical Mpc/h
-
rho0_c
(c, z)[source]¶ Computes density normalization as a function of concentration parameter.
Parameters: - c – concentration
- z – redshift
Returns: density normalization in h^2/Mpc^3 (physical)
-
rhoc
= 277536627000.0¶
-
static