Lensing¶
The following functions/classes provide functionality for simulating gravitational lensing by the matter distribution in the universe.
Iterative lensing¶
- class glass.MultiPlaneConvergence(cosmo)[source]¶
Compute convergence fields iteratively from multiple matter planes.
- Parameters:
cosmo (Cosmology)
- glass.multi_plane_matrix(shells, cosmo)[source]¶
Compute the matrix of lensing contributions from each shell.
- Parameters:
shells (
Sequence
[RadialWindow
]) – The shells of the mass distribution.cosmo (
Cosmology
) – Cosmology instance.
- Returns:
The matrix of lensing contributions.
- Return type:
- glass.multi_plane_weights(weights, shells, cosmo)[source]¶
Compute effective weights for multi-plane convergence.
Converts an array weights of relative weights for each shell into the equivalent array of relative lensing weights.
This is the discretised version of the integral that turns a redshift distribution \(n(z)\) into the lensing efficiency sometimes denoted \(g(z)\) or \(q(z)\).
- Parameters:
weights (
ndarray
[Any
,dtype
[float64
]]) – Relative weight of each shell. The first axis must broadcast against the number of shells, and is normalised internally.shells (
Sequence
[RadialWindow
]) – Window functions of the shells.cosmo (
Cosmology
) – Cosmology instance.
- Returns:
The relative lensing weight of each shell.
- Raises:
ValueError – If the shape of weights does not match the number of shells.
- Return type:
Lensing fields¶
- glass.from_convergence(kappa, lmax=None, *, potential=False, deflection=False, shear=False, discretized=True)[source]¶
Compute other weak lensing maps from the convergence.
Takes a weak lensing convergence map and returns one or more of deflection potential, deflection, and shear maps. The maps are computed via spherical harmonic transforms.
- Parameters:
kappa (
ndarray
[Any
,dtype
[float64
]]) – HEALPix map of the convergence field.lmax (
int
|None
) – Maximum angular mode number to use in the transform.potential (
bool
) – Which lensing maps to return.deflection (
bool
) – Which lensing maps to return.shear (
bool
) – Which lensing maps to return.discretized (
bool
) – Correct the pixel window function in output maps.
- Returns:
psi – Map of the lensing (or deflection) potential. Only returned if
potential
is true.alpha – Map of the deflection (complex). Only returned if
deflection
if true.gamma – Map of the shear (complex). Only returned if
shear
is true.
- Return type:
Notes
The weak lensing fields are computed from the convergence or deflection potential in the following way. [1]
Define the spin-raising and spin-lowering operators of the spin-weighted spherical harmonics as
\[\eth {}_sY_{lm} = +\sqrt{(l-s)(l+s+1)} \, {}_{s+1}Y_{lm} \;, \\ \bar{\eth} {}_sY_{lm} = -\sqrt{(l+s)(l-s+1)} \, {}_{s-1}Y_{lm} \;.\]The convergence field \(\kappa\) is related to the deflection potential field \(\psi\) by the Poisson equation,
\[2 \kappa = \eth\bar{\eth} \, \psi = \bar{\eth}\eth \, \psi \;.\]The convergence modes \(\kappa_{lm}\) are hence related to the deflection potential modes \(\psi_{lm}\) as
\[2 \kappa_{lm} = -l \, (l+1) \, \psi_{lm} \;.\]The deflection \(\alpha\) is the gradient of the deflection potential \(\psi\). On the sphere, this is
\[\alpha = \eth \, \psi \;.\]The deflection field has spin weight \(1\) in the HEALPix convention, in order for points to be deflected towards regions of positive convergence. The modes \(\alpha_{lm}\) of the deflection field are hence
\[\alpha_{lm} = \sqrt{l \, (l+1)} \, \psi_{lm} \;.\]The shear field \(\gamma\) is related to the deflection potential \(\psi\) and deflection \(\alpha\) as
\[2 \gamma = \eth\eth \, \psi = \eth \, \alpha \;,\]and thus has spin weight \(2\). The shear modes \(\gamma_{lm}\) are related to the deflection potential modes as
\[2 \gamma_{lm} = \sqrt{(l+2) \, (l+1) \, l \, (l-1)} \, \psi_{lm} \;.\]References
- [1] Tessore N., et al., OJAp, 6, 11 (2023).
doi:10.21105/astro.2302.01942
- glass.shear_from_convergence(kappa, lmax=None, *, discretized=True)[source]¶
Weak lensing shear from convergence.
Computes the shear from the convergence using a spherical harmonic transform.
Deprecated since version 2023.6: Use the more general
from_convergence()
function instead.
Applying lensing¶
- glass.deflect(lon, lat, alpha)[source]¶
Apply deflections to positions.
Takes an array of deflection values and applies them to the given positions.
- Parameters:
lon (
float
|ndarray
[Any
,dtype
[float64
]]) – Longitudes to be deflected.lat (
float
|ndarray
[Any
,dtype
[float64
]]) – Latitudes to be deflected.alpha (
complex
|list
[float
] |ndarray
[Any
,dtype
[complex128
]] |ndarray
[Any
,dtype
[float64
]]) – Deflection values. Must be complex-valued or have a leading axis of size 2 for the real and imaginary component.
- Returns:
The longitudes and latitudes after deflection.
- Return type:
tuple
[ndarray
[Any
,dtype
[float64
]],ndarray
[Any
,dtype
[float64
]]]
Notes
Deflections on the sphere are defined as follows: The complex deflection \(\alpha\) transports a point on the sphere an angular distance \(|\alpha|\) along the geodesic with bearing \(\arg\alpha\) in the original point.
In the language of differential geometry, this function is the exponential map.