subroutine s2hat_alm2map( precompute_plms, pixelization, scan, nlmax, nmmax, nmvals, mvals, nmaps, nstokes, first_ring, last_ring, map_size, local_map, lda, local_alm, nplm, local_plm, numprocs, myid, comm)
Purpose:
memory distributed direct spherical harmonic transform: synthesizes a map (or maps) of the sky patch as defined in
scan using the pixelization as defined by pixelization. The output, stored in local_map,
is distributed over all procs (as defined by first_ring and last_ring) and calculated
given input spherical harmonic coefficients, stored in local_alm, distributed over all procs (as defined by mvals).
All procs of the communicator comm have to make the call. Blocking.
Parameter description:
input:
- precompute_plms(*)-- [integer(i4b)], = 0,1,2 for no precomputed,
precomputed scalar and
precomputed both, scalar and spin-2, associated Legendre functions, Plm, respectively.
- pixelization(*)-- [pixeltype], an S2HAT structure defining the pixelization to be used;
- scan(*)-- [scandef], an S2HAT structure defining the sky coverage;
- nlmax(*)-- [integer(4b)], maximum considered value of l (polar number) (included);
- nmmax(*)-- [integer(4b)], maximum considered value of m (azimuthal number) (included);
- nmvals-- [integer(4b)], a number of m modes assigned to this proc;
- mvals-- [integer (4b)], a vector of nmvals values of m number to be processed by
a given proc;
- nmaps(*)-- [integer(4b)], a number of input sets of alm coefficients each having nstokes
components. All the sets have to have the same properties and be distributed identically over the processors.
On the input they are stored in a single array local_alm.
This number also determines a number of the output maps (each made of nstokes maps)
which are to be produced. All maps have to have the same properties and be distributed
over the processors in the same way.
- nstokes(*)-- [integer(4b)], a number of Stokes parameters characterizing the maps to be computed:
Similarly, the nstokes value also determines a number of types of the input alm coefficients required for the run.
- first_ring, last_ring-- [integer(4b)], numbers of a first and last ring (both included)
of the Northern hemisphere which together with their Southern hemisphere counterparts define a sub map which will
be computed by a given proc;
- map_size-- [integer(4b)], a number of pixels corresponding to the ring interval defined
by first_ring and last_ring;
- lda(*)-- [integer(4b)], a leading dimension of local_alm array. Equal either to
= nstokes or = nlmax;
- local_alm-- [complex(dp)], an allocated pointer to a 4dim array of dimensions:
storing relevant alm coefficients for all m values processed on a given proc and all l
and all nstokes components for each of the nmaps maps which are to be processed simultanously;
For nstokes == 1 the coefficients are those for total intensity only. For nstokes == 3 the
component corresponds to the total intensity, E-mode and B-mode polarizations, respectively;
- nplm-- [integer(8b)], defines one of the two dimensions of the local_plm array.
It is equal to a product of a number of all rings (Northern hemisphere plus equator) times a number of all l,m
processed by a given proc. The value of this parameter is irrelevant if precompute_plms == 0;
- local_plm-- [real(dp)], an allocated (if needed) pointer to a 2dim array.
It is relevant only if precompute_plms is not zero, and then:
if precommpute_plms == 1, it contains only scalar Plm precomputed for
every ring of the map and subset of m values as defined by
a vector mvals. If precompute_plms == 2, it stores a scalar and two spin-2
harmonics (the conventions and normalizations of which are as defined in Healpix.
The layout of the array depends on the value of lda parameter as for the local_alm array
(and must agree with it),
where
- numprocs(*)-- [integer(4b)], total number of procs used;
- myid-- [integer(4b)], proc id (0,...,numprocs-1);
- comm(*)-- [integer(4b)], an mpi communicator for the numprocs procs;
(*)The value must be the same on all procs of the communicator comm.
output:
- local_map-- [real(dp)], an allocated pointer of a 3dim array of a size (0:map_size-1,1:nstokes,1:nmaps)
which on each processor stores all the maps of nstokes Stokes parameters ( I if nstokes == 1 and
I,Q,U if nstokes == 3) computed for a given range of rings.
Top of the page
Comments:
a subset of m values to be processed by the routine (mvals) is arbitrary,
though can not be empty. The subsets assigned to any two procs must be disjoint. The load balanced distribution of
m modes over procs is calculated by the routines nummvalues and
find_mvalues or, alternately, by the wrapper routine get_local_data_sizes
.
a subset of a map to be computed by any proc will be composed of
a complete iso-latitude rings of the Northern hemisphere, with rings from the interval
defined by first_ring and last_ring as well as of the corresponding rings
of the Southern hemisphere.
All rings are included in full as long as at least are partially observed.
Thus only fully unobserved rings from the interval will be omitted. Consequently, if a final
map consists of partially observed rings, the unobserved pixels of those rings may need to be
marked as such in post processing. (See submap description for
more info.)
The ring intervals can be arbitrary though can not be empty on any processor and can not overlap.
Moreover, it is assumed that
the proc 0 gets the first set of rings, proc 1 next etc.
(These assumptions are left unchanged from Healpix 2.1).
The load balanced distribution of rings can be calculated using the routine find_ring_range or
the wrapper routine get_local_data_sizes.
local_plm can be computed using the routine plm_mvals_gen.
On each processor only modes with m values as defined in mvals are stored.
The values correspond to each Northern hemisphere ring and for all m (as they appear
in the vector mvals) and only for those l modes which have non-zero values (i.e., l >= m ).
the parameter map_size is in fact redundant and can be computed from first_ring
and last_ring using the function local_map_size for any
predefined pixelization.
See also comments to the s2hat_map2alm routine.
Top of the page
|