grale.lenses

This module contains classes to work with different kinds of gravitational lens models.

The core class is GravitationalLens, but this should not be instantiated directly. Instead allocate a class derived from this; currently available are

LensException

exception grale.lenses.LensException

This exception is raised if something goes wrong in the GravitationalLens derived classes.

Helper functions

grale.lenses.getCriticalDensity(Dd, Ds, Dds)

Returns the critical density corresponding to angular diameter distances Dd (observer to lens), Ds (observer to source) and Dds (lens to source).

This is defined as:

Sigma_crit = c^2/(4*pi*G*Dd) * (Ds/Dds)
grale.lenses.createLensFromLenstoolFile(inputData, mirrorX=False)

Based on a LensTool model, a corresponding lens model is constructed. The function returns a tuple consisting of the lens model, the lens redshift and the cosmological model. An illustration can be found in the notebook lenstooltest.ipynb

Note that this is preliminary code, and currently only the PIEMD model (LensTool model type 81) is handled.

Arguments:
  • inputData: the data from a LensTool file (typically with ‘.par’ extension), either as a file name, a file object or a string.

  • mirrorX: if True, the model will be mirrored along the X-axis.

grale.lenses.createEquivalentPotentialGridLens(lens, bottomLeft, topRight, NX, NY, maskRegions, potentialGradientWeight, densityGradientWeight, densityWeight, pixelEnlargements=2, enlargeDiagonally=False, circleToPolygonPoints=10000, feedbackObject='default', qpsolver='scs', laplacianKernel=array([[0., 0., 1., 0., 0.], [0., 1., 2., 1., 0.], [1., 2., -16., 2., 1.], [0., 1., 2., 1., 0.], [0., 0., 1., 0., 0.]]), maxDensityConstraints=[], exactDensityConstraints=[], exactDeflectionConstraints=[], exactDeflectionTolerance=0, ignorePixelMismatch=False, ignorePositiveDensityConstraint=False, exceptionOnFail=True)

This uses a quadratic programming approach to extrapolate the lens potential values in certain regions (typically covering the images in a lensing system), thereby creating a lens that has the same effect (because the lens potential is the same in the image regions).

Examples can be found in the msdexample-equivlenstests.ipynb and potentialextrap_multisheet.ipynb notebooks.

Arguments:
  • lens: the procedure will start from the lens potential values of this GravitationalLens model. It will sample the lens potential on a grid, and keep some values fixed, determined by a mask. The other lens potential values of the new lens model will be extrapolated.

  • bottomLeft, topRight, NX, NY: these values determine the grid on which the lens potential will be sampled from the original lens, and which will be used to define the new lens (a PotentialGridLens). See createThetaGrid.

  • maskRegions: this is a list of region descriptions that will be combined to create a binary mask of NX by NY that indicates which lens potential values should be kept. This is passed to the createThetaGridAndImagesMask function.

  • potentialGradientWeight, densityGradientWeight, densityWeight: the quadratic programming problem tries to optimize a combination of three parts: one for the gradient of the lens potential, one for the gradient of the resulting mass density, and one for the mass density itself. These weights are used to specify their respective contributions.

  • pixelEnlargements, enlargeDiagonally, circleToPolygonPoints: these are passed on to the similarly named arguments of createThetaGridAndImagesMask.

  • feedbackObject: can be used to specify a particular feedback mechanism.

  • qpsolver: for the quadratic programming optimization, the qpsolvers module is used, which itself allows for different solver implementations to be used. This argument specifies the name of the solver, and is passed as the solver argument to qpsolvers.solve_qp.

  • laplacianKernel: to go from the lens potential values to the density, a convolution with this 2D kernel is performed. The default is

    \[\begin{split}\left[ \begin{array}{cccc} 0 & 0 & 1 & 0 & 0 \\ 0 & 1 & 2 & 1 & 0 \\ 1 & 2 & -16 & 2 & 1 \\ 0 & 1 & 2 & 1 & 0 \\ 0 & 0 & 1 & 0 & 0 \end{array} \right]\end{split}\]
  • maxDensityConstraints: since this kernel allows the density to be calculated at the grid points, constraints can be specified for those. This is a list of these constraints where each entry is a dictionary with the following keys:

    • maskRegions: a list of regions to which the specified density applies. Is passed on to the createThetaGridAndImagesMask again, without any pixel enlargements.

    • density: can either be a scalar value or a grid with the same dimensions as the lens potential grid.

    • upperlimit (optional): if omitted, this is interpreted as True, and means that the specified density in the specified region should not be exceeded. Set this to False for a lower bound on the density instead.

  • exactDensityConstraints: similar to maxDensityConstraints, but requests the exact mass density values in certain regions. This is again a list of similar dictionaries, but in this case only with maskRegions and density keys.

  • exactDeflectionConstraints: instead of providing constraints on the density, constraints on the deflection field (the gradient of the lens potential) can be set as well. This is a list of dictionaries with entries

    • maskRegions: again a list of regions, this time for which the specified deflection field applies.

    • lens: get the desired deflection field at the grid points from this lens.

    Alternatively, instead of the lens key, ax and ay keys may be present, each describing one component of the deflection field. In case these are specified manually, it will be necessary to take into account that the kernel used (e.g. [-1, 1]) gives an approximation in between grid points.

  • exactDeflectionTolerance: using the constraints above as exact constraints rarely works. Instead it is possible to request that this deflection field is obtained within some tolerance, e.g. 0.1 arcsec, on each side of the exact value.

  • ignorePixelMismatch: if the specified grid parameters does not yield the same spacing between grid points in x- and y-direction, the routine will abort by default. You can override this and continue anyway by setting this flag.

  • ignorePositiveDensityConstraint: by default, a constraint is added that requires the density resulting from the lens potential to be positive everywhere. To ignore this condition, set this flag.

  • exceptionOnFail: by default, if the quadratic programming solver is unable to come up with a solution, an exception is thrown. Setting this to False disables this and returns the dictionary mentioned below anyway, but with the resulting solution and resulting lens set to None.

The return value is a dictionary with the following keys:

  • philens_equiv: this the the main result of the optimization routine, a PotentialGridLens that’s based on lens potential values defined on the grid that you specified, where some values were taken from the input lens, and others were optimized according to the weights and constraints that were specified.

  • philens_orig: as an intermediate step, the input lens is approximated by sampling the lens potential on the grid points, and this is the approximate lens model. You can compare this to the original model to see how much sense the approximation makes.

  • mask: the main maskRegions argument, together with other parameters (e.g. pixelEnlargements) leads to this binary mask. Where it’s True, the lens potential values are fixed, the others will be optimized.

  • masksMaxDens: each entry of maxDensityConstraints has its own maskRegions argument, leading to a binary mask again. This entry contains a list of these masks.

  • masksExactDens: similar, but for the exactDensityConstraints.

  • masksDeflection: similar, but for the exactDeflectionConstraints.

  • P, q, G, h, A, b: these are the matrices that are created for the quadratic programming problem

  • initvals: the initial values that are passed to the solver, these are the lens potential values of the input lens.

  • x: this contains the solution to the quadratic programming problem.

GravitationalLens

class grale.lenses.GravitationalLens

A base class for a gravitational lens.

static fromBytes(bytes b)

Instantiates a gravitational lens object from the byte array previously generated by toBytes().

getAlphaVector(thetas)

Returns the deflection angles at the positions in ‘thetas’.

getAlphaVectorDerivatives(thetas)

Returns the derivatives of the deflection angle at the positions specified in ‘thetas’. For each position (thetax, thetay), three values (axx, ayy, axy) are returned, where

\[\begin{split}\begin{array}{l} {\rm axx} = \frac{\partial \hat{\alpha}_x}{\partial \theta_x} \\ {\rm ayy} = \frac{\partial \hat{\alpha}_y}{\partial \theta_y} \\ {\rm axy} = \frac{\partial \hat{\alpha}_x}{\partial \theta_y} = \frac{\partial \hat{\alpha}_y}{\partial \theta_x} \end{array}\end{split}\]
getCLLensProgram(derivatives=True, potential=True)

Returns a concatenation of getCLLensQuantitiesStructure() and getCLProgram()

getCLLensQuantitiesStructure(derivatives=True, potential=True)

The OpenCL routine returned by getCLProgram() returns a structure, for which the code can be queried by this function.

getCLParameterCounts()

Returns a tuple containing the number of integer parameters and floating point parameters this lens’s OpenCL implementation needs.

getCLParameters(deflectionscale, potentialscale)

Returns a tuple of integer and floating point parameters to use this lens in OpenCL, expressed using the specified deflection and potential scales.

getCLProgram(derivatives=True, potential=True)

Returns an OpenCL program to calculate the deflection and optionally deflection derivatives and lensing potential, for this lens model.

getCriticalDensity(Ds, Dds)

Calls getCriticalDensity() where Dd is the angular diameter distance of this gravitational lens (see getLensDistance())

getInverseMagnification(Ds, Dds, thetas)

Returns the inverse magnification values at the positions in ‘thetas’, for a source plane with angular diameter distances Ds and Dds.

getLensDistance()

Returns the angular diameter distance Dd that was specified when creating this gravitational lens object.

getLensParameters()

Obtain the parameters (or equivalent) that were used to create this gravitational lens instance.

getProjectedPotential(Ds, Dds, thetas)

Returns the projected potential for a source plane with angular diameter distances Ds and Dds, measured at the positions in ‘thetas’.

getRadialDensityProfile(thetaRadii)

If the lens is a symmetric lens, this function returns the density profile at the specified radii.

getRadialMassProfile(thetaRadii)

If the lens is a symmetric lens, this function returns the radial mass profile at the specified radii.

getSuggestedScales()

Returns a dictionary with suggested deflection and potential scales, which can be used in the OpenCL routines to be able to calculate in single precision floating point.

getSurfaceMassDensity(thetas)

Returns the surface mass density of the gravitational lens at the positions in ‘thetas’.

getSurfaceMassDensityMap(bottomLeft, topRight, numX, numY, feedbackObject='default', renderer='default', reduceToPixels=True)

Returns the surface density map for the area between the bottomLeft and topRight coordinates. The surface mass density is sampled at numX and numY points and these values are returned directly if reduceToPixels is False. The grid can also be seen as a grid of pixels, and if reduceToPixels is True, the values will be averaged to obtain (numX-1)*(numY-1) pixel values.

The feedbackObject specifies the method which is used for feedback when calculating the values as the grid points. See grale.feedback for more information.

If renderer is None (the default), calculations will occur in the same process, on a single CPU core. For other possible values, see the grale.renderers module.

getTimeDelay(zd, Ds, Dds, thetas, beta)

Returns the time delay (up to an unknown constant) for this gravitational lens with redshift ‘zd’, for a source position ‘beta’ in a source plane specified by the angular diameter distances ‘Ds’ and ‘Dds’, and for an image plane positions ‘thetas’.

Note that for this function, the ‘thetas’ positions do not necessarily need to trace back to ‘beta’. This allows, for example, to calculate the time delays for a grid of image plane points (and one fixed ‘beta’) and identify the locations where images are seen as the points where the time delay reaches a local optimum.

static load(fileName)

Loads the gravitational lens that was previously saved (see save()) to the file with name ‘fileName’.

save(fileName)

Write this gravitational lens object to the file with name ‘fileName’. Can be loaded again using load().

setDerivativeAngularDistanceScale(scale)

In case the getAlphaVectorDerivatives() function is not provided for a specific lens, the derivatives will be calculated numerically using points that are this distance apart.

setLensDistance(Dd)

This function changes the angular diameter distance to the lens to Dd.

toBytes()

Converts the lens object to a byte array, which can be converted to a real lens object again using fromBytes()

traceTheta(Ds, Dds, thetas)

Traces the image plane positions in ‘thetas’ to the source plane that corresponds to the angular diameter distances Ds and Dds.

GaussLens

class grale.lenses.GaussLens

A gravitational lens with a gaussian profile

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘mass’: the total mass of the lens

    • ‘width’: the (angular) width of the lens

MultiplePlummerLens

class grale.lenses.MultiplePlummerLens

A gravitational lenses consisting of a number of Plummer basis functions

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a list of which each entry is a dictionary with the following entries:

    • ‘mass’: the total mass of this Plummer basis function

    • ‘width’: the (angular) width of this basis function

    • ‘x’: the X-position of the center of this basis function

    • ‘y’: the Y-position of the center of this basis function

PlummerLens

class grale.lenses.PlummerLens

A gravitational lens corresponding to a projected Plummer sphere

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘mass’: the total mass of this Plummer lens

    • ‘width’: the (angular) width of the Plummer profile

PointmassLens

class grale.lenses.PointmassLens

A point mass gravitational lens

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘mass’: the total mass of this point mass lens

SISLens

class grale.lenses.SISLens

A singular isothermal sphere (SIS) gravitational lens

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘velocityDispersion’: the velocity dispersion of the SIS lens

NSIELens

class grale.lenses.NSIELens

A non-singular isothermal ellipse (NSIE) gravitational lens

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘velocityDispersion’: the velocity dispersion of the NSIE lens

    • ‘ellipticity’: the ellipticity parameter for this lens

    • ‘coreRadius’: the (angular) radius of the core of the NSIE profile

NSISLens

class grale.lenses.NSISLens

A non-singular isothermal sphere (NSIS) gravitational lens

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘velocityDispersion’: the velocity dispersion of the NSIS lens

    • ‘coreRadius’: the (angular) radius of the core of the NSIS profile

SIELens

class grale.lenses.SIELens

A singular isothermal ellipse (SIE) gravitational lens

\[ \begin{align}\begin{aligned}\Sigma(\vec{\theta}) = \frac{\sigma_v^2\sqrt{f}}{2 G D_d\sqrt{\theta_x^2+f^2\theta_y^2}}\\\vec{\hat{\alpha}}(\vec{\theta}) = \frac{4 \pi \sigma_v^2}{c^2} \frac{\sqrt{f}}{\sqrt{1-f^2}} \left[ \mathrm{asinh}\left(\frac{\sqrt{1-f^2}}{f}\frac{\theta_x}{|\vec{\theta}|}\right)\vec{e}_x + \mathrm{asin}\left(\sqrt{1-f^2}\frac{\theta_y}{|\vec{\theta}|}\right)\vec{e}_y \right]\end{aligned}\end{align} \]

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘velocityDispersion’: the velocity dispersion of the SIE lens

    • ‘ellipticity’: the ellipticity parameter for this lens

SquareLens

class grale.lenses.SquareLens

The gravitational lens effect of a square shaped region of contant surface mass density

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘mass’: the total mass of this square shaped lens

    • ‘width’: the (angular) width of the lens

MultipleSquareLens

class grale.lenses.MultipleSquareLens

A gravitational lenses consisting of a number of square shaped basis functions

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a list of which each entry is a dictionary with the following entries:

    • ‘mass’: the total mass of this square basis function

    • ‘width’: the (angular) width of this basis function

    • ‘x’: the X-position of the center of this basis function

    • ‘y’: the Y-position of the center of this basis function

MultipleGaussLens

class grale.lenses.MultipleGaussLens

A gravitational lenses consisting of a number of gaussian basis functions

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a list of which each entry is a dictionary with the following entries:

    • ‘mass’: the total mass of this gaussian basis function

    • ‘width’: the (angular) width of this basis function

    • ‘x’: the X-position of the center of this basis function

    • ‘y’: the Y-position of the center of this basis function

MassSheetLens

class grale.lenses.MassSheetLens

Gravitational lens effect of “n infinitely large mass sheet of constant surface mass density

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have either the following entry:

    • ‘density’: the surface mass density of the lens

    or the following entries:

    • ‘Ds’ : the angular diameter distance between observer and source plane

    • ‘Dds’ : the angular diameter distance between lens and source plane

    In the second case, the critical density corresponding to Dd, Dds and Ds is used.

CompositeLens

class grale.lenses.CompositeLens

This models a gravitational lens that consists of other known lens models

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a list of which each entry is a dictionary with the following entries:

    • ‘lens’: the gravitational lens object of this component

    • ‘factor’: the effect of this component can be amplified or reduced by setting this factor. Set to 1 to use the component as-is.

    • ‘angle’: the angle in degrees that the component should be rotated (counter clockwise)

    • ‘x’: the center of the component is now placed at this X-position

    • ‘y’: the center of the component is now placed at this Y-position

Warning! No check is done to make sure that the same Dd parameter from the constructor is also used in the individual components. Make sure that this is the case, or deal with undefined behaviour!

findCLSubroutines(derivatives, potential)

Analyzes the current CompositeLens instance, returns a tuple of:

  • the recursion level needed (how many other CompositeLens levels it contains)

  • a dictionary where the keys are subroutine names needed by the full program (e.g. a plummer lens function) and the values are the actual OpenCL programs for these subroutines

  • an array where each entry is either the subroutine name for a sublens or None if not needed.

The recursion level and array can be fed into the getCompositeCLProgram() function.

static getCompositeCLProgram(otherRoutineNames, maxRecursion, derivatives, potential)

Returns the name of the OpenCL function as well a the OpenCL program itself for a CompositeLens where maxRecursion levels are needed, and calls to the subroutine names in otherRoutineNames need to be present. These two parameters are returned by the findCLSubroutines() function. The code for these other subroutines is not included.

MassDiskLens

class grale.lenses.MassDiskLens

This models the gravitational lens effect of a disk with a constant density that’s centered on the origin.

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have either the following entry:

    • ‘density’: the surface mass density of the lens

    • ‘radius’: the angular radius of the disk

    or the following entries:

    • ‘Ds’ : the angular diameter distance between observer and source plane

    • ‘Dds’ : the angular diameter distance between lens and source plane

    • ‘radius’: the angular radius of the disk

    In the second case, the critical density corresponding to Dd, Dds and Ds is used.

PolynomialMassProfileLens

class grale.lenses.PolynomialMassProfileLens

With this gravitational lens model, a circularly symmetric lens centered on the origin can be modeled based on a specified enclosed mass profile. This profile can be specified by different pieces, each one being a polynomial.

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a list of which each entry is a dictionary with the following entries:

    • ‘xoffset’: see below

    • ‘yoffset’: see below

    • ‘xscale’: see below

    • ‘yscale’: see below

    • ‘xend’: see below

    • ‘coeffs’: see below

    The ‘xend’ values specify the x values (the \(\theta\) radii) where each polynomial part is valid: the current profile part is used if x is between the previous list entry’s ‘xend’ value (or 0 if it doesn’t exist) and the current ‘xend’ value.

    If it is used, the mass profile’s value M(x) will be calculated as:

    M(x) = yscale * ( ... + coeffs[k] * ((x-xoffset)/scale)^k + ... ) + yoffset
    

TimeDelayAdjustLens

class grale.lenses.TimeDelayAdjustLens

A circularly symmetric lens with which the timedelay in the central part can be changed by a specific amount. An example can be found in the notebook timedelayadjust.ipynb

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘theta1’: angular radius withing which the lensing potential has a constant value, thereby adjusting the time delay.

    • ‘theta2’: the angular radius beyond which the lensing potential is zero. In between the potential will change smoothly.

    • ‘dt’: the amount with which the time delay in the central region should be modified, determining the value of the lensing potential within ‘theta1’.

    • ‘z’: the redshift of the lens.

ZeroMassLens

class grale.lenses.ZeroMassLens

This produces a symmetric gravitational lens with zero total mass, as described in the article below.

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have the following entries:

    • ‘density’: the surface mass density at the center of the lens

    • ‘radius’: the angular radius beyond which the total mass of the lens will be zero

    • ‘zeropoint’: a value between 0 and 1 that specified the position (as a fraction of the radius) where the surface mass density becomes zero and the negative density begins.

MassDiskLensSmoothed

class grale.lenses.MassDiskLensSmoothed

This is the same a a mass disk (MassDiskLens), but with a smoother edge.

The density profile is:

\[\begin{split}\Sigma(\theta) = \left\{\begin{array}{lr} \Sigma_0 & \theta < \theta_1 \\ \Sigma_0 \left[ 2\left(\frac{\theta-\theta_1}{\theta_2-\theta_1}\right)^3 -3\left(\frac{\theta-\theta_1}{\theta_2-\theta_1}\right)^2 +1 \right] & \theta \in [\theta_1, \theta_2] \\ 0 & \theta > \theta_2 \end{array}\right.\end{split}\]
__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params is a dictionary that should have either the following entry:

    • ‘density’: the surface mass density of the lens

    • ‘radius’: the angular radius of the disk itself

    • ‘endradius’: the angular radius where the density becomes zero

    or the following entries:

    • ‘Ds’ : the angular diameter distance between observer and source plane

    • ‘Dds’ : the angular diameter distance between lens and source plane

    • ‘radius’: the angular radius of the disk itself

    • ‘endradius’: the angular radius where the density becomes zero

    In the second case, the critical density corresponding to Dd, Dds and Ds is used.

MultipleWendlandLens

class grale.lenses.MultipleWendlandLens

This creates a lens that’s composed of the Wendland basis functions described in the LensPerfect article. You can either specify the basis functions yourself, or you can specify a list of points at which certain deflection angles should be present in the resulting lens.

An example notebook can be found here: wendland.ipynb, the lens file used in the notebook is the following: reallens_nosheet.lensdata

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: should be one of two cases:

    1. a dictionary with entries ‘phix’ and ‘phiy’. Each of those should be a list of dictionaries, each describing a Wendland basis function that’s oriented along the x- or y-axis. These dictionaries should have the entries:

      • ‘weight’: a weight for the basis function

      • ‘scale’: an angular scale for the basis function

      • ‘position’: the angular position of the basis function

    2. a dictionary containing entries

      • ‘points’: a list of angular positions at which the deflection angle will be specified

      • ‘angles’: a list of equal length as the ‘points’ list, describing the deflection angle the lens should have at that point.

      • ‘scale’: an angular scale of the basis functions

      If this method is used, the basis functions will be determined automatically in such a way that they generate the specified deflection angles at the specified points, where each basis function has the angular scale given by ‘scale’.

DeflectionGridLens

class grale.lenses.DeflectionGridLens

Creates a lens based on provided deflection angles, provided on a regular grid. When a deflection angle is later retrieved using e.g. getAlphaVector, each component will be interpolated bilinearly based on the provided deflections. This type of lens is what is created when LensPlane.createDeflectionGridLens is called.

Note that this is only an approximation and should be used with some limitations in mind. For one thing, it will not have a valid lens potential associated with it. Also, as the deflection angles rely on interpolation, they will in general not be curl-free.

An example notebook can be found here: deflectiongridlens.ipynb, the lens file used in the notebook is the following: reallens_nosheet.lensdata

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘angles’: a NumY*NumX*2 numpy array containing the deflection angles.

    • ‘bottomleft’: the bottom-left coordinate, which will have the deflection angle stored at position (0, 0) in the ‘angles’ grid.

    • ‘topright’: the top-right coordinat, which will have the deflection angle stored at position (NumY-1, NumX-1) in the ‘angles’ grid.

NFWLens

class grale.lenses.NFWLens

A lens model for a 2D symmetric lens profile, based on a 3D symmetric Navarro-Frenk-White (NFW) mass distribution.

We’ll be using the following definitions of the helper functions \(F(x)\), \(G(x)\), \(H(x)\), \(A(x)\) and \(B(x)\) as follows

\[\begin{split}F(x) = \left\{ \begin{array}{ll} \frac{1}{\sqrt{1-x^2}}\textrm{atanh}{\sqrt{1-x^2}} & x < 1 \\ 1 & x = 1 \\ \frac{1}{\sqrt{x^2-1}}\textrm{atan}{\sqrt{x^2-1}} & x > 1 \\ \end{array} \right. \Rightarrow \frac{d F}{d x} = \frac{1-x^2 F(x)}{x(x^2-1)}\end{split}\]
\[G(x) = \frac{1-F(x)}{x^2-1} \qquad H(x) = \textrm{ln}\left(\frac{x}{2}\right) + F(x) \Rightarrow \frac{d H}{d x} = G(x) x\]
\[\begin{split}B(x) = \left\{ \begin{array}{ll} \textrm{ln}^2\left(\frac{x}{2}\right) - \textrm{atanh}^2\sqrt{1-x^2} & x < 1 \\ \textrm{ln}^2 2 & x = 1 \\ \textrm{ln}^2\left(\frac{x}{2}\right) + \textrm{atan}^2\sqrt{x^2-1} & x > 1 \\ \end{array} \right. \qquad A(x) = \frac{d H}{d x} \Rightarrow \frac{1}{2} \frac{d B}{d x} = A(x)\end{split}\]

The following relations then hold for the 2D surface mass density \(\Sigma(\theta)\) and the 2D integrated mass profile \(\textrm{M}(\theta)\)

\[\Sigma(\theta) = 2 r_s \rho_s G\left(\frac{\theta}{\theta_s}\right) \qquad \textrm{M}(\theta) = 4\pi r_s^3\rho_s H\left(\frac{\theta}{\theta_s}\right)\]

Being a symmetric lens, this \(\textrm{M}(\theta)\) then leads to the deflection angles

\[\vec{\hat{\alpha}}\left(\vec{\theta}\right) = \frac{16 \pi G r_s^2 \rho_s}{c^2} A\left(\frac{\theta}{\theta_s}\right) \frac{\vec{\theta}}{\theta}\]

while the (unscaled) lensing potential becomes

\[\psi(\theta) = \frac{D_{ds}}{D_s} \frac{8 \pi G r_s^2 \rho_s \theta_s}{c^2} B\left(\frac{\theta}{\theta_s}\right)\]

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘rho_s’: a density scale for the NFW profile, corresponding to \(\rho_s\) in the equations above

    • ‘theta_s’: an angular scale factor for the NFW profile, corresponding to \(\theta_s\) in the equations above

EllipticNFWLens

class grale.lenses.EllipticNFWLens

An elliptic generalization of the NFWLens lens model. The relations between the circularly symmetric profile and the elliptic generalization are described in the reference.

The degree of ellipticity is encoded by a parameter \(q\), so that for the 2D mass density the relation holds:

\[\Sigma(\theta_x, \theta_y) = \Sigma_{\rm circular}\left(\theta_x, \frac{\theta_y}{q}\right)\]

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘rho_s’: a density scale for the NFW profile, corresponding to \(\rho_s\) in the equations in NFWLens

    • ‘theta_s’: an angular scale factor for the NFW profile, corresponding to \(\theta_s\) in the equations in NFWLens

    • ‘q’: a measure of the ellipticity of the lens, where 1 describes again the circular profile.

SersicLens

class grale.lenses.SersicLens

A lens model with a 2D projected mass density based on the Sérsic profile. This mass density is the following:

\[\Sigma(\theta) = \Sigma_{\rm central}\exp\left[\left(\frac{\theta}{\theta_s}\right)^{\frac{1}{n}}\right]\]
__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘centraldensity’: the value of the 2D density at the center of the lens corresponding to \(\Sigma_{\rm central}\) in the equation above

    • ‘scale’: an angular scale of the lens, corresponding to \(\theta_s\) in the equation above

    • ‘index’: the so-called Sérsic index, corresponding to \(n\) in the equation above

EllipticSersicLens

class grale.lenses.EllipticSersicLens

An elliptic generalization of the SersicLens lens model. The relations between the circularly symmetric profile and the elliptic generalization are described in the reference.

The degree of ellipticity is encoded by a parameter \(q\), so that for the 2D mass density the relation holds:

\[\Sigma(\theta_x, \theta_y) = \Sigma_{\rm circular}\left(\theta_x, \frac{\theta_y}{q}\right)\]

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘centraldensity’: the value of the 2D density at the center of the lens corresponding to \(\Sigma_{\rm central}\) in the equation in SersicLens

    • ‘scale’: an angular scale of the lens, corresponding to \(\theta_s\) in the equation in SersicLens

    • ‘index’: the so-called Sérsic index, corresponding to \(n\) in the equation in SersicLens

    • ‘q’: a measure of the ellipticity of the lens, where 1 describes again the circular profile.

ProfileLens

class grale.lenses.ProfileLens

A circularly symmetric lens based on a discrete density profile. In between the points of the profile, the density is linearly interpolated.

An example notebook can be found here: profilelens.ipynb.

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘radius’: the profile will contain density values for evenly spaced points between zero and this end radius. Beyond this radius, the model does not contain any mass.

    • ‘profile’: a list of density values, the first of which corresponds to the central density of the circularly symmetric lens, and the last corresponds to the density at the specified end radius.

PIEMDLens

class grale.lenses.PIEMDLens

A lens with a PIEMD profile with core and cut radius. Although this is usually referred to as a PIEMD, it is not a PIEMD as described in the article of Kassiola and Kovner. In Eliasdottir et al, is is called a dPIE to make the distinction more clear.

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘centraldensity’: the central projected density \(\Sigma_0\)

    • ‘coreradius’: the core radius \(a\) of the profile

    • ‘scaleradius’: the scale radius \(s\) of the profile

    • ‘epsilon’: a measure for the ellipticity of the lens

PIMDLens

class grale.lenses.PIMDLens

The circularly symmetric version of the PIEMD lens (PIEMDLens), corresponding to en ‘epsilon’ parameter of zero.

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘centraldensity’: the central projected density \(\Sigma_0\)

    • ‘coreradius’: the core radius \(a\) of the profile

    • ‘scaleradius’: the scale radius \(s\) of the profile

AlphaPotLens

class grale.lenses.AlphaPotLens

This lens is based on the ‘alphapot’ softened power law lens model from lenstool/gravlens. The unscaled (\(D_{ds}/D_s = 1\)) lensing potential is taken to be

\[\psi(\vec{\theta}) = b\left(s^2 + \theta_x^2 + \frac{\theta_y^2}{q^2} + K^2\theta_x\theta_y\right)^{\frac{\alpha}{2}}\]

Note that because \(\kappa(\vec{\theta}) = \frac{1}{2}\left(\frac{\partial^2 \psi}{\partial \theta_x^2} + \frac{\partial^2 \psi}{\partial \theta_y^2}\right)\) you’ll need to do an additional scaling (e.g. by adjusting \(b\) or using a CompositeLens) by \((1 \textrm{ arcsec})^{\alpha-2}\) to get the same mass density as lenstool (there, 1 arcsec is used as angular unit).

References

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘b’: the value of \(b\) in lensing potential

    • ‘s’: the value of \(s\) in lensing potential

    • ‘q’: the value of \(q\) in lensing potential

    • ‘K2’: the value of \(K^2\) in lensing potential

    • ‘alpha’: the value of \(\alpha\) in the lensing potential

HarmonicLens

class grale.lenses.HarmonicLens

Harmonic 2D mass density, similar to what’s used in the JPEG article.

\[\Sigma(\vec{\theta}) = \Sigma_0 \cos\left(\frac{k}{2}\theta_x + \phi_x\right) \cos\left(\frac{l}{2}\theta_y + \phi_y\right)\]
References:
__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary containing the following entries:

    • ‘sigma0’: Maximum density \(\Sigma_0\)

    • ‘k’: value of ‘k’ in the formula

    • ‘l’: value of ‘l’ in the formula

    • ‘phi_x’: value of \(\phi_x\) in the formula

    • ‘phi_y’: value of \(\phi_y\) in the formula

PotentialGridLens

class grale.lenses.PotentialGridLens

Create a lens based on the values of the projected potential (for \(D_{ds}/D_s = 1\)) defined on a grid. In between the grid points, bicubic interpolation is used using routines from GSL.

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary consisting of the following entries

    • ‘bottomleft’: the coordinates of the bottom-left corner of the grid

    • ‘topright’: the top-right corner of the grid

    • ‘values’: 2D numpy array containing the potential values, sampled at the grid points. If the shape is (numY, numX) (first index describes y-direction) then point [0, 0] corresponds to the bottom-left corner, and [numY-1, numX-1] to the top-right corner.

CircularPiecesLens

class grale.lenses.CircularPiecesLens

Create a lens based on the projected potentials of existing lenses. These potentials will be used exactly in ring-shaped regions (possibly with a different offset, and rescaled), and interpolated in between. If \(\psi_i\) is the projected potential for one such component, the scaled/offset version is

\[\psi_i^s = \mathrm{scale}\times(\psi_i + \mathrm{offset})\]

In between such regions, the potentials are merged in the following way:

\[\psi_i^{\mathrm{interp}}(\vec{\theta}) = f(x) \psi_i^s(\vec{\theta}) + (1-f(x)) \psi_{i+1}^s(\vec{\theta})\]

where

\[x = \frac{|\vec{\theta}| - \theta_i^{\rm end}}{\theta_{i+1}^{\rm start} - \theta_i^{\rm end}}\]

and \(\theta_i^{\rm end}\) is the outer radius for the ring-like region of potential \(i\), and \(\theta_{i+1}^{\rm start}\) is the inner radius of the next ring-like region.

__init__(Dd, params)
Parameters:
  • Dd is the angular diameter distance to the lens.

  • params: a dictionary with two entries

    • ‘coeffs’: (optional) a list of numbers that represent the coefficients of the interpolation function

      \[f(x) = \sum_k a_k x^k\]

      The default is [ 1, 0, 0, -10, 15, -6  ], similar to a smoothstep-like function.

    • ‘pieces’: a list of dictionaries with the following entries, representing the existing lens potentials and their ring-like regions:

      • ‘lens’: lens model

      • ‘r0’: start radius (zero if not present)

      • ‘r1’: end radius (infinity if not present)

      • ‘potentialOffset’: offset to the lens potential for this lens (zero if not present)

      • ‘potentialScale’: scale factor for the lens potential for this lens (one if not present)

MultiPlaneContainer

class grale.lenses.MultiPlaneContainer

This ‘lens’ can be used as a container for multiple lenses as different redshifts

__init__(Dd, params)
Parameters:
  • Dd is not used in this container, must be set to 0

  • params is a list of which each entry is a dictionary with the following entries:

    • ‘lens’: the gravitational lens object at this redshift

    • ‘z’: the redshift for this lens

Warning! No check is done to make sure that the lens distance parameter matches with the specified redshift (there’s no notion of a cosmological model here)