-
Notifications
You must be signed in to change notification settings - Fork 103
/
shexpanddh.1
142 lines (142 loc) · 5.75 KB
/
shexpanddh.1
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
.\" Automatically generated by Pandoc 2.17.1.1
.\"
.\" Define V font for inline verbatim, using C font in formats
.\" that render this, and otherwise B font.
.ie "\f[CB]x\f[]"x" \{\
. ftr V B
. ftr VI BI
. ftr VB B
. ftr VBI BI
.\}
.el \{\
. ftr V CR
. ftr VI CI
. ftr VB CB
. ftr VBI CBI
.\}
.TH "shexpanddh" "1" "2021-02-15" "Fortran 95" "SHTOOLS 4.10"
.hy
.SH SHExpandDH
.PP
Expand an equally sampled or equally spaced grid into spherical
harmonics using Driscoll and Healy\[cq]s (1994) sampling theorem.
.SH Usage
.PP
call SHExpandDH (\f[V]griddh\f[R], \f[V]n\f[R], \f[V]cilm\f[R],
\f[V]lmax\f[R], \f[V]norm\f[R], \f[V]sampling\f[R], \f[V]csphase\f[R],
\f[V]lmax_calc\f[R], \f[V]exitstatus\f[R])
.SH Parameters
.TP
\f[V]griddh\f[R] : input, real(dp), dimension (\f[V]n\f[R], \f[V]n\f[R]) or (\f[V]n\f[R], 2*\f[V]n\f[R])
A 2D equally sampled (default) or equally spaced grid that conforms to
the sampling theorem of Driscoll and Healy (1994).
The first latitudinal band corresponds to 90 N, the latitudinal band for
90 S is not included, and the latitudinal sampling interval is
180/\f[V]n\f[R] degrees.
The first longitudinal band is 0 E, the longitude band for 360 E is not
included, and the longitudinal sampling interval is 360/\f[V]n\f[R] for
an equally and 180/\f[V]n\f[R] for an equally spaced grid, respectively.
.TP
\f[V]n\f[R] : input, integer(int32)
The number of samples in latitude of \f[V]griddh\f[R].
If \f[V]sampling\f[R] is 1 (default) then the number of samples in
longitude is \f[V]n\f[R].
If \f[V]sampling\f[R] is 2 then the number of longitudinal samples is
\f[V]2n\f[R].
\f[V]n\f[R] must be even.
.TP
\f[V]cilm\f[R] : output, real(dp), dimension (2, \f[V]n\f[R]/2, \f[V]n\f[R]/2) or (2, \f[V]lmax_calc\f[R]+1, \f[V]lmax_calc\f[R]+1)
The real spherical harmonic coefficients of the function.
These will be exact if the function is bandlimited to degree
\f[V]lmax=n/2-1\f[R].
The coefficients \f[V]c1lm\f[R] and \f[V]c2lm\f[R] refer to the cosine
(\f[V]clm\f[R]) and sine (\f[V]slm\f[R]) coefficients, respectively,
with \f[V]clm=cilm(1,l+1,m+1)\f[R] and \f[V]slm=cilm(2,l+1,m+1)\f[R].
.TP
\f[V]lmax\f[R] : output, integer(int32)
The maximum spherical harmonic bandwidth of the input grid, which is
\f[V]n/2-1\f[R].
If the optional parameter \f[V]lmax_calc\f[R] is not specified, this
corresponds to the maximum spherical harmonic degree of the output
coefficients \f[V]cilm\f[R].
.TP
\f[V]norm\f[R] : input, optional, integer(int32), default = 1
1 (default) = 4-pi (geodesy) normalized harmonics; 2 = Schmidt
semi-normalized harmonics; 3 = unnormalized harmonics; 4 = orthonormal
harmonics.
.TP
\f[V]sampling\f[R] : input, optional, integer(int32), default = 1
If 1 (default) the input grid is equally sampled (\f[V]n\f[R] by
\f[V]n\f[R]).
If 2, the grid is equally spaced (\f[V]n\f[R] by \f[V]2n\f[R]).
.TP
\f[V]csphase\f[R] : input, optional, integer(int32), default = 1
1 (default) = do not apply the Condon-Shortley phase factor to the
associated Legendre functions; -1 = append the Condon-Shortley phase
factor of (-1)\[ha]m to the associated Legendre functions.
.TP
\f[V]lmax_calc\f[R] : input, optional, integer(int32), default = \f[V]lmax\f[R]
The maximum spherical harmonic degree calculated in the spherical
harmonic expansion.
.TP
\f[V]exitstatus\f[R] : output, optional, integer(int32)
If present, instead of executing a STOP when an error is encountered,
the variable exitstatus will be returned describing the error.
0 = No errors; 1 = Improper dimensions of input array; 2 = Improper
bounds for input variable; 3 = Error allocating memory; 4 = File IO
error.
.SH Description
.PP
\f[V]SHExpandDH\f[R] will expand an equally sampled (\f[V]n\f[R] by
\f[V]n\f[R]) or equally spaced grid (\f[V]n\f[R] by \f[V]2n\f[R]) into
spherical harmonics using the sampling theorem of Driscoll and Healy
(1994).
The number of latitudinal samples, \f[V]n\f[R], must be even, and the
transform is exact if the function is bandlimited to spherical harmonic
degree \f[V]n/2-1\f[R].
The inverse transform is given by the routine \f[V]MakeGridDH\f[R].
If the optional parameter \f[V]lmax_calc\f[R] is specified, the
spherical harmonic coefficients will only be calculated to this degree
instead of \f[V]n/2-1\f[R].
The algorithm is based on performing FFTs in longitude and then
integrating over latitude using an exact quadrature rule.
.PP
The default is to use an input grid that is equally sampled (\f[V]n\f[R]
by \f[V]n\f[R]), but this can be changed to use an equally spaced grid
(\f[V]n\f[R] by \f[V]2n\f[R]) by the optional argument
\f[V]sampling\f[R].
When using an equally spaced grid, the Fourier components corresponding
to degrees greater than \f[V]n/2-1\f[R] are simply discarded; this is
done to prevent aliasing that would occur if an equally sampled grid was
constructed from an equally spaced grid by discarding every other column
of the input grid.
.PP
The employed spherical harmonic normalization and Condon-Shortley phase
convention can be set by the optional arguments \f[V]norm\f[R] and
\f[V]csphase\f[R]; if not set, the default is to use geodesy 4-pi
normalized harmonics that exclude the Condon-Shortley phase of
(-1)\[ha]m.
The normalized legendre functions are calculated in this routine using
the scaling algorithm of Holmes and Featherstone (2002), which are
accurate to about degree 2800.
The unnormalized functions are accurate only to about degree 15.
.SH References
.PP
Driscoll, J.R.
and D.M.
Healy, Computing Fourier transforms and convolutions on the 2-sphere,
Adv.
Appl.
Math., 15, 202-250, 1994.
.PP
Holmes, S.
A., and W.
E.
Featherstone, A unified approach to the Clenshaw summation and the
recursive computation of very high degree and order normalised
associated Legendre functions, J.
Geodesy, 76, 279-299, 2002.
.SH See also
.PP
makegriddh, makegriddhc, shexpanddhc, makegridglq, shexpandglq,
makegridglqc, shexpandglqc, shexpandlsq