plancklens.n0s

This module provides convenience functions to calculate curved-sky responses and reconstruction noise curve for lensing or other estimators In plancklens a QE is often described by a short string. For example ‘ptt’ stands for lensing (or lensing gradient mode) from temperature x temperature. Anisotropy source keys are a one-letter string including

‘p’ (lensing gradient) ‘x’ (lensing curl) ‘s’ (point sources) ‘f’ (modulation field) ‘a’ (polarization rotation)

Typical keys include then:

‘ptt’, ‘xtt’, ‘stt’, ‘ftt’ for the corresponding QEs from temperature only ‘p_p’, ‘x_p’, ‘f_p’, ‘a_p’ for the corresponding QEs from polarization only (combining EE EB and BB if relevant) ‘p’, ‘x’, ‘f’, ‘a’, ‘f’ … for the MV (or GMV) combination ‘p_eb’, … for the EB estimator (this is the symmetrized version (‘peb’ + ‘pbe’) / 2 so that E and B appear each once on the gradient and inverse-variance filtered leg)

Bias-hardening can be included by inserting ‘_bh_’.

E.g. ‘ptt_bh_s’ is the lensing TT QE bias-hardened against point source contamination using the ‘stt’ estimator

Responses method takes as input the QE weights (typically the lensed CMB spectra) and the filtering cls (‘fals’) which describes the filtering applied to the maps (the \((C + N)^{-1}\) operation) get_N0_iter calculates an estimate of the N0s for iterative lensing estimator beyond the QE

plancklens.n0s.get_N0(beam_fwhm=1.4, nlev_t: float = 5.0, nlev_p: array | None = None, lmax_CMB: dict = 3000, lmin_CMB=100, lmax_out=None, cls_filt: dict | None = None, cls_len: dict | None = None, cls_weight: dict | None = None, cls_sky: dict | None = None, joint_TP=True, ksource='p', wfleg_Tcut=None)[source]

Example function to calculates reconstruction noise levels for a bunch of quadratic estimators

Parameters:
  • beam_fwhm – beam fwhm in arcmin, assuming gaussian shape.

  • nlev_t – T white noise level in uK-arcmin (an array of size lmax_CMB can be passed for scale-dependent noise level)

  • nlev_p – P white noise level in uK-arcmin (defaults to root(2) nlevt, if None) (can also be a (non-white) array of shape (lmax_CMB), or with either polarization or E and B noise separately, in which case the array is of shape (2,lmax_CMB)), and E noise expected at first

  • lmax_CMB – max. CMB multipole used in the QE (use a dict with ‘t’ ‘e’ ‘b’ keys instead of int to set different CMB lmaxes)

  • lmin_CMB – min. CMB multipole used in the QE

  • lmax_out – max lensing ‘L’ multipole calculated

  • cls_filt – set of spectra used for the filtering, defaults to FFP10 lensed CMB spectra

  • cls_len – CMB spectra entering the sky response to the anisotropy (defaults to FFP10 lensed CMB spectra, in general should be the lensed gradient spectra)

  • cls_sky – actual spectra of the measured CMB (without noise); defaults to FFP10 lensed CMB spectra

  • cls_weight – CMB spectra entering the QE weights (defaults to FFP10 lensed CMB spectra; optimally, gradient spectra)

  • joint_TP – if True include calculation of the N0s for the GMV estimator (incl. joint T and P filtering)

  • ksource – anisotropy source to consider (defaults to ‘p’, lensing)

  • wfleg_Tcut – high-l cut on the T gradient leg if set

Returns:

N0s array for the lensing gradient and curl modes for the T-only, P-onl and (G)MV estimators

Prompted by AL

plancklens.n0s.get_N0_iter(qe_key: str, nlev_t: float, nlev_p: float, beam_fwhm: float, cls_unl_fid: dict, lmin_cmb: int, lmax_cmb: int, itermax, cls_unl_dat=None, lmax_qlm=None, ret_delcls=False, datnoise_cls: dict | None = None, ret_curl=False, rho_sqd_ext: float = 0, filter_E=False)[source]

Iterative lensing-N0 estimate

This calculates iteratively partially lensed spectra and lensing noise levels. At each iteration this takes out the resolved part of the lenses (ignoring N1) and recomputes a N0

Parameters:
  • qe_key – QE estimator key

  • nlev_t – temperature noise level (in :math:`mu `K-arcmin) (an array can be passed for scale-dependent noise level)

  • nlev_p – polarisation noise level (in :math:`mu `K-arcmin)(an array can be passed for scale-dependent noise level) can also be for E and B noise separately, in which case the array is of shape (2,lmax_CMB)), and E noise expected at first

  • beam_fwhm – Gaussian beam full width half maximum in arcmin

  • cls_unl_fid (dict) – unlensed CMB power spectra

  • lmin_cmb – minimal CMB multipole used in the QE

  • lmax_cmb – maximal CMB multipole used in the QE

  • itermax – number of iterations to perform

  • lmax_qlm (optional) – maximum lensing multipole to consider. Defaults to 2 lmax_ivf

  • ret_delcls (optional) – returns the partially delensed CMB cls as well if set

  • datnoise_cls (optional) – feeds in custom noise spectra to the data. The nlevs and beam only apply to the filtering in this case

  • ret_curl (optional) – also returns lensing curl N0s

  • rho_sqd_ext (optional scalar or array) – option cross-correlation coefficent squared of external tracer to use for partial delensing

  • filter_E (optional) – do linear delensing by substractng B estimate constructed by partial lensing of Wiener-filtered unlensed E

Returns

tuple of arrays of shape (itermax + 1, lmax_qlm + 1) with all iterated N0s. First entry is standard N0. See code below.

Note

this requires camb python package for the lensed spectra calc.