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 $\ {_s}Y_{\ell m}$

$\displaystyle {_s}S(p) = \sum_{\ell m} {_s}a_{\ell m}\ \ {_s}Y_{\ell m}(p)$ (12)

for $\ell \ge \vert m\vert, \ell \ge \vert s\vert$, and is such that ${_s}S^* = {_{-s}}S$.
The usual phase convention for the spin weighted harmonics is ${_s}Y_{\ell m}^* = (-1)^{s+m} {_{-s}}Y_{\ell -m}$ and therefore ${_s}a_{\ell m}^* = (-1)^{s+m} {_{-s}}a_{\ell -m}$. The two (real) input maps for map2alm_spin* are defined respectively as

$\displaystyle {_{\vert s\vert}}S^+$ $\displaystyle = ({_{\vert s\vert}}S + {_{-\vert s\vert}}S)/2,$ (13)
$\displaystyle {_{\vert s\vert}}S^-$ $\displaystyle = ({_{\vert s\vert}}S - {_{-\vert s\vert}}S)/(2i).$ (14)

map2alm_spin* outputs the alm coefficients defined as

$\displaystyle {_{\vert s\vert}}a^{+}_{\ell m}$ $\displaystyle = - ( {_{\vert s\vert}}a_{\ell m} + (-1)^s {_{-\vert s\vert}}a_{\ell m} )/2$ (15)
$\displaystyle {_{\vert s\vert}}a^{-}_{\ell m}$ $\displaystyle = - ( {_{\vert s\vert}}a_{\ell m} - (-1)^s {_{-\vert s\vert}}a_{\ell m} )/(2i),$ (16)

for $m\ge 0$, knowing that, just as for spin 0 maps, the coefficients for $m<0$ are given by

$\displaystyle {_{\vert s\vert}}a^{+}_{\ell-m}$ $\displaystyle = (-1)^m {_{\vert s\vert}}a^{+*}_{\ell m},$ (17)
$\displaystyle {_{\vert s\vert}}a^{-}_{\ell-m}$ $\displaystyle = (-1)^m {_{\vert s\vert}}a^{-*}_{\ell m}.$ (18)

With these definitions, ${_2}a^{+}$, ${_2}a^{-}$, ${_2}S^+$ and ${_2}S^-$ match HEALPix polarization $a^E, a^B, Q$ and $U$ respectively. However, for $s=0$, $\ _{0}a^+_{\ell m} = -a^T_{\ell m}$, $\ _{0}a^-_{\ell m} = 0$, $\ {_0}S^+ = T$, $\ {_0}S^- = 0.$

When dealing only with scalar quantities, like temperature or intensity maps, having a spin $s=0$, it is highly recommended, and much more memory-efficient, to use directly the routine map2alm, rather then setting spin$=0$ in map2alm_spin*.

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 $N_{\mathrm{side}}$ value of the map to analyse.
nlmax I4B IN the maximum $\ell$ 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 ${_{\vert s\vert}}S^+$ and ${_{\vert s\vert}}S^-$ input maps
alm(1:2, 0:nlmax, 0:nmmax) SPC/ DPC OUT The ${_{\vert s\vert}}a^+_{\ell m}$ and ${_{\vert s\vert}}a^-_{\ell m}$ output values.
zbounds(1:2), OPTIONAL DP IN section of the map on which to perform the $a_{\ell m}$ analysis, expressed in terms of $z=\sin(\mathrm{latitude}) =
\cos(\theta).$ If zbounds(1)$<$zbounds(2), it is performed on the strip zbounds(1)$<z<$zbounds(2); if not, it is performed outside the strip zbounds(2)$\le z \le$zbounds(1). If absent, the whole map is processed.
w8ring_TQU(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 $N_{\mathrm{side}}$ of 256, and the analysis is performed up to 512 in $\ell$ and $m$. The resulting $a_{\ell m}$ 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_{\ell m}$ recursion
misc_util
module, containing:
assert_alloc
routine to print error message when an array is not properly allocated
Note: Starting with version 3.80, some libsharp routines will be called for any $\vert s\vert$ value.


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

Version 3.82, 2022-07-28