map2alm_iterative
![]() |
(10) |
a(0) = A.w.m. | (11) |
call map2alm_iterative( nsmax, nlmax, nmmax, iter_order, map_TQU, alm_TGC[, zbounds, w8ring_TQU, plm, mask] )
name & dimensionality | kind | in/out | description |
---|---|---|---|
nsmax | I4B | IN | the Nside value of the map to analyse. |
nlmax | I4B | IN | the maximum ![]() ![]() |
nmmax | I4B | IN | the maximum m value for the analysis. |
iter_order | I4B | IN | the order of Jacobi iteration. Increasing that order
improves the accuracy of the final ![]() ![]() |
map_TQU(0:12*nsmax**2-1, 1:p) | SP/ DP | INOUT | input map. p is 1 or 3 depending if temperature (T) only or temperature and polarisation (T, Q, U) are to be analysed. It will be altered on output if a mask is provided and/or if iter_order>0 and zbounds is provided. |
alm_TGC(1:p, 0:nlmax, 0:nmmax) | SPC/ DPC | OUT | The ![]() |
zbounds(1:2), OPTIONAL | DP | IN | section of the map on which to perform the ![]() ![]() ![]() |
Analyses temperature and polarisation signals in the first 10000 pixels of map (as determined by mask). The map has an Nside of 256, and the analysis is supposed to be performed up to 512 inand m. The resulting
coefficients for temperature and polarisation are returned in alm. Uniform weights are assumed. In order to improve the
accuracy, 2 Jacobi iterations are performed.
Version 3.50, 2018-12-10