mirror of
https://bitbucket.org/cosmicvoids/vide_public.git
synced 2025-07-05 07:41:11 +00:00
103 lines
4.5 KiB
Text
103 lines
4.5 KiB
Text
/*! \mainpage libpsht documentation
|
|
<ul>
|
|
<li>\ref introduction "Introduction"
|
|
<li><a href="modules.html">Programming interface</a>
|
|
</ul>
|
|
*/
|
|
|
|
/*! \page introduction Introduction to libpsht
|
|
|
|
"PSHT" is an acronym for <i>Performant Spherical Harmonic Transforms</i>.
|
|
All user-visible data types and functions in this library start with
|
|
the prefix "psht_", or with "pshts_" and "pshtd_" for single- and
|
|
double precision variants, respectively.
|
|
|
|
<i>libpsht</i>'s main functionality is the conversion between <i>maps</i>
|
|
on the sphere and <i>spherical harmonic coefficients</i> (or <i>a_lm</i>).
|
|
A map is defined as a set of <i>rings</i>, which in turn consist of
|
|
individual pixels that
|
|
<ul>
|
|
<li>all have the same colatitude and</li>
|
|
<li>are uniformly spaced in azimuthal direction.</li>
|
|
</ul>
|
|
Consequently, a ring is completely defined by
|
|
<ul>
|
|
<li>its colatitute (in radians)</li>
|
|
<li>the number of pixels it contains</li>
|
|
<li>the azimuth (in radians) of the first pixel in the ring</li>
|
|
<li>the weight that must be multiplied to every pixel during a map
|
|
analysis (typically the solid angle of a pixel in the ring) </li>
|
|
<li>the offset of the first ring pixel in the <i>map array</i></li>
|
|
<li>the stride between consecutive pixels in the ring.</li>
|
|
</ul>
|
|
The map array is a one-dimensional array of type <i>float</i> or
|
|
<i>double</i>, which contains the values of all map pixels. It is assumed
|
|
that the pixels of every ring are stored inside this array in order of
|
|
increasing azimuth and with the specified stride. Note however that the rings
|
|
themselves can be stored in any order inside the array.
|
|
|
|
The a_lm array is a one-dimensional array of type <i>pshts_cmplx</i> or
|
|
<i>pshtd_cmplx</i>, which contains all spherical harmonic coefficients
|
|
for 0<=m<=mmax and m<=l<=lmax. There is only one constraint on the
|
|
internal structure of the array, which is:
|
|
|
|
<code>Index[a_l+1,m] = Index[a_l,m] + stride</code>
|
|
|
|
That means that coefficients with identical <i>m</i> but different <i>l</i>
|
|
can be interpreted as a one-dimensional array in <i>l</i> with a unique
|
|
stride.
|
|
|
|
Several functions are provided for efficient index computation in this array;
|
|
they are documented \ref almgroup "here".
|
|
|
|
Information about a pixelisation of the sphere is stored in objects of
|
|
type psht_geom_info. It is possible to create such an object for any
|
|
supported pixelisation by using the function psht_make_geometry_info();
|
|
however, several easier-to-use functions are \ref geominfogroup "supplied"
|
|
for generating often-used pixelisations like ECP grids, Gaussian grids,
|
|
and Healpix grids.
|
|
|
|
A distinctive feature of PSHT is the ability to execute several transforms
|
|
simultaneously (as long as they have the same geometry information and
|
|
a_lm structure); this has the advantage that the spherical
|
|
harmonics only have to be evaluated once for all of the transforms, saving
|
|
a substantial amount of CPU time. The set of functions dealing with such
|
|
<i>job lists</i> is implemented twice, both for \ref sjoblistgroup "single"
|
|
and \ref djoblistgroup "double" precision transforms.
|
|
|
|
Currently, PSHT supports the following kinds of transforms:
|
|
<ul>
|
|
<li>scalar a_lm to map</li>
|
|
<li>scalar map to a_lm</li>
|
|
<li>polarised a_lm to map</li>
|
|
<li>polarised map to a_lm</li>
|
|
<li>spin(1-100) a_lm to map</li>
|
|
<li>spin(1-100) map to a_lm</li>
|
|
<li>scalar a_lm to maps of first derivatives</li>
|
|
</ul>
|
|
|
|
All these different kinds can be mixed freely within a job list; i.e., it is
|
|
possible to perform an a_lm->map and a map->a_lm transform simultaneously.
|
|
|
|
PSHT supports shared-memory parallelisation via OpenMP; this feature will
|
|
be automatically enabled if the compiler supports it.
|
|
|
|
PSHT will also make use of SSE2 instructions when compiled for a platform
|
|
known to support them (basically all Intel/AMD processors running a 64bit
|
|
operating system).
|
|
|
|
The spherical harmonic transforms can ce executed on double-precision and
|
|
single-precision maps and a_lm, but for accuracy reasons the computations
|
|
will always be performed in double precision. As a consequence,
|
|
single-precision transforms will most likely not be faster than their
|
|
double-precision counterparts, but they will require significantly less
|
|
memory.
|
|
|
|
Two example and benchmark programs are distributed with PSHT:
|
|
<ul>
|
|
<li>psht_test.c checks the accuracy of the (iterative) map analysis
|
|
algorithm</li>
|
|
<li>psht_perftest.c checks the performance of single or multiple
|
|
transforms</li>
|
|
</ul>
|
|
*/
|