NASA - Jet Propulsion Laboratory
    + View the NASA Portal
Search JPL
Jet Propulsion Laboratory Home Earth Solar System Stars & Galaxies Technology
Introduction Background Software Links


map2alm_spin*


This routine extracts the alm coefficients out of maps of spin s and -s. A (complex) map S of spin s is a linear combination of the spin weighted harmonics sYlm
\begin{displaymath}
{_s}S(p) = \sum_{lm} {_s}a_{lm}\ \ {_s}Y_{lm}(p)
\end{displaymath} (8)

for $l \ge \vert m\vert, l \ge \vert s\vert$, and is such that sS* = -sS.
The two (real) input maps for map2alm_spin* are defined respectively as
|s|S+ $\textstyle \myequal$ (|s|S + -|s|S)/2 (9)
|s|S- $\textstyle \myequal$ (|s|S - -|s|S)/(2i). (10)

map2alm_spin* outputs the alm coefficients defined as
htmlcommentmark>362 |s|a+lm $\textstyle \myequal$ - ( |s|alm + (-1)s -|s|alm )/2 (11)
|s|a-lm $\textstyle \myequal$ - ( |s|alm - (-1)s -|s|alm )/(2i) (12)

for $m\ge 0$, knowing that, just as for spin 0 maps, the coefficients for m<0 are given by
|s|a+l-m $\textstyle \myequal$ (-1)m |s|a+*lm, (13)
|s|a-l-m $\textstyle \myequal$ (-1)m |s|a-*lm. (14)

With these definitions, 2a+, 2a-, 2S+ and 2S- match HEALPix polarization aE, aB, Q and U respectively. However, for s=0, $\ _{0}a^+_{lm} = -a^T_{lm}$, $\ _{0}a^-_{lm} = 0$, $\ {_0}S^+ = T$, $\ {_0}S^- = 0.$


Location in HEALPix directory tree: src/f90/mod/alm_tools.f90


FORMAT

call map2alm_spin*( nsmax, nlmax, nmmax, spin, map, alm [, zbounds, w8ring_TQU] )


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.
spin I4B IN the spin s of the maps to be analysed (only its absolute value is relevant).
map(0:12*nsmax**2-1, 1:2) SP/ DP IN |s|S+ and |s|S- input maps
alm(1:2, 0:nlmax, 0:nmmax) SPC/ DPC OUT The |s|a+lm and |s|a-lm output values.
zbounds(1:2), OPTIONAL DP IN section of the map on which to perform the alm analysis, expressed in terms of $z=\sin({\rm latitude}) =
\cos(\theta).$ 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(1:2*nsmax,1:2), OPTIONAL DP IN ring weights for quadrature corrections. If ring weights are not used, this array should be 1 everywhere.


EXAMPLE:

use healpix_types
use alm_tools
use pix_tools
integer(i4b) :: nside, lmax, spin
real(sp), allocatable, dimension(:,:) :: map
complex(spc), allocatable, dimension(:,:,:) :: alm

nside = 256
lmax = 512
spin = 5
allocate(map(0:nside2npix(nside)-1,1:2))
allocate(alm(1:2, 0:lmax, 0:lmax)
...
call map2alm_spin(nside, lmax, lmax, spin, map, alm)

Analyses spin 5 and -5 maps. The maps have an Nside of 256, and the analysis is performed up to 512 in l and m. The resulting alm coefficients for are returned in alm.


MODULES & ROUTINES

This section lists the modules and routines used by map2alm_spin*.

ring_analysis
Performs FFT for the ring analysis.
compute_lam_mm, get_pixel_layout,
gen_lamfac_der, gen_mfac,
gen_recfac, init_rescale, l_min_ylm
Ancillary routines used for $\ {_s}Y_{lm}$ recursion
misc_util
module, containing:
assert_alloc
routine to print error message when an array is not properly allocated


RELATED ROUTINES

This section lists the routines related to map2alm_spin*


alm2map_spin
routine performing the inverse transform of map2alm_spin*.
map2alm
routine analyzing temperature and polarization maps

Eric Hivon 2010-06-18
Privacy / Copyright
FIRST GOV Contact: NASA Home Page Site Manager:
Webmaster:

CL 03-2650