Up: Preliminaries    Previous: Parameters/structures


S2HAT data objects

  • pointers and arrays;
  • maps;
  • submaps;
  • harmonic coefficients;
  • full-size maps and alm arrays;
  • associated legendre functions, plm.

    pointers and arrays:
    the library routines expect on the input, and provide on the output, pointers to the multidimensional arrays. As such they need to follow the indexing schemes as described below and their indices need to conform with the assumed index limits. That can be ensured either by directly allocating the pointer with the required index bounds, i.e.,

    real(DP), dimension(:,:,:), pointer :: pmap

    allocate(pmap(0:local_npix-1, 1:ncomp, 1:nmaps))

    or allocate an array predefined with a target keyword and associate a pointer with it, i.e.,

    real(DP), dimension(:,:,:), allocatable, target :: map
    real(DP), dimension(:,:,:), pointer :: pmap

    allocate(map(0:local_npix-1, 1:ncomp, 1:nmaps))

    pmap => map.


    Refer to a Fortran 90 manual for more details.

    Top of the page

    maps -- a pointer to a double precision three-dimensional array: map[0:npix-1, 1:ncomp, 1:nmaps],
        where npix is a number of pixels of a map (or a submap as assigned to a given processor);
    ncomp is a number of map components (e.g., Stokes parameters) assuming values of 1, 2, or 3 corresponding to:

    1 - a total intensity map;

    2 - two real maps corresponding to a single complex field of (non-zero, but arbitrary) spin s; the convention is as that: map[ 0:npix-1, 1, 1:nmaps] + i (-SPIN_CONV_SIGN) map[ 0:npix-1, 2, 1:nmaps] is a complex field with spin s;
    (where SPIN_CONV_SIGN is fixed by default to -1 (can be changed, see s2hat_types.h and s2hat_types.f90)).

        NB. for a single real field spin s = 0 it may be more efficient to use the scalar (total intensity) transforms. However, not the difference in the sign conventions.

    3 - three (real) maps for each of three Stokes parameters (I,Q,U).

    nmaps -- a number of maps to be processed simultanously.


    • The values in the map are always assumed to be ordered ring by ring from one pole (usually identified as North Pole) to the other. The pixel numbering scheme follows the same assumptions (i.e., corresponds to the Healpix ring ordering).
    • In the language used hereafter, the object map[ 0:npix-1, 1:ncomp, 1:nmaps] contains nmaps maps (or submaps) composed of ncomp components (e.g., maps of different Stokes parameters), each containing npix pixels.
    • Note that the full-size, non-distributed maps are stored only as 2dim arrays with the 3rd index (nmaps) suppressed. (With an exception of a "do-it-all" driver distribute_local_data_objects_map2alm.) See full-size arrays for more info.


    Top of the page

    submaps ("local_map") -- all maps considered here are divided into submaps which are distributed over the processors. The submap arrays are of the same form as all the maps considered here as described above. Each submap consists of a number of complete iso-latitude rings and includes the corresponding rings (i.e., symmetric with respect to the equator) of the Northern and Southern hemispheres. Moreover, for each hemisphere the submap is composed of the consecutive rings, however with the completely unobserved rings excised. Consequently, the submaps are defined by two integer numbers standing for a first and last ring of the Northern hemisphere (plus the equator, if present) denoted hereafter as first_ring and last_ring. Note that both of the limits are always assumed to belong to the submap.

    For cut sky experiments, the actual sky considered in calculation of the transforms is in fact composed of all rings crossing the observed patch, i.e., includes all the rings which have at least one observed pixel. The inverse transforms (map -> alm) thus utilizes the map values in all the pixels of all observed rings (including those in fact unobserved). Hence, the map values for the unobserved pixels need to be set explicitly to zero prior to calling the routines. For the direct transforms the extra unobserved values will be also computed and thus need to be rejected (or marked as unobserved) in post-processing.

    Note that the extra computational cost due to the patch expansion in the azimuthal direction is in fact negligible. The memory overhead can be however occasionally significant (e.g., for the narrow sky patches elongated along the meridian) that can be overcome by simple distributing the run over more processors and/or rotating the pixelization coordinate system. There seems not to be any general way around, performing well for the nearly full sky and small sky coverages. The solution adopted here seems to be the "lesser evil".

    Top of the page


    harmonic coefficients ("alm") -- a pointer to a double precision complex four-dimensional array. Two different layouts are supported to provide some level of compatibility with the serial Healpix routines. (Healpix alm arrays are 3-dimensional as they do not allow for a multiple map option. How to combine Healpix serial and S2HAT libraries is suggested here. They are:

    S2HAT native layout: alm[ 0:lmax, 0:mvals-1, 1:ncomp, 1:nmaps]
    HEALPIX-like layout: alm[ 1:ncomp, 0:lmax, 0:mvals-1, 1:nmaps]

    where

    lmax -- a maximal l value (harmonic polar index) to be recovered;
    nmvals -- a number of different m (harmonic azimuthal index) values as assigned to a given proc. The actual m associated with a given proc are stored in a vector mvals;
    ncomp -- a number of components;
    nmaps -- a number of the "alm" sets (i.e., maps) to be processed simultaneously.


    The choice of the input matrix layout is decided via an extra input parameter, lda (4-byte integer), defining the leading dimension of the array, which is either equal to lmax -- the S2HAT convention or ncomp for the HEALPIX convention.

    • In the language used hereafter the array, alm defined here contains nmaps sets of alm coefficients each including ncomp components (e.g., corresponding to I, E, B ).
    • Note that full-size, non-distributed alm array appear explicitly only as 3 dimensional array with the fourth index (nmaps) suppressed. (With an exception of a "do-it-all" driver distribute_local_data_objects_alm2map.) See full-size arrays for more info.


    Top of the page


    full-size maps and alm arrays -- in the S2HAT library full-size, non-distributed alm arrays and full-size non-distributed maps appear explicitly only as 3 and 2 dimensional arrays with the last index (nmaps) suppressed in comparison with their distributed versions (*). This is dictated by the assumption that storing multiple complete sets of the alm coefficients in a memory of a single processor will be generally prohibitively expensive.
    Note the full-size arrays of alm or maps are only used as by some of the gather/scatter routines and provided, for example, to facilitate interface with the I/O part of a user code. The underlying assumption of the S2HAT library is that all operations can be performed on the memory-distributed objects bypassing any memory bottlenecks.

    The 3dim full-size arrays provide also a level of compatibility with serial HEALpix routines. (See Compatibility with Healpix).

    (*)With an exception of two "do-it-all" drivers distribute_local_data_objects_map2alm and distribute_local_data_objects_alm2map which are provided only for completeness.

    Top of the page


    associated Legendre functions -- a pointer to a double precision two-dimensional array. Two different layouts are supported, corresponding to the S2HAT native (plm[0:nplm-1,1:ncomp]) and HEALPIX (plm[1:ncomp,0:nplm-1]) choices.

    Here
    nplm -- a number of (l,m) modes assigned to a given processor multiplied by a number of constant elevation rings as also assigned to it.
    (This number can be estimated by the get_local_data_sizes or, more directly, by the nummmodes routines.)
    ncomp -- a number of map components in this case meaning specifically Stokes parameters and thus equal either to 1 or 3.


    For every component the values are assigned to the plm array as follows:

      for every ring (from 1 to nringsall)
        for every m ( m = mvals[0],...,mvals[nmvals-1])
          for every l ( l = m, ..., nlmax)


        NB. The layout choice has to be consistent for all objects of any specific routine.

        Top of the page


  • Up: Preliminaries    Previous: Parameters/structures

    radek stompor 2007-09-15