subroutine s2hat_map2alm_spin( pixelization, scan, spin, nlmax, nmmax, nmvals, mvals, nmaps, first_ring, last_ring, local_w8ring, map_size, local_map, lda, local_alm, numprocs, myid, comm)
Purpose:
memory distributed inverse spherical harmonic transform of a field with an arbitrary spin value ( = spin):
calculates a subset of all alm harmonic coefficients for an input field map (or maps) of the sky patch as defined in
scan using the pixelization as defined by pixelization. The output, stored in local_alm,
is distributed over all procs in a way determined by mvals. The input maps are
distributed over all procs as defined by first_ring and last_ring.
Parameter description:
input:
- pixelization-- [pixeltype], an S2HAT structure defining the pixelization to be used;
- scan-- [scandef], an S2HAT structure defining the sky coverage;
- spin-- [integer(4b)], a value of a spin of the field to be synthesized;
- 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 the input maps each made of 2 component maps corresponding
to the real and imaginary parts of the spin field.
All maps have to have the same properties and be distributed over the processors in the same way.
This is also a number of the alm sets computed on the output for each map and each composed of
two sets of Elm and Blm multipoles. All the sets have to have the same properties and be
distributed identically over the processors. On the output they are stored in a single array local_alm.
- first_ring, last_ring-- [integer(4b)],
numbers of a first and last ring (both included) of the Northern hemisphere,
which together with their Southern counterparts define a sub map to be computed by a given proc;
- local_w8ring-- [real(dp)], a pointer to a 2dim array of a size
(1:last_ring-first_ring+1,1:2) 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 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:2,1:nmaps)
which stores all the maps and stokes parameters computed for a given subset of pixels.
If the input spin value is positive then on the input the real part of the complex field should be stored first and the imaginary
one - second. If the spin value is negative, the real part should be stored first while the second part is (-1) x the imaginary part of the field.
That corresponds to a convention that for two real fields local_map[ :,1,:] and local_map[:,2,:] the complex field
of spin spin is given as:
and SPIN_CONV_SIGN is an S2HAT parameter set by default to -1.
- lda-- [integer(4b)], a leading dimension of the local_alm array. Equal either to
= 2 or = nlmax;
- 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:
On the output it stores relevant alm coefficients, both Elm and Blm
for all m values processed on a given proc, all l ([0,nlmax])
and all nmaps which are processed simultanously.
Top of the page
Comments:
a subset of mvalues 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 must be composed of
a complete iso-latitude rings of the northern hemisphere.
The rings of the southern hemisphere are added automatically if needed.
The rings creating the submap are from the interval defined by first_ring and last_ring.
Note that not all rings from the interval have to belong to the map. It is however assumed that the map is partitioned in
a way that the proc 0 gets the first set of rings, proc 1 next etc. No overlaps are allowed.
(These assumptions are left unchanged from Healpix 2.1).
The distribution of rings can be calculated using the routine find_ring_range or
the wrapper routine get_local_data_sizes.
Top of the page
|