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