subroutine s2hat_map2alm( precompute_plms, pixelization, scan, nlmax, nmmax, nmvals, mvals, nmaps, nstokes, first_ring, last_ring, local_w8ring, map_size, local_map, lda, local_alm, nplm, local_plm, numprocs, myid, comm)
Purpose:
memory distributed inverse spherical harmonic transform: calculates a subset of alm coefficients
for the sky patch (or patches) defined in
scan assuming the pixelization as defined by pixelization. The input is stored in
local_map and is assumed to be distributed over all procs (as defined by first_ring
and last_ring). The harmonic coefficients are stored in local_alm and distributed over
all procs (as defined by mvals).
This routine has to be called by all procs of the communicator comm.
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 value of l (polar number) (included);
- nmmax(*)-- [integer(4b)], maximum value of m (azimuthal number) (included);
- nmvals-- [integer(4b)], a number of m modes assigned to a given proc;
- mvals-- [integer (4b)], a vector of nmvals values of m
assigned to a given proc;
- nmaps(*)-- [integer(4b)], a number of input maps (each made of nstokes components corresponding to Stokes parameters) to be processed simultaneously.
All maps need to be identical and be distributed over the processors in the same way. They are all stored in the same array, local_map. Similarly,
on the output nmaps determines a number of calculated sets of the alm coefficients.
- nstokes(*)-- [integer(4b)], a number of Stokes parameters characterizing the maps to be computed.
Simmilarly, the nstokes value also determines a number of types of the output alm coefficients computed for
each of the maps.
- first_ring, last_ring-- [integer(4b)], numbers of a first and last ring (both included) which together with their Southern hemisphere counterparts define
a sub map a transform of which is to be computed. If some of the rings are observed partially the corresponding map values have to be set to zero on the input.
- local_w8ring-- [real(dp)], a pointer to a 2dim array of a size
(1:last_ring-first_ring+1,1:nstokes) containing a subset of the map weights to be applied
to a sub-map to be processed on a given proc;
- map_size-- [integer(4b)], a number of pixels in a map to be computed on a given proc and corresponding to a set
of rings defined by first_ring and last_ring above;
- local_map-- [real(dp)], an allocated pointer of a 3dim array of a size (0:map_size-1,1:nstokes,1:nmaps)
which stores all the maps of all required Stokes parameters for a required subset of pixels (i.e., as defined by first_ring and
last_ring);
- lda(*)-- [integer(4b)], a leading dimension of local_alm array. Equal either to
= nstokes or = nlmax;
- nplm-- [integer(8b)], a total number of precomputed values of
Plm for one of the Stokes params.
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. Its value is irrelevant if precompute_plms == 0;
- local_plm-- [real(dp)], an allocated (if needed) pointer to a 2dim of a layout,
where
Relevant only if precompute_plms is not zero.
If nstokes == 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 nstokes == 3, it stores a scalar and two spin-2 harmonics (the conventions and normalizations of which are as
defined in Healpix.
- 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;
Output:
- local_alm- [complex(dp)], an allocated pointer to a 4dim array of dimensions:
storing computed relevant alm coefficients for all m values processed on a given proc and all l
and of all requested types (I if nstokes == 1; and I, E, B if nstokes == 3) for all nmaps
which are processed simultanously.
(*)The value must be the same on all procs of the communicator comm.
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 processed by any proc has to be composed of
complete iso-latitude rings of the Northern hemisphere
defined by first_ring and last_ring as well as corresponding rings from the South.
If some of the rings are observed only partially the unobserved pixels need to be set to zero on the input.
The fully unobserved rings (as defined in the scan structure) need not to be included.
The 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 a given proc only modes with m values as defined in nmvals are stored.
The Plm values are stored for each ring of the Northern hemisphere
and for all m (as they appear
in the vector mvals) and only those l values which Plm have non-zero values (i.e., l >= m ).
See also comments to the s2hat_alm2map routine.
Top of the page
|