lensit.ffs_deflect

lensit.ffs_deflect.ffs_deflect.displacement_fromplm(lib_plm, plm, **kwargs)[source]

Returns deflection-field instance from lensing gradient potential

Parameters:
  • lib_plmffs_alm or ffs_alm_pyFFTW instance describing the plm array
  • plm – flat-sky lensing gradient potential alm array
Returns:

Deflection field as ffs_displacement instance

class lensit.ffs_deflect.ffs_deflect.ffs_displacement(dx, dy, lsides, LD_res=(11, 11), verbose=False, NR_iter=3, lib_dir=None, cache_magn=False)[source]

Flat-sky deflection-field class

Used to perform lensing on maps and to obtain the inverse deflection for iterative lensing estimation

Parameters:
  • dx – deflection field, x-component \(\alpha_x\) (2d-array or path to array on disk)
  • dy – deflection field, y-component \(\alpha_y\) (2d-array or path to array on disk)
  • lsides (tuple) – physical size in radians of the flat-sky patch.
  • LD_res (optional) – to perform inversion or lensing large maps are split into chunks sized these powers of two
  • verbose (optional) – various prints out
  • NR_iter (optional) – Number of Newton-Raphson iterations in deflection inversion. Default works very well for LCDM-like deflection fields
  • cache_magn (optional) – optionally caches magnification determinant matrix when needed
  • lib_dir (optional) – required only if cache_magn is set
alm2lenmap(lib_alm, alm, use_Pool=0, crude=0)[source]

Return deflected position-space map from its unlensed input harmonic coeffients.

Parameters:
  • lib_almlensit.ffs_covs.ell_mat.ffs_alm instance adapted to input alm array
  • alm – input unlensed flat-sky alm array
Returns:

position space map of shape lib_alm.shape

get_Omega()[source]

Lensing (curl) potential map

get_det_magn()[source]

Returns magnification determinant map

\(|M| = \det \begin{pmatrix} 1 + \frac{\partial \alpha_x}{\partial x} & \frac{\partial \alpha_x}{\partial y} \\ \frac{\partial \alpha_y}{\partial x} & 1 + \frac{\partial \alpha_y}{\partial y} \end{pmatrix}\)

get_inverse(NR_iter=None, use_Pool=0, crude=0, HD_res=None)[source]

Builds deflection field inverse

The inverse deflection is defined by the condition that forward-deflected points remap to the original points under the inverse displacement.

Parameters:
  • NR_iter (optional) – Number of Newton-Raphson iterations. Superseeds the instance NR_iter argument if set
  • use_Pool (optional) – Send the calculation to the GPU if negative
  • crude (optional) – Uses some faster but crude method for the inverse (check get_inverse_crude)
  • HD_res (optional) – Set this to perform the inversion at higher resolution (tuple of powers of two)

Returns: inverse deflection field as ffs_displacement instance

get_inverse_crude(crude)[source]

Crude inversions of the displacement field

Parameters:crude – method key

Supported now is only crude = 1, for which the inverse is approximated as the negative deflection

get_kappa()[source]

Convergence map.

Returns:\(\kappa = -\frac 12 \left( \frac{\partial \alpha_x}{\partial x} + \frac{\partial \alpha_y}{\partial y} \right)\)
get_noisefreemf(lib_qlm)[source]

Deflection induced mean-field estimate for a noisefree, isotropic experimental setting

get_omega()[source]

Field rotation map

Returns:\(\omega = \frac 12 \left( \frac{\partial \alpha_x}{\partial y} - \frac{\partial \alpha_y}{\partial x} \right)\)
get_phi()[source]

Lensing (gradient) potential map

lens_alm(lib_alm, alm, lib_alm_out=None, mult_magn=False, use_Pool=0)[source]

Returns lensed harmonic coefficients from the unlensed input coefficients

Parameters:
  • lib_almlensit.ffs_covs.ell_mat.ffs_alm instance adapted to input alm array
  • alm – input unlensed flat-sky alm array
  • lib_alm_out (optional) – output ffs_alm instance if difference from input
  • mult_magn (optional) – optionally multiplies the real-space lensed map with the magnification det. if set.
  • use_Pool (optional) – calculations are performed on the GPU if negative.
Returns:

lensed alm array

lens_map(m, use_Pool=0, crude=0, do_not_prefilter=False)[source]

Lens the input flat-sky map, using a bicubic spline interpolation algorithm

The task is split in chunks (of typically (2048 * 2048) or specified by the LD_res parameters) with a buffer size to ensure the junctions are properly performed.

Parameters:
  • m – real-space map to deflect. numpy array, or the path to an array on disk.
  • use_Pool (optional) – set this to < 0 to perform the operation on the GPU
  • crude (optional) – uses an alternative crude approximation to lensing if set (check lens_map_crude)
  • do_not_prefilter (optional) – sidesteps the bicubic interpolation prefiltering step. Only use this if you know what you are doing

Returns: deflected real-space map (array)

lens_map_crude(m, crude)[source]

Performs crude approximations to lens operation

Parameters:
  • m – map to lens (2d-array of the right shape)
  • crude – approximation method key

Now supported are crude = 1 (nearest pixel rouding) or 2 (first order series expansion in deflection)

class lensit.ffs_deflect.ffs_deflect.ffs_id_displacement(shape, lsides)[source]

Displacement instance where there is actually no displacement.

Mostly useful for testing purposes