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
get_Hessian(it, key)[source]

Build the L-BFGS Hessian at iteration it

get_Plm(it, key)[source]

Loads solution at iteration it

get_gradPpri(it, key, cache_only=False)[source]

Builds prior gradient at iteration it

how_many_iter_done(key)[source]

Returns the number of points already calculated. 0th is the qest.

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

load_gradquad(it, key)[source]

Loads likelihood quadratic piece gradient at iteration it

Gradient must have already been calculated

load_total_grad(it, key)[source]

Load the total gradient at iteration it.

All gradients 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.
calc_gradplikpdet(it, key)[source]

Calculates the likelihood gradient (quadratic and mean-field parts)

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.
calc_gradplikpdet(it, key)[source]

Calculates the likelihood gradient (quadratic and mean-field parts)

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.
build_pha(it)[source]

Builds sims for the mean-field evaluation at iter it

calc_gradplikpdet(it, key, callback='default_callback')[source]

Caches the det term for iter via MC sims, together with the data one, with MPI maximal //isation.