Welcome to lensit’s documentation!¶
pip install -e . [—user]
Flat-sky python quadratic and iterative CMB lensing simulation package by Julien Carron.
lensit¶
This lensit package contains some convenience functions in its __init__.py for quick startup.
-
lensit.
get_fidcls
(ellmax_sky=6000)[source]¶ Returns lensit fiducial CMB spectra (Planck 2015 cosmology)
Parameters: ellmax_sky – optionally reduces outputs spectra \(\ell_{\rm max}\) Returns: unlensed and lensed CMB spectra (dicts)
-
lensit.
get_lencmbs_lib
(res=14, cache_sims=True, nsims=120, num_threads=1)[source]¶ Default lensed CMB simulation library
Lensing is always performed at resolution of \(0.75\) arcmin
Parameters: - res – lensed CMBs are generated on a square box with of physical size \(\sim 0.74 \cdot 2^{\rm res}\) arcmin
- cache_sims – saves the lensed CMBs when produced for the first time
- nsims – number of simulations in the library
- num_threads – number of threads used by the pyFFTW fft-engine.
Note
All simulations random phases will be generated at the very first call if not performed previously; this might take some time
-
lensit.
get_maps_lib
(exp, LDres, HDres=14, cache_lenalms=True, cache_maps=False, nsims=120, num_threads=1)[source]¶ Default CMB data maps simulation library
Parameters: - exp – experimental configuration (see get_config)
- LDres – the data is generated on a square patch with :math:` 2^{rm LDres}` pixels on a side
- HDres – The physical size of the path is \(\sim 0.74 \cdot 2^{\rm HDres}\) arcmin
- cache_lenalms – saves the lensed CMBs when produced for the first time (defaults to True)
- cache_maps – saves the data maps when produced for the first time (defaults to False)
- nsims – number of simulations in the library
- num_threads – number of threads used by the pyFFTW fft-engine.
Note
All simulations random phases (CMB sky and noise) will be generated at the very first call if not performed previously; this might take some time
lensit.ffs_iterators¶
-
class
lensit.ffs_iterators.ffs_iterator.
ffs_iterator
(lib_dir, typ, filt, dat_maps, lib_qlm, Plm0, H0, cpp_prior, use_Pool_lens=0, use_Pool_inverse=0, chain_descr=None, opfilt=None, soltn0=None, cache_magn=False, no_deglensing=False, NR_method=100, tidy=10, verbose=True, maxcgiter=150, PBSSIZE=None, PBSRANK=None, **kwargs)[source]¶ Flat-sky iterator template class
Parameters: - lib_dir – many things will be written there
- typ – ‘T’, ‘QU’ or ‘TQU’ for estimation on temperature data, polarization data or jointly
- filt – inverse-variance filtering instance (e.g. lensit.qcinv.ffs_ninv_filt )
- dat_maps – data maps or path to maps.
- lib_qlm – lib_alm (lensit.ffs_covs.ell_mat.ffs_alm) instance describing the lensing estimate Fourier arrays
- Plm0 – Starting point for the iterative search. alm array consistent with lib_qlm
- H0 – initial isotropic likelihood curvature approximation (roughly, inverse lensing noise bias \(N^{(0)}_L\))
- cpp_prior – fiducial lensing power spectrum, used for the prior part of the posterior density.
- chain_descr – multigrid conjugate gradient inversion chain description
-
build_incr
(it, key, gradn)[source]¶ Search direction
BGFS method with ‘self.NR method’ BFGS updates to the Hessian. Initial Hessian are built from N0s. It must be rank 0 here.
Parameters: - it – current iteration level. Will produce the increment to phi_{k-1}, from gradient est. g_{k-1} phi_{k_1} + output = phi_k
- key – ‘p’ or ‘o’
- gradn – current estimate of the gradient (alm array)
Returns: increment for next iteration (alm array)
-
calc_gradplikpdet
(it, key)[source]¶ Calculates the likelihood gradient (quadratic and mean-field parts)
-
get_Gaussnoisesample
(it, key, plm_noisephas, real_space=False, verbose=False)[source]¶ Produce a Gaussian random field from the approximate BFGS covariance
Parameters: - it – iteration index
- key – ‘p’ or ‘o’ for lensing gradient or curl iteration
- plm_noisepha – unit spectra random phases of the right shape
- real_space – produces random field in real space if set, otherwise alm array
-
iterate
(it, key, cache_only=False)[source]¶ Performs an iteration
This builds the gradients at iteration it, and the potential estimate, and saves the it + 1 estimate.
-
load_graddet
(it, key)[source]¶ Loads mean-field gradient at iteration it
Gradient must have already been calculated
-
load_gradpri
(it, key)[source]¶ Loads prior gradient at iteration it
Gradient must have already been calculated
-
class
lensit.ffs_iterators.ffs_iterator.
ffs_iterator_cstMF
(lib_dir, typ, filt, dat_maps, lib_qlm, Plm0, H0, MF_qlms, cpp_prior, **kwargs)[source]¶ Iterator instance, that uses fixed, input mean-field at each step.
Parameters: - lib_dir – many things will be written there
- typ – ‘T’, ‘QU’ or ‘TQU’ for estimation on temperature data, polarization data or jointly
- filt – inverse-variance filtering instance (e.g. lensit.qcinv.ffs_ninv_filt )
- dat_maps – data maps or path to maps.
- lib_qlm – lib_alm (lensit.ffs_covs.ell_mat.ffs_alm) instance describing the lensing estimate Fourier arrays
- Plm0 – Starting point for the iterative search. alm array consistent with lib_qlm
- H0 – initial isotropic likelihood curvature approximation (roughly, inverse lensing noise bias \(N^{(0)}_L\))
- MF_qlms – mean-field alm array (also desribed by lib_qlm)
- cpp_prior – fiducial lensing power spectrum, used for the prior part of the posterior density.
-
class
lensit.ffs_iterators.ffs_iterator.
ffs_iterator_pertMF
(lib_dir, typ, filt, dat_maps, lib_qlm, Plm0, H0, cpp_prior, init_rank=0, init_barrier=<function <lambda>>, **kwargs)[source]¶ Iterator instance, that uses the deflection-perturbative prediction for the mean-field at each step.
Parameters: - lib_dir – many things will be written there
- typ – ‘T’, ‘QU’ or ‘TQU’ for estimation on temperature data, polarization data or jointly
- filt – inverse-variance filtering instance (e.g. lensit.qcinv.ffs_ninv_filt )
- dat_maps – data maps or path to maps.
- lib_qlm – lib_alm (lensit.ffs_covs.ell_mat.ffs_alm) instance describing the lensing estimate Fourier arrays
- Plm0 – Starting point for the iterative search. alm array consistent with lib_qlm
- H0 – initial isotropic likelihood curvature approximation (roughly, inverse lensing noise bias \(N^{(0)}_L\))
- cpp_prior – fiducial lensing power spectrum, used for the prior part of the posterior density.
-
class
lensit.ffs_iterators.ffs_iterator.
ffs_iterator_simMF
(lib_dir, typ, MFkey, nsims, filt, dat_maps, lib_qlm, Plm0, H0, cpp_prior, **kwargs)[source]¶ Iterator instance, that estimate the mean-field at each steps from Monte-Carlos.
Parameters: - lib_dir – many things will be written there
- typ – ‘T’, ‘QU’ or ‘TQU’ for estimation on temperature data, polarization data or jointly
- MFkey – mean-field estimator key
- nsims – number of sims to use at each step
- filt – inverse-variance filtering instance (e.g. lensit.qcinv.ffs_ninv_filt )
- dat_maps – data maps or path to maps.
- lib_qlm – lib_alm (lensit.ffs_covs.ell_mat.ffs_alm) instance describing the lensing estimate Fourier arrays
- Plm0 – Starting point for the iterative search. alm array consistent with lib_qlm
- H0 – initial isotropic likelihood curvature approximation (roughly, inverse lensing noise bias \(N^{(0)}_L\))
- cpp_prior – fiducial lensing power spectrum, used for the prior part of the posterior density.
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_plm – ffs_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_alm – lensit.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_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)\)
-
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_alm – lensit.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)
lensit.ffs_covs.ell_mat¶
-
class
lensit.ffs_covs.ell_mat.
ell_mat
(lib_dir, shape, lsides, cache=1)[source]¶ Library helping with flat-sky patch discretization and harmonic mode structure.
This handles Fourier mode structure on the flat sky, at the given resolution and size of specified rectangular box.
Parameters: - lib_dir – various things might be cached there (unlens cache is set to 0, in which case only a hashdict is put there at instantiation)
- shape (2-tuple) – pair of int defining the number of pixels on each side of the box
- lsides (2-tuple) – physical size (in radians) of the box sides
- cache (optional) – if non-zero, a bunch of matrices might be cached to speed up some calculations.
-
alm2rfftmap
(alm, filt_func=<function ell_mat.<lambda>>)[source]¶ Returns a ffs_alm array from a rfftmap.
Pseudo-inverse to fftmap2alm
-
bin_inell
(rfftmap)[source]¶ Bin in a rfftmap according to multipole \(\ell\)
Input must have the same shape then self.rshape
-
degrade
(LDshape, lib_dir=None)[source]¶ Returns an ell_mat instance on the same physical flat-sky patch but a degraded sampling resolution
-
get_ellmat
(ellmax=None)[source]¶ Returns the matrix containing the multipole ell assigned to 2d-frequency \(k = (k_x,k_y)\)
-
get_phasemat
(ellmax=None)[source]¶ Returns the rfft array containing the phase \(\phi\) where 2d-frequencies are \(k = |k| e^{i \phi}\)
-
get_pixwinmat
()[source]¶ Pixel window function rfft array \(\sin(k_x l_{\rm cell, x} / 2) \sin (k_y l_{\rm cell, y} / 2 )\)
-
get_unique_ells
(return_index=False, return_inverse=False, return_counts=False)[source]¶ Returns list of multipoles \(\ell\) present in the flat-sky patch
-
static
k2ell
(k)[source]¶ Mapping of 2d-frequency \(k\) to multipole \(\ell\)
\(\ell = \rm{int}\left(|k| - \frac 12 \right)\)
-
class
lensit.ffs_covs.ell_mat.
ffs_alm
(ellmat, filt_func=<function ffs_alm.<lambda>>, kxfilt_func=None, kyfilt_func=None)[source]¶ Library to facilitate operations on flat-sky alms in the ffs scheme.
Methods includes alm2map and map2alm, among other things. The instance contains info on the sky patch size and discretization.
Parameters: - ellmat – lensit.ffs_covs.ell_mat.ell_mat instance carrying flat-sky patch definitions
- filt_func (callable, optional) – mutlipole filtering. Set this to discard multipoles in map2alm Defaults to lambda ell : ell > 0
- kxfilt_func (callable, optional) – filtering on x-component of the 2d-frequencies
- kyfilt_func (callable, optional) – filtering on y-component of the 2d-frequencies
-
alm2Pk_minimal
(alm)[source]¶ Power spectrum estimation on the grid, with minimal binning (only exact same frequencies).
Many bins will have only number count 4 or something.
Returns: list with integer array of squared frequencies, integer array number counts, and Pk estimates. (only non zero counts frequencies)
-
almxfl
(alm, fl, inplace=False)[source]¶ Multiply :flat-sky math:a_{ell lm} array with isotropic function \(f_\ell\)
-
bicubic_prefilter
(alm)[source]¶ Refilter the map for use in bicubic spline prefiltering.
Mostly useful for lensing of maps.
-
bin_realpart_inell
(alm, ellmax=None)[source]¶ Bin real part of a rfftmap according to multipole \(\ell\)
Input must have the same shape then self.rshape
lensit.ffs_covs.ffs_cov¶
-
class
lensit.ffs_covs.ffs_cov.
ffs_diagcov_alm
(lib_dir, lib_datalm, cls_unl, cls_len, cl_transf, cls_noise, lib_skyalm=None, init_rank=0, init_barrier=<function <lambda>>)[source]¶ Library for flat-sky calculations of various lensing biases, responses, etc. in a idealized, isotropic case
Parameters: - lib_dir – many things will be saved there
- lib_datalm – lib_alm instance (see lensit.ffs_covs.ell_mat containing mode filtering and flat-sky patch info
- cls_unl (dict) – unlensed CMB cls
- cls_len (dict) – lensed CMB cls
- cl_transf – instrument transfer function
- cls_noise (dict) – ‘t’, ‘q’ and ‘u’ noise arrays
- lib_skyalm (optional) – lib_alm instance describing the sky mode. Irrelevant with some exceptions. Defaults to lib_datalm
-
get_N0cls
(typ, lib_qlm, use_cls_len=True, cls_obs=None, cls_obs2=None, cls_weights=None, cls_filt=None)[source]¶ Lensing Gaussian bias \(N^{(0)}_L\)
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance decribing the lensing alm arrays
- cls_obs (dict, optional) – empirical observed data spectra, for a realization-dependent estimate
- cls_weights (dict, optional) – CMB cls entering the QE weights (numerator)
- cls_filt (dict, optional) – CMB cls entering the inverse-variance filtering step (denominator of QE weights)
- use_cls_len (optional) – Uses lensed or unlensed CMB cls (when not superseeded byt the other keywords)
Returns: Gradient and curl lensing mode noise \(N^{(0)}_L\)
-
get_RDdelensingcorrbias
(typ, lib_qlm, alwfcl, clsobs_deconv, clsobs_deconv2=None, cls_weights=None)[source]¶ Calculate delensing bias given a reconstructed potential map spectrum
Same as get_delensingcorrbias using empirical input CMB spectra
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- alwfcl – normalization of the Wiener-filter quadratic estimate (i.e. Wiener-filter times inverse response)
- clsobs_deconv (dict) – emprical beam-deconvolved CMB data spectra.
Returns: (3, 3, lmax +1) array with bias \(C_\ell^{ij}, \textrm{ with }i,j \in (T,E,B)\)
-
get_delensingcorrbias
(typ, lib_qlm, alwfcl, CMBonly=False)[source]¶ Calculate delensing bias given a reconstructed potential map spectrum
Crudely, the delensing bias is defined as \(C_\ell^{\rm lensed} - C_\ell^{\rm delensed}\) on Gaussian CMB maps (with lensed power spectrum). More precisely this method calculates the 4-point disconnected contribution (see https://arxiv.org/abs/1701.01712) from the statistical dependence of the lensing tracer reconstruction noise to the CMB maps
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- alwfcl – normalization of the Wiener-filter quadratic estimate (i.e. Wiener-filter times inverse response)
- CMBonly – do not include noise spectra if set (defaults to False)
Returns: (3, 3, lmax +1) array with bias \(C_\ell^{ij}, \textrm{ with }i,j \in (T,E,B)\)
-
get_delensinguncorrbias
(lib_qlm, clpp_rec, wNoise=True, wCMB=True, recache=False, use_cls_len=True)[source]¶ Calculate delensing bias given a reconstructed potential map spectrum
Crudely, the delensing bias is defined as \(C_\ell^{\rm lensed} - C_\ell^{\rm delensed}\) on Gaussian CMB maps (with lensed power spectrum). More precisely this method calculates the 4-point disconnected contribution (see https://arxiv.org/abs/1701.01712) if no statistical dependence of the lensing tracer reconstruction noise to the CMB maps
Parameters: - lib_qlm – ffs_alm instance describing the lensing alm arrays
- clpp_rec – lensing tracer power. For the true bias calculation clpp_rec should be reconstruction noise power only. For perturbative lensing calculation, include the signal power in clpp_rec
- wNoise – include noise spectra in total observed CMB spectra (defaults to True)
- wCMB – include signal CMB spectra in total observed CMB spectra (defaults to True)
- use_cls_len – use lensed CMB cls if set, unlensed if not (defaults to True)
Returns: (3, 3, lmax +1) array with bias \(C_\ell^{ij}, \textrm{ with }i,j \in (T,E,B)\)
-
get_dfishertrace
(typ, cmb_dcls, lib_qlm, K1=None, K2=None, dK1=None, dK2=None, use_cls_len=True, recache=False)[source]¶ Variation of Fisher-trace get_fishertrace with respect to input CMB variation
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis cmb_dcls (dict): CMB spectra variation
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- use_cls_len – use lensed or unlensed cls in QE weights (numerator), defaults to lensed cls
-
get_dlndetcurv
(typ, cmb_dcls, lib_qlm, K=None, dK=None, use_cls_len=False)[source]¶ Variation of get_lndetcurv with respect to CMB spectra variation
get_lndetcurv is the harmonic transform of \(\frac 1 2 \rm{Tr} A \frac{\delta^2 \rm{Cov}}{\delta \alpha(x) \delta \alpha(y)}\) for vanishing deflection
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- cmb_dcls (dict) – CMB spectra variation
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- use_cls_len – use lensed or unlensed cls in QE weights (numerator), defaults to lensed cls
Returns: gradient and curl Fourier components
-
get_dmfrespcls
(typ, cmb_dcls, lib_qlm, use_cls_len=False)[source]¶ Variation of the linear mean-field response \(R_L\) of the unnormalized quadratic estimators
Variation of get_mfrespcls given input CMB spectra variation
This deflection-induced mean-field linear response is
\(\frac{\delta g^{\rm MF}(x)}{\delta \alpha(y)} = \frac 12 \frac{\delta \ln \rm{Cov}}{\delta \alpha(x) \delta \alpha(y)} = R(x- y)\)
See https://arxiv.org/abs/1704.08230
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- cmb_dcls (dict) – CMB spectra variations
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- use_cls_len – use lensed CMB cls if true, unlensed if not (default)
-
get_dplmRDlikcurvcls
(typ, cmb_dcls, datcls_obs, lib_qlm, use_cls_len=True, recache=False, dat_only=False)[source]¶ derivative of plmRDlikcurvcls (data held fixed) Finite difference test OK (+ much faster)
-
get_fishertrace
(typ, lib_qlm, get_A1=None, get_A2=None, use_cls_len=True)[source]¶ Fisher-matrix like trace
\(\frac 1 2 \textrm{Tr}\: A_1 \:\frac{\delta\rm{Cov}}{\delta\phi} \:A_2 \:\frac{\delta\rm{Cov}}{\delta\phi}\)
if operators \(A_1\) or \(A_2\) not set, they reduce to \(\rm{Cov}^{-1}\), and the result is the Fisher info for vanishing deflection
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- get_A1 (optional, callable) – 1st trace operator
- get_A2 (optional, callable) – 2nd trace operator
- use_cls_len – use lensed or unlensed cls in QE weights (numerator), defaults to lensed cls
-
get_iblms
(typ, datalms, use_cls_len=True, use_Pool=0, **kwargs)[source]¶ Inverse-variance filters input CMB maps
Produces \(B^t \rm{Cov}^{-1} X^{\rm dat}\) (inputs to quadratic estimator routines) This instance applies isotropic filtering, with no deflection field
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- datalms – data alms array
- use_cls_len – use lensed cls in filtering (defaults to True)
Returns: inverse-variance filtered alms (T and/or Q, Ulms)
-
get_lndetcurv
(typ, lib_qlm, get_A=None, use_cls_len=False)[source]¶ (Part of) covariance log-det second variation
Harmonic transform of \(\frac 1 2 \rm{Tr} A \frac{\delta^2 \rm{Cov}}{\delta \alpha(x) \delta \alpha(y)}\) for vanishing deflection
if get_A argument is not set, output becomes \(\frac 1 2 \rm{Tr} \rm{Cov}^{-1} \delta^2 \rm{Cov}\)
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- get_A (optional, callable) – 1st trace operator. Defaults to \(\rm{Cov}^{-1}\)
- use_cls_len – use lensed or unlensed cls in QE weights (numerator), defaults to lensed cls
Returns: gradient and curl Fourier components
-
get_mfrespcls
(typ, lib_qlm, use_cls_len=True)[source]¶ Linear mean-field response \(R_L\) of the unnormalized quadratic estimators (gradient and curl components)
This deflection-induced mean-field linear response is
\(\frac{\delta g^{\rm MF}(x)}{\delta \alpha(y)} = \frac 12 \frac{\delta \ln \rm{Cov}}{\delta \alpha(x) \delta \alpha(y)} = R(x- y)\)
See https://arxiv.org/abs/1704.08230
The expected normalised MF spectrum is then \(\frac{R_L^2}{\mathcal R_L^2} C_L^{\hat g\hat g}\)
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- use_cls_len – use lensed CMB cls if true, unlensed if not (default)
-
get_mllms
(typ, datmaps, use_cls_len=True, use_Pool=0, **kwargs)[source]¶ Returns maximum likelihood sky CMB modes.
This instance uses isotropic filtering
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- datmaps – data real-space maps array
- use_cls_len – use lensed cls in filtering (defaults to True)
Returns: T,E,B alm array
-
get_plmRDlikcurvcls
(typ, datcls_obs, lib_qlm, use_cls_len=True, use_cls_len_D=None, recache=False, dat_only=False)[source]¶ Returns realization-dependent second variation (curvature) of lensing deflection likelihood
Second variation of
\(\frac 12 X^{\rm dat} \rm{Cov_\alpha}^{-1} X^{\rm dat} + \frac 1 2 \ln \det \rm{Cov_\alpha}\)
with respect to deflection for custom data spectra. See https://arxiv.org/abs/1808.10349
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- datcls_obs (dict) – data (beam-convovled) CMB and noise spectra
- lib_qlm – ffs_alm instance describing the lensing alm’s arrays
- use_cls_len – use lensed cls in filtering (defaults to True)
-
get_plmlikcurvcls
(typ, datcmb_cls, lib_qlm, use_cls_len=True, recache=False, dat_only=False)[source]¶ Returns realization-dependent second variation (curvature) of lensing deflection likelihood
Second variation of
\(\frac 12 X^{\rm dat} \rm{Cov_\alpha}^{-1} X^{\rm dat} + \frac 1 2 \ln \det \rm{Cov_\alpha}\)
with respect to deflection for custom data CMB spectra. See https://arxiv.org/abs/1808.10349
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- datcmb_cls (dict) – data CMB spectra (if they do match lensed cls attribute, the result should be inverse lensing Gaussian noise bias)
- lib_qlm – ffs_alm instance describing the lensing alm’s arrays.
- use_cls_len – use lensed cls in filtering (defaults to True)
-
get_qlms
(typ, iblms, lib_qlm, use_cls_len=True, **kwargs)[source]¶ Unormalized quadratic estimates (potential and curl).
Note
the output differs by a factor of two from the standard QE. This is because this was written initially as gradient function of the CMB likelihood w.r.t. real and imaginary parts. So to use this as QE’s, the normalization is 1/2 the standard normalization the inverse response. The lensit.ffs_qlms.qlms.py module contains methods of the QE’s with standard normalizations which you may want to use instead.
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- iblms – inverse-variance filtered CMB alm arrays
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- use_cls_len – use lensed or unlensed cls in QE weights (numerator), defaults to lensed cls
-
get_response
(typ, lib_qlm, cls_weights=None, cls_filt=None, cls_cmb=None, use_cls_len=True)[source]¶ Lensing quadratic estimator gradient and curl response functions.
Parameters: - cls_filt – CMB spectra used in the filtering procedure ( i.e. those entering \(\rm{Cov}^{-1}\)). Defaults to self.cls_len if set else self.cls_unl
- cls_weights (dict) – CMB spectra used in the QE weights (those entering the numerator in the usual Okamoto & Hu formulae) Defaults to self.cls_len if set else self.cls_unl
- cls_cmb (dict) – CMB spectra of the sky entering the response contractions (in principle lensed cls or grad-lensed cls)
-
iterateN0cls
(typ, lib_qlm, itmax, return_delcls=False, _it=0, _cpp=None)[source]¶ Iterative flat-sky \(N^{(0)}_L\) calculation to estimate the noise levels of the iterative estimator.
This uses perturbative approach in Wiener-filtered displacement, consistent with box shape and mode structure.
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- itmax – Number of iterations to performs
- return_delcls – optionally return partially delensed cmb cls as well
-
class
lensit.ffs_covs.ffs_cov.
ffs_lencov_alm
(lib_dir, lib_datalm, lib_skyalm, cls_unl, cls_len, cl_transf, cls_noise, f, fi, init_rank=0, init_barrier=<function <lambda>>)[source]¶ Extended ffs_diagcov_alm sub-class, adding a deflection field and its inverse
Parameters: - lib_dir – many things will be saved there
- lib_datalm – ffs_alm instance containing mode filtering and flat-sky patch info
- lib_skyalm – ffs_alm instance describing the sky modes.
- cls_unl (dict) – unlensed CMB cls
- cls_len (dict) – lensed CMB cls
- cl_transf – instrument transfer function
- cls_noise (dict) – ‘t’, ‘q’ and ‘u’ noise arrays
- f – deflection field, forward operation (e.g. ffs_displacement instance)
- fi – deflection field, backward operation (ideally the inverse f-deflection, f.get_inverse())
-
degrade
(LD_shape, no_lensing=False, ellmax=None, ellmin=None, libtodegrade='sky', lib_dir=None)[source]¶ Degrades covariance matrix to some lower resolution.
-
eval_mf
(typ, mfkey, xlms_sky, xlms_dat, lib_qlm, use_Pool=0, **kwargs)[source]¶ Delfection-induced mean-field estimation from input random phases
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- mfkey – mean-field estimator key
- xlms_sky – i.i.d. sky modes random phases (if relevant)
- xlms_dat – i.i.d. data modes random phases (if relevant)
- lib_qlm – ffs_alm instance describing the lensing alm arrays
Returns: mean-field estimate (gradient and curl component)
-
get_iblms
(typ, datalms, use_cls_len=False, use_Pool=0, **kwargs)[source]¶ Inverse-variance filters input CMB maps
Produces \(B^t \rm{Cov_\alpha}^{-1}X^{\rm dat}\) (inputs to quadratic estimator routines)
The covariance matrix here includes the lensing deflection field \(\alpha\) as given by the self.f and self.fi attributes
-
get_mllms
(typ, datmaps, use_Pool=0, use_cls_len=False, **kwargs)[source]¶ Returns maximum likelihood sky CMB modes.
This instance uses anisotropic filtering, with lensing deflections as defined by self.f and self.fiv
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- datmaps – data real-space maps array
- use_cls_len – use lensed cls in filtering (defaults to True)
Returns: T,E,B alm array
-
get_qlms
(typ, iblms, lib_qlm, use_Pool=0, use_cls_len=False, **kwargs)[source]¶ Unormalized quadratic estimates (potential and curl), including current lensing deflection estimate
Note
the output differs by a factor of two from the standard QE. This is because this was written initially as gradient function of the CMB likelihood w.r.t. real and imaginary parts. So to use this as QE’s, the normalization is 1/2 the standard normalization the inverse response. The lensit.ffs_qlms.qlms.py module contains methods of the QE’s with standard normalizations which you may want to use instead.
Parameters: - typ – ‘T’, ‘QU’, ‘TQU’ for temperature-only, polarization-only or joint analysis
- iblms – inverse-variance filtered CMB alm arrays
- lib_qlm – ffs_alm instance describing the lensing alm arrays
- use_cls_len – use lensed or unlensed cls in QE weights (numerator), defaults to lensed cls