charm_sha#
Module to perform surface spherical harmonic analysis using:
point data values,
mean data values.
Note
This documentation is written for double precision version of CHarm.
Analysis of point values
Functions to perform spherical harmonic analysis of point values.
-
void charm_sha_point(const charm_point *pnt, const double *f, unsigned long nmax, charm_shc *shcs, charm_err *err)#
Performs the exact surface spherical harmonic analysis based on point data values:
\[\begin{split}\left.\begin{aligned} \bar{C}_{nm} \\ \bar{S}_{nm} \end{aligned}\right\} = \dfrac{1}{4 \pi} \, \dfrac{R}{\mu} \, \left( \dfrac{r}{R} \right)^n \displaystyle\iint_{\sigma} f(r, \varphi, \lambda) \, \bar{Y}_{nm}(\varphi,\lambda) \, \mathrm{d}\sigma {.}\end{split}\]The coefficients are automatically normalized to the constants \(\mu\) and \(R\), which are taken from
shcs->mu
andshcs->r
, respectively. If \(r \neq R\) (pnt->r[i] != shcs->r
), the coefficients are also rescaled from the data sphere with the radius \(r\) to the reference sphere with the radius \(R\) usingcharm_shc_rescale()
. If \(r = R = \mu = 1\) (pnt->r[i] == shcs->r == shcs->mu == 1.0
), one gets the well-known relation for spherical harmonic analysis on the unit sphere (no normalization, no rescaling).The integrals are evaluated by exact quadratures, either the Gauss-Legendre or the Driscoll-Healy one. The type of the quadrature to be applied is determined by
pnt->type
.The function is parallelized using OpenMP.
The function exploits the symmetry property of Legendre functions with respect to the equator, thereby improving its performance. For further details on the grid requirements, see
charm_shs_point()
.References:
Driscoll, J. R., Healy, D. M. (1994) Computing Fourier transforms and convolutions on the 2-sphere. Advances in Applied Mathematics 15:202-250
Sneeuw, N. (1994) Global spherical harmonic analysis by least-squares and numerical quadrature methods in historical perspective. Geophysical Journal International 118:707-716
Heiskannen, W. A., Moritz, H. (1967) Physical Geodesy. W. H. Freeman and Company, San Francisco, 364 pp
Fukushima, T. (2012) Numerical computation of spherical harmonics of arbitrary degree and order by extending exponent of floating point numbers. Journal of Geodesy 86:271-285
Note
The function modifies
shcs->c
andshcs->s
, but it does not touch the remaining members ofshcs
. Ifshcs->nmax > nmax
, the coefficients beyondnmax
are set to zero.- Parameters:
pnt – [in] Structure defining grid points and integration weights (see
charm_point
). Thepnt
structure must be prepared by calling thecharm_crd_point_gl()
,charm_crd_point_dh1()
orcharm_crd_point_dh2()
with somenmax_signal
. Importantly, it is assumed that no harmonics (frequencies) are present in the input signal beyondnmax_signal
. Otherwise, the output spherical harmonic coefficients will suffer from the alliasing effect.f – [in] Input signal with point data given at the nodes of the
pnt
structure and containing harmonics up to degreenmax_signal
at most. The pointerf
must have an access topnt->nlat * pnt->nlon
array elements. The value off
atpnt->lat[i]
andpnt->lon[j]
must be stored asf[i * pnt->nlon + j]
withi = 0, 1, ..., pnt->nlat - 1
andj = 0, 1, ..., pnt->nlon - 1
.nmax – [in] Maximum harmonic degree of the output spherical harmonic coefficients. If
f
is properly given at the nodes of the Gauss-Legendre or the Driscoll-Healy grid fornmax_signal
and no harmonics are present beyondnmax_signal
inf
, than the output coefficients are computed exactly for allnmax <= nmax_signal
.shcs – [inout] Output spherical harmonic coefficients (see
charm_shc
).err – [out] Error reported by the function (if any).
Analysis of mean values
Functions to perform spherical harmonic analysis of mean values.
-
void charm_sha_cell(const charm_cell *cell, const double *f, unsigned long nmax, int method, charm_shc *shcs, charm_err *err)#
Performs an approximate surface spherical harmonic analysis based on mean data values (area means, block-means) given over grid cells:
\[\begin{split}\left.\begin{aligned} \bar{C}_{nm} \\ \bar{S}_{nm} \end{aligned} \right\} =& \dfrac{1}{4 \pi} \, \dfrac{R}{\mu} \, \left( \dfrac{r}{R} \right)^n \, \sum_{i,j} \tilde{f}(\varphi_i^{\min}, \varphi_i^{\max}, \lambda_j^{\min}, \lambda_j^{\max}) \\ &\times \int\limits_{\varphi_i^\mathrm{\min}}^{\varphi_i^\mathrm{\max}} \bar{P}_{nm}(\sin\varphi) \, \cos\varphi \, \mathrm{d}\varphi \, \int\limits_{\lambda_j^\mathrm{\min}}^{\lambda_j^\mathrm{\max}} \begin{Bmatrix} \cos m\lambda \\ \sin m\lambda \end{Bmatrix} \mathrm{d}\lambda{.}\end{split}\]See also the comments from
charm_sha_point()
on the normalization and rescaling of the coefficients with \(R\), \(r\) and \(\mu\).Applied is the approximate numerical quadrature due to Colombo (1981).
The function is parallelized using OpenMP.
A note on grid cells
Since the quadrature is approximate, the grid restrictions are more relaxed when compared with
charm_sha_point()
. Yet, some conditions has to be fulfilled:at least
nmax + 1
data values inf
are required in the latitudinal direction,at least
2 * nmax + 1
data values inf
are required in the longitudinal direction,cell->latmax[0]
has to be equal toPI / 2.0
(within thecharm_glob_threshold
limit, see charm_glob),cell->latmin[cell->nlat - 1]
has to be equal to-PI / 2.0
(again, within thecharm_glob_threshold
limit),cell->lonmin[0]
has to be equal to0.0
(still within thecharm_glob_threshold
limit),cell->lonmax[cell->nlon - 1]
has to be equal to2.0 * PI
(still within thecharm_glob_threshold
limit),the maximum latitude of one cell has to be equal to the minimum latitude of the next cell that follows (the same holds also for the longitudes; both within the
charm_glob_threshold
limit).
If possible, the function exploits the symmetry property of Legendre functions and their integrals with respect to the equator, thereby improving its performance (see
charm_shs_cell()
).References:
Colombo, O. L. (1981) Numerical methods for harmonic analysis on the sphere. Reports of the Department of Geodetic Science, Report No. 310, 140 pp
Heiskannen, W. A., Moritz, H. (1967) Physical Geodesy. W. H. Freeman and Company, San Francisco, 364 pp
Jekeli, C., Lee, J. K., Kwon, J. H. (2007) On the computation and approximation of ultra-high-degree spherical harmonic series. Journal of Geodesy 81:603-615
Fukushima, T. (2012) Numerical computation of spherical harmonics of arbitrary degree and order by extending exponent of floating point numbers. Journal of Geodesy 86:271-285.
Balmino, Vales, Bonvalot, Briais, 2012. Spherical harmonic modelling to ultra-high degree of Bouguer and isostatic anomalies. Journal of Geodesy 86:499-520, doi: 10.1007/s00190-011-0533-4
Fukushima, T. (2014) Numerical computation of spherical harmonics of arbitrary degree and order by extending exponent of floating point numbers: III integral. Computers & Geosciences 63:17-21
Note
The function modifies
shcs->c
andshcs->s
, but it does not touch the remaining members ofshcs
. Ifshcs->nmax > nmax
, the coefficients beyondnmax
are set to zero.- Parameters:
cell – [in] Structure defining grid cells (see
charm_cell
).f – [in] Input signal with block mean data given at the nodes of the
cell
structure. The pointerf
must have an access tocell->nlat * cell->nlon
array elements. The mean value of the signal at thei
-th cell in the latitudinal direction and thej
-th cell in the longitudinal direction must be stored asf[i * cell->nlon + j]
withi = 0, 1, ..., cell->nlat 1
andj = 0, 1, ..., pnt->nlon 1
.nmax – [in] Maximum harmonic degree of the output spherical harmonic coefficients.
method – [in] Method of spherical harmonic analysis of cell data. The only method currently supported is the approximate quadrature
CHARM_SHA_CELL_AQ
.shcs – [inout] Output spherical harmonic coefficients (see
charm_shc
).err – [out] Error reported by the function (if any).