alm2map*
This routine is a wrapper to 10 other routines: alm2map_sc_X,
alm2map_sc_pre_X, alm2map_pol_X, alm2map_pol_pre1_X,
alm2map_pol_pre2_X, where X stands for either s or d. These routines
synthesize a HEALPix RING ordered temperature map (and if specified, polarisation maps)
from input
(and if specified
and
) values.
The different routines are called dependent on what parameters are passed.
Some routines synthesize maps with or without precomputed harmonics (note that
since HEALPix v2.20 precomputed harmonics most likely won't speed up computation)
and some with or without polarisation.
The routines accept both single and double precision arrays for alm_TGC and
map_TQU. The precision of these arrays should match.
Location in HEALPix directory tree: src/f90/mod/alm_tools.F90
FORMAT
ARGUMENTS
name & dimensionality |
kind |
in/out |
description |
|
|
|
|
nsmax |
I4B |
IN |
the
value of the map to synthesize. |
nlmax |
I4B |
IN |
the maximum value used for the
. |
nmmax |
I4B |
IN |
the maximum value used for the
. |
alm_TGC(1:p, 0:nlmax, 0:nmmax) |
SPC or DPC |
IN |
The
values to make
the map from. is 3 or 1 depending on wether polarisation is
respectively included or not. In the former case, the first
index runs from 1 to 3 corresponding to (T,E,B). |
map_TQU(0:12*nsmax**2-1) |
SP or DP |
OUT |
if only a temperature map is
to be synthesized, the map-array should be passed with this rank. |
map_TQU(0:12*nsmax**2-1, 1:3) |
SP or DP |
OUT |
if both temperature an
polarisation maps are to be synthesized, the map array should have this rank,
where the second index is (1,2,3) corresponding to (T,Q,U). |
plm(0:n_plm-1), OPTIONAL |
DP |
IN |
If this optional matrix is passed with
this rank, precomputed
are used instead of recursion. (
n_plm = nsmax*(nmmax+1)*(2*nlmax-nmmax+2). |
plm(0:n_plm-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. (n_plm =
nsmax*(nmmax+1)*(2*nlmax-nmmax+2). |
zbounds(1:2), OPTIONAL |
DP |
IN |
section of the sphere on which to perform the map synthesis, expressed in terms of
If zbounds(1)zbounds(2), it is
performed on the strip zbounds(1)zbounds(2); if not,
it is performed outside the strip
zbounds(2)zbounds(1). If absent, the whole map is processed.
Currently, zbounds and plm can not be used together. |
EXAMPLE:
use healpix_types
use pix_tools, only : nside2npix
use alm_tools, only : alm2map
integer(I4B) :: nside, lmax, mmax, npix
real(SP), dimension(:,:), allocatable :: map
complex(SPC), dimension(:,:,:), allocatable :: alm
real(DP), dimension(1:2) :: zrange
...
nside=256 ; lmax=512 ; mmax=lmax
npix=nside2npix(nside)
allocate(alm(1:3,0:lmax,0:mmax))
allocate(map(0:npix-1,1:3))
...
zrange =(/ 0.0_dp, 0.5_dp /)
call alm2map(nside, lmax, mmax, alm, map, zbounds=zrange)
Make temperature and polarisation maps from the scalar and tensor
passed in alm. The maps have
of 256, and are constructed from
values up to 512 in and . In order to save time,
the maps are only generated on the range
(leaving the other pixels to 0)
even though the input
are those of a full sky map.
MODULES & ROUTINES
This section lists the modules and routines used by alm2map*.
-
ring_synthesis
- Performs FFT over for synthesis of the rings.
-
compute_lam_mm, get_pixel_layout,
-
gen_lamfac,gen_mfac, gen_normpol,
-
gen_recfac, init_rescale, l_min_ylm
- Ancillary routines used
for
recursion
-
misc_utils
- module, containing:
-
assert_alloc
- routine to print error message, when an array can not be
allocated properly
Note: Starting with version 3.10, libsharp routines will be called when precomputed
are not provided.
RELATED ROUTINES
This section lists the routines related to alm2map*
-
alm2map_der
- routine generating a map and
its derivatives from its
-
alm2map_spin
- routine generating maps of
arbitrary spin from their
-
smoothing
- executable using alm2map* to smooth maps
-
synfast
- executable using alm2map* to synthesize maps.
-
map2alm
- routine performing the inverse transform
of alm2map*.
-
create_alm
- routine to generate randomly
distributed
coefficients according to a given power spectrum
-
pixel_window,
generate_beam
- return the -space HEALPix -pixel and beam window function respectively
-
alter_alm
- modifies
to emulate effect
of real space filtering
Version 3.83, 2024-11-13