|
map2alm*
This routine is a wrapper to 5 internal routines:map2alm_sc,
map2alm_sc_pre, map2alm_pol, map2alm_pol_pre1,
map2alm_pol_pre2. These routines analyse a HEALPix RING ordered map and return
almT (and if specified almE and almB) values up to
the desired order in l (maximum 3*Nside). The different
routines are called depending on what parameters are passed. Some
routines analyse with or without precomputed harmonics and some with
or without polarisation. Location in HEALPix directory tree:
src/f90/mod/alm_tools.f90
FORMAT call map2alm*(
nsmax, nlmax, nmmax, map_TQU, alm_TGC, zbounds, w8ring_TQU [, plm]
)
ARGUMENTS
name & dimensionality |
kind |
in/out |
description |
|
|
|
|
nsmax |
I4B |
IN |
the Nside value of the map to analyse. |
nlmax |
I4B |
IN |
the maximum l value for the analysis. |
nmmax |
I4B |
IN |
the maximum m value for the analysis. |
map_TQU(0:12*nsmax**2-1) |
SP/ DP |
IN |
if only the temperature map is to be analyse, the map-array should be passed with this rank. |
map_TQU(0:12*nsmax**2-1, 1:3) |
SP/ DP |
IN |
if both temperature an polarisation maps are to be analysed, the map array should have this rank, where the second index is (1,2,3) corresponding to (T,Q,U). |
alm_TGC(1:p, 0:nlmax, 0:nmmax) |
SPC/ DPC |
OUT |
The alm values output from the analysis. p is 1 or 3 dependent on wether polarisation is included or not. In the former case, the first index is (1,2,3) corresponding to (T,E,B). |
zbounds(1:2) |
DP |
IN |
section of the map on which to perform the alm
analysis, expressed in terms of
If zbounds(1)<zbounds(2), the analysis is
performed on the strip zbounds(1)<z<zbounds(2); if not,
it is performed outside of the strip
zbounds(2)<z<zbounds(1). |
w8ring_TQU(1:2*nsmax, 1:p) |
DP |
IN |
ring weights for quadrature corrections. If ring weights are not used, this array should be 1 everywhere. p is 1 for a temperature analysis and 3 for (T,Q,U). |
plm(0:(nlmax+1)(nlmax+2)nsmax-1), OPTIONAL |
DP |
IN |
If this optional matrix is passed with this rank, precomputed
are used instead of recursion. |
plm(0:(nlmax+1)(nlmax+2)nsmax-1,1:3), OPTIONAL |
DP |
IN |
If this optional matrix is passed with this rank, precomputed
AND precomputed tensor harmonics are used instead of recursion. |
EXAMPLE:
use healpix_types
use alm_tools
use pix_tools
integer(i4b) :: nside, lmax
real(dp), allocatable, dimension(:,:) :: dw8
real(dp), dimension(2) :: z
real(sp), allocatable, dimension(:,:) :: map
complex(spc), allocatable, dimension(:,:,:) :: alm
nside = 256
lmax = 512
allocate(dw8(1:2*nside, 1:3))
allocate(map(0:nside2npix(nside)-1,1:3))
allocate(alm(1:3, 0:lmax, 0:lmax)
dw8 = 1.0_dp
z = sin(10.0_dp * DEG2RAD)
call map2alm(nside, lmax, lmax, map, alm, ( z, -z ) , dw8, plm(0:(lmax+1)*(lmax+2)*nside-1))
Analyses temperature and polarisation maps passed in map. The map has
an Nside of 256, and the analysis is performed up
to 512 in l and m. The resulting alm coefficients for
temperature and polarisation are returned in alm. A cut on
each side of the equator is applied. Uniform weights are used. Since
the optional plm array is provided with rank one, precomputed scalar
are
used while tensor harmonics are computed with a recursion.
MODULES & ROUTINESThis section lists the modules and routines used by map2alm*.
-
ring_analysis
- Performs FFT for the ring analysis.
-
misc_util
- module, containing:
-
assert_alloc
- routine to print error message when an array is not
properly allocated
RELATED ROUTINESThis section lists the routines related to map2alm*
-
anafast
- executable using map2alm*to analyse maps.
-
alm2map
- routine performing the inverse transform
of map2alm*.
-
map2alm_iterative
- similar to
map2alm* with iterative scheme.
Eric Hivon
2010-06-18
|
|