Content-type: text/html
void forwardSHT( void *map, int inputIsComplex, coordStruct coords, int lmax, fftw_complex *alm );
#include ccSHTmpi.h
void forwardSHTmpi( void *map, int inputIsComplex, coordStruct coords, int lmax, fftw_complex *alm, int outputRoot, MPI_Comm theComm );
These functions are the members of the ccSHT library which perform the forward discrete spherical harmonic transform (SHT). forwardSHT() performs the computation in serial and forwardSHTmpi() performs the computation in parallel. For some background, and more information about the library see the man page for ccSHT.
forwardSHT() and forwardSHTmpi perform a band limited discrete SHT, where the band limit is set by the input parameter called lmax. More exactly, the transform computes the spherical harmonic coefficients (a = a_{l,m}) for all l in the set {0, 1, 2, ..., lmax} and m in the set {-l, -l+1, -l+2, ..., l}. The spherical harmonic coefficients are defined to be
/
| _
a = | f * Y dA
|
/
where
a = a_{l,m}
f = f(theta,phi)
_ _
Y = Y_{l,m}(theta, phi)
and the integral is over the surface of the unit sphere. In the above formulation, Y bar are the complex conjugates of the spherical harmonics. More explicitly, the complex conjugates of the spherical harmonics are defined
________________
_ / (2*l+1)*(l-m)!
Y = / ----------------*P(cos(theta))*exp(-i*m*phi)
\/ 4*pi*(l+m)!
where P(cos(theta)) = P_{l,m}(cos(theta)) are the associated Legendre polynomials, and i is the imaginary number. The above formulation is only valid for non-negative m since P is only defined for non-negative m. Y_{l,m} is related to Y_{l,-m} as follows:
m _
Y_{l,-m} = (-1) *Y_{l,m}
so we can infer the values of Y for negative m.
The input to the transform is a vector of evaluations of f at a set of positions on the sphere. The location of the positions at which the function f is evaluated are stored in the structure called coords (see man page ccSHT for more details about this structure). Let us refer to a vector of values of the spherical harmonics for a fixed l and m over all of the pixels on the sphere as a spherical harmonic vector. The integral is then calculated as a dot product in pixel space of the function vector and the spherical harmonic vector with a constant quadrature weight (coords.pixSize). If a more complicated quadrature scheme is preferred, simply set coords.pixSize to one and pass f multiplied by the preferred quadrature weighting to forwardSHT() in place of f.
For more information about the calculation of the associated Legendre polynomials see the man page for calculateQlm(). The transform does not actually call calculateQlm() for optimization purposes. Instead most of the coding required to do the calculation of the associated Legendre polynomials is done in the functions calculateSmallQlm() and lmCalculations() which the transform does call.
Below is a table describing the parameters passed to forwardSHT() and forwardSHTmpi().
Copyright (C) 2003 Christopher M. Cantalupo
forwardSHT is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
Christopher Cantalupo <cmc@nersc.gov>
Send bug reports or comments to the above address.