/
makegravgradgriddh.1
226 lines (226 loc) · 7.4 KB
/
makegravgradgriddh.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
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
.TH "makegravgradgriddh" "1" "2015\-04\-26" "Fortran 95" "SHTOOLS 3.1"
.SH MakeGravGradGridDH
.PP
Create 2D cylindrical maps on a flattened ellipsoid of the components of
the gravity "gradient" tensor in a local north\-oriented reference
frame.
.SH Usage
.PP
call MakeGravGradGridDH (\f[C]cilm\f[], \f[C]lmax\f[], \f[C]gm\f[],
\f[C]r0\f[], \f[C]a\f[], \f[C]f\f[], \f[C]vxx\f[], \f[C]vyy\f[],
\f[C]vzz\f[], \f[C]vxy\f[], \f[C]vxz\f[], \f[C]vyz\f[], \f[C]n\f[],
\f[C]sampling\f[], \f[C]lmax_calc\f[])
.SH Parameters
.TP
.B \f[C]cilm\f[] : input, real*8, dimension (2, \f[C]lmax\f[]+1, \f[C]lmax\f[]+1)
The real gravitational potential spherical harmonic coefficients.
The coefficients \f[C]c1lm\f[] and \f[C]c2lm\f[] refer to the cosine and
sine coefficients, respectively, with \f[C]c1lm=cilm(1,l+1,m+1)\f[] and
\f[C]c2lm=cilm(2,l+1,m+1)\f[].
.RS
.RE
.TP
.B \f[C]lmax\f[] : input, integer
The maximum spherical harmonic degree of the coefficients \f[C]cilm\f[].
This determines the number of samples of the output grids,
\f[C]n=2lmax+2\f[], and the latitudinal sampling interval,
\f[C]90/(lmax+1)\f[].
.RS
.RE
.TP
.B \f[C]gm\f[] : input, real*8
The gravitational constant multiplied by the mass of the planet.
.RS
.RE
.TP
.B \f[C]r0\f[]: input, real*8
The reference radius of the spherical harmonic coefficients.
.RS
.RE
.TP
.B \f[C]a\f[] : input, real*8
The semi\-major axis of the flattened ellipsoid on which the field is
computed.
.RS
.RE
.TP
.B \f[C]f\f[] : input, real*8
The flattening of the reference ellipsoid:
\f[C]f=(R_equator\-R_pole)/R_equator\f[].
.RS
.RE
.TP
.B \f[C]vxx\f[] : output, real*8, dimension (2*\f[C]lmax\f[]+2, \f[C]sampling\f[]*(2*\f[C]lmax\f[]+2))
A 2D equally sampled (\f[C]n\f[] by \f[C]n\f[]) or equally spaced
(\f[C]n\f[] by 2\f[C]n\f[]) grid of the \f[C]xx\f[] component of the
gravity tensor.
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[C]n\f[] degrees.
The first longitudinal band is 0 E, the longitudinal band for 360 E is
not included, and the longitudinal sampling interval is 360/\f[C]n\f[]
for an equally sampled and 180/\f[C]n\f[] for an equally spaced grid,
respectively.
.RS
.RE
.TP
.B \f[C]vyy\f[] : output, real*8, dimension (2*\f[C]lmax\f[]+2, \f[C]sampling\f[]*(2*\f[C]lmax\f[]+2))
A 2D equally sampled or equally spaced grid of the \f[C]yy\f[] component
of the gravity tensor.
.RS
.RE
.TP
.B \f[C]vzz\f[] : output, real*8, dimension (2*\f[C]lmax\f[]+2, \f[C]sampling\f[]*(2*\f[C]lmax\f[]+2))
A 2D equally sampled or equally spaced grid of the \f[C]zz\f[] component
of the gravity tensor.
.RS
.RE
.TP
.B \f[C]vxy\f[] : output, real*8, dimension (2*\f[C]lmax\f[]+2, \f[C]sampling\f[]*(2*\f[C]lmax\f[]+2))
A 2D equally sampled or equally spaced grid of the \f[C]xy\f[] component
of the gravity tensor.
.RS
.RE
.TP
.B \f[C]vxz\f[] : output, real*8, dimension (2*\f[C]lmax\f[]+2, \f[C]sampling\f[]*(2*\f[C]lmax\f[]+2))
A 2D equally sampled or equally spaced grid of the \f[C]xz\f[] component
of the gravity tensor.
.RS
.RE
.TP
.B \f[C]vyz\f[] : output, real*8, dimension (2*\f[C]lmax\f[]+2, \f[C]sampling\f[]*(2*\f[C]lmax\f[]+2))
A 2D equally sampled or equally spaced grid of the YZ component of the
gravity tensor.
.RS
.RE
.TP
.B \f[C]n\f[] : output, integer
The number of samples in latitude of the output grids.
This is equal to \f[C]2lmax+2\f[].
.RS
.RE
.TP
.B \f[C]sampling\f[] : optional, input, integer, default = 1
If 1 (default) the output grids are equally sampled (\f[C]n\f[] by
\f[C]n\f[]).
If 2, the grids are equally spaced (\f[C]n\f[] by 2\f[C]n\f[]).
.RS
.RE
.TP
.B \f[C]lmax_calc\f[] : optional, input, integer
The maximum spherical harmonic degree used in evaluating the functions.
This must be less than or equal to \f[C]lmax\f[].
.RS
.RE
.SH Description
.PP
\f[C]MakeGravGradGridDH\f[] will create 2\-dimensional cylindrical maps
from the spherical harmonic coefficients \f[C]cilm\f[], equally sampled
(\f[C]n\f[] by \f[C]n\f[]) or equally spaced (\f[C]n\f[] by 2\f[C]n\f[])
in latitude and longitude, for six components of the gravity "gradient"
tensor (all using geocentric coordinates):
.PP
\f[C](Vxx,\ \ Vxy,\ \ \ \ Vxz)\f[]
.PD 0
.P
.PD
\f[C](Vyx,\ \ Vyy,\ \ \ \ Vyz)\f[]
.PD 0
.P
.PD
\f[C](Vzx,\ \ Vzy,\ \ \ \ Vzz)\f[]
.PP
The reference frame is north\-oriented, where \f[C]x\f[] points north,
\f[C]y\f[] points west, and \f[C]z\f[] points upward (all tangent or
perpendicular to a sphere of radius r).
The gravitational potential is defined as
.PP
\f[C]V\ =\ GM/r\ Sum_{l=0}^lmax\ (r0/r)^l\ Sum_{m=\-l}^l\ C_{lm}\ Y_{lm}\f[],
.PP
where \f[C]r0\f[] is the reference radius of the spherical harmonic
coefficients \f[C]Clm\f[], and the gravitational acceleration is
.PP
\f[C]B\ =\ Grad\ V\f[].
.PP
The gravity tensor is symmetric, and satisfies \f[C]Vxx+Vyy+Vzz=0\f[],
though all three diagonal elements are calculated independently in this
routine.
.PP
The components of the gravity tensor are calculated according to eq.
1 in Petrovskaya and Vershkov (2006), which is based on eq.
3.28 in Reed (1973) (noting that Reed\[aq]s equations are in terms of
latitude and that the \f[C]y\f[] axis points east):
.PP
\f[C]Vzz\ =\ Vrr\f[]
.PD 0
.P
.PD
\f[C]Vxx\ =\ 1/r\ Vr\ +\ 1/r^2\ Vtt\f[]
.PD 0
.P
.PD
\f[C]Vyy\ =\ 1/r\ Vr\ +\ 1/r^2\ /tan(t)\ Vt\ +\ 1/r^2\ /sin(t)^2\ Vpp\f[]
.PD 0
.P
.PD
\f[C]Vxy\ =\ 1/r^2\ /sin(t)\ Vtp\ \-\ cos(t)/sin(t)^2\ /r^2\ Vp\f[]
.PD 0
.P
.PD
\f[C]Vxz\ =\ 1/r^2\ Vt\ \-\ 1/r\ Vrt\f[]
.PD 0
.P
.PD
\f[C]Vyz\ =\ 1/r^2\ /sin(t)\ Vp\ \-\ 1/r\ /sin(t)\ Vrp\f[]
.PP
where \f[C]r\f[], \f[C]t\f[], \f[C]p\f[] stand for radius, theta, and
phi, respectively, and subscripts on \f[C]V\f[] denote partial
derivatives.
.PP
The output grid are in units of s^\-2 and are cacluated on a flattened
ellipsoid with semi\-major axis \f[C]a\f[] and flattening \f[C]f\f[].
To obtain units of Eotvos (10^\-9 s^\-2), multiply the output by 10^9.
The calculated values should be considered exact only when the radii on
the ellipsoid are less than the maximum radius of the planet (the
potential coefficients are simply downward continued in the spectral
domain).
.PP
The default is to calculate grids for use in the Driscoll and Healy
(1994) routines that are equally sampled (\f[C]n\f[] by \f[C]n\f[]), but
this can be changed to calculate equally spaced grids (\f[C]n\f[] by
2\f[C]n\f[]) by setting the optional argument \f[C]sampling\f[] to 2.
The input value of \f[C]lmax\f[] determines the number of samples,
\f[C]n=2lmax+2\f[], and the latitudinal sampling interval,
90/(\f[C]lmax\f[]+1).
The first latitudinal band of the grid corresponds to 90 N, the
latitudinal band for 90 S is not calculated, and the latitudinal
sampling interval is 180/\f[C]n\f[] degrees.
The first longitudinal band is 0 E, the longitudinal band for 360 E is
not calculated, and the longitudinal sampling interval is 360/\f[C]n\f[]
for equally sampled and 180/\f[C]n\f[] for equally spaced grids,
respectively.
.SH References
.PP
Reed, G.B., Application of kinematical geodesy for determining the short
wave length components of the gravity field by satellite gradiometry,
Ohio State University, Dept.
of Geod.
Sciences, Rep.
No.
201, Columbus, Ohio, 1973.
.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
Petrovskaya, M.S.
and A.N.
Vershkov, Non\-singular expressions for the gravity gradients in the
local north\-oriented and orbital reference frames, J.
Geod., 80, 117\-127, 2006.
.SH See also
.PP
makegravgriddh, makegeoidgrid, makegriddh