alm2map_der*

This routine is a wrapper to four other routines that synthesize a HEALPix temperature (and polarisation) map(s), its (their) first derivatives, and optionally its (their) second derivatives. The routines accept both single and double precision arrays for alm, map, der1 and der2. The precision of these arrays should match. All maps produced are RING ordered.

See ”Fortran Facilities” Appendix for a note on a bug affecting the calculation of polarisation derivatives on past versions of this routine.

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


FORMAT

call alm2map_der*( nsmax, nlmax, nmmax, alm, map, der1[, der2=, zbounds=] )


ARGUMENTS

name & dimensionality kind in/out description
       
nsmax I4B IN the $N_{\mathrm{side}}$ value of the map to synthesize.
nlmax I4B IN the maximum $\ell$ value used for the $a_{\ell m}$.
nmmax I4B IN the maximum $m$ value used for the $a_{\ell m}$.
alm(1:p, 0:nlmax, 0:nmmax) SPC/ DPC IN The $a_{\ell m}$ values to make the map from. p is either 1 (temperature only) or 3 (temperature+polarisation).
map(0:12*nsmax**2-1) or (0:12*nsmax**2-1,1:3) SP/ DP OUT temperature map $T(p)$ or temperature + polarisation maps $T(p)$, $Q(p)$, $U(p)$ to be synthesized.
der1(0:12*nsmax**2-1, 1:2*p) SP/ DP OUT contains on output the first derivatives of T: $\left({\partial T}/{\partial \theta}, {\partial T}/{\partial \phi}/\sin\theta \right)
$ or the interleaved derivatives of T, Q, and U: $\left({\partial T}/{\partial
\theta}, {\partial Q}/{\partial \theta}, {\partial U}/{\partial \theta};\right. $ $\left.{\partial T}/{\partial \phi}/\sin\theta, \ldots \right)
$
der2(0:12*nsmax**2-1,1:3*p), OPTIONAL SP/ DP OUT If this optional matrix is passed with this rank, it will contain on output the second derivatives $\left({\partial^2 T}/{\partial \theta^2}, {\partial^2 T}/{\partial
\theta\partial\phi}/\sin\theta,\right. $ $\left.{\partial^2 T}/{\partial \phi^2}/\sin^2\theta \right) $ or $\left({\partial^2 T}/{\partial \theta^2}, {\partial^2 Q}/{\partial \theta^2},
{\partial^2 Q}/{\partial \theta^2}, \ldots \right) $
zbounds(1:2), OPTIONAL DP IN section of the sphere on which to perform the map synthesis, 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.


EXAMPLE:

use healpix_types
use pix_tools, only : nside2npix
use alm_tools, only : alm2map_der
integer(I4B) :: nside, lmax, mmax, npix
real(SP), dimension(:), allocatable :: map
real(SP), dimension(:,:), allocatable :: der1, der2
complex(SPC), dimension(:,:,:), allocatable :: alm
...
nside=256 ; lmax=512 ; mmax=lmax
npix=nside2npix(nside)
allocate(alm(1:1,0:lmax,0:mmax))
allocate(map(0:npix-1))
allocate(der1(0:npix-1,1:2), der2(0:npix-1,1:3))
...
call alm2map_der(nside, lmax, mmax, alm, map, der1, der2=der2)
Make temperature maps and its derivatives from the $a_{\ell m}$ passed in alm. The maps have $N_{\mathrm{side}}$ of 256, and are constructed from $a_{\ell m}$ values up to 512 in $\ell$ and $m$.


MODULES & ROUTINES

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

ring_synthesis
Performs FFT over $m$ for synthesis of the rings.
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_utils
module, containing:
assert_alloc
routine to print error message, when an array can not be allocated properly


RELATED ROUTINES

This section lists the routines related to alm2map_der*

alm2map
routine generating maps of temperature and polarisation from their $a_{\ell m}$
alm2map_spin
routine generating maps of arbitrary spin from their ${_s}a_{\ell m}$
synfast
executable using alm2map_der* to synthesize maps.
create_alm
routine to generate randomly distributed $a_{\ell m}$ coefficients according to a given power spectrum

Version 3.83, 2024-11-13