Next: distribute_local_data_objects_alm2map    Up: Gather/scatter routines


distribute_local_data_objects_map2alm

subroutine distribute_local_data_objects_map2alm( int precompute_plms, pixeltype pixelization, scandef scan, int nlmax, int nmmax, int nmaps, int ncomp, double *map, int map_size, double *local_map, int nmvals, int *mvals, int lda, long int nplm, double *local_plm, int first_ring, int last_ring, double *local_w8ring, double *w8ring, int myid, int numprocs, int root, MPI_Comm comm)

Purpose:

A general "do-it-all" driver which distributes the map and ring weights initially stored on a proc root, precomputes a distribution of the m-modes of processor and, on each proc, computes a relevant set of Plm if requested. To be used together and after a call to get_local_data_size for a workload balanced distribution and prior to a call to one of the map2alm transform rouines.
This routine is to be called by all procs of the communicator comm. Blocking.

Parameter description:

input:
precompute_plms(*)-- [int], (= 0,1,2) for no precomputation precomputation of only scalar, and of both scalar and spin associated Legendre functions in a distributed form, respectively, if ncomp == 2 (i.e., spin transforms are to be performed) then precompute_plms is forced to be 0;

pixelization(*)-- [pixeltype], an S2HAT structure containing the info about the pixelization;

scan(*)-- [ scandef], an S2HAT structure containing the info about the sky coverage;

nlmax(*)-- [int], maximum l value considered (included);

nmmax(*)-- [int], maximum m value considered (included);

nmaps(*)-- [int], a total number of maps (each composed of ncomp component maps) to be distributed over all processors;

ncomp(*)-- [int], a number of components of every map to be distributed (=1,2,3);

map_size-- [int], a number of pixels in a assumed pixelization corresponding to a set of rings defined by first_ring and last_ring as assigned to each proc;

map-- [double], a vector of doubles storing a 3dim array of dimensions (0:npixall-1, 1:ncomp, 1:nmaps) in the column-wise order. It contains all the maps on the proc root. This needs to be allocated only on the proc root;

first_ring, last_ring -- [int], numbers of a first and last ring which define a sub map (both included) as assigned to each proc;

w8ring-- [double], a vector of doubles corresponding to a 2-dim array of a size (1:nringsall,1:ncomp) stored in the column-wise order. The vector stores on a proc root all the ring weights. It needs to be allocated only on the proc root;

nmvals-- [int], a number of m values stored on a given proc;

lda(*)-- [int], a leading dimension of local_plm array. Equal either to = ncomp or = nlmax. Meaningfull only for precompute_plms /= 0;

nplm-- [long int (8bytes)], a total number of precomputed Plm of one kind, (irrelevant if precompute_plms == 0);

myid-- [int], proc rank (0,...,numprocs-1);

numprocs(*)-- [int], total number of procs used;

root(*)-- [int], defines a proc root on which the non-distributed data are initially stored;

comm(*)-- [MPI_Comm], an mpi communicator for the numprocs procs.
(*)The value must be the same on all procs of the communicator comm.

output:

mvals-- [int], a vector of a size nmvals containing m values to be processed by a given proc;

local_map-- [double], a vector of doubles corresponding to a 3-dim array of a size (0:map_size-1,1:ncomp,1:nmaps) stored in the column-wise order. On each proc it stores its part of the input map map initially stored only on the proc root;

local_w8ring-- [double], a vector of doubles storing a 2-dim array of a size (1:last_ring-first_ring+1,1:ncomp) in the column-wise order. It contains a subset of the ring weights corresponding to the sub-map stored in the memory of a given proc.

local_plm-- [double], a vector of doubles corresponding to a 2dim array stored in the column-wise order. The layout of the 2dim array depends on the value of lda parameter as for the local_alm array (and must agree with it),

(1:nprec,0:nplm-1), if lda == ncomp;      (HEALpix convention);
(0:nplm-1,1:nprec), if lda == nlmax;      (S2HAT convention);
where
nprec = 1 if precompute_plms == 1, i.e., only scalar Plm are precomputed,
nprec = 3 if precompute_plms == 2, i.e., scalar and spin-2 Plm are precomputed.

It is relevant only if precompute_plms is not zero, and then: if precompute_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.
Comments:

  • memory inefficient as on the root a full map and weights need to be stored at the same time as all distributed objects.
  • all the output (distributed) objects need to be allocated beforehand.
  • This routine exceptionally requires a full map to be defined as a 3 dimensional object storing all the info for all components and maps.


  • Next: distribute_local_data_objects_alm2map Up: Gather/scatter routines

    radek stompor 2009-10-15