-
Notifications
You must be signed in to change notification settings - Fork 106
/
makemaggradgriddh.1
209 lines (209 loc) · 7.79 KB
/
makemaggradgriddh.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
.\" Automatically generated by Pandoc 2.9.2
.\"
.TH "makemaggradgriddh" "1" "2020-01-17" "Fortran 95" "SHTOOLS 4.6"
.hy
.SH MakeMagGradGridDH
.PP
Create 2D cylindrical maps on a flattened ellipsoid of the components of
the magnetic field tensor in a local north-oriented reference frame.
.SH Usage
.PP
call MakeMagGradGridDH (\f[C]cilm\f[R], \f[C]lmax\f[R], \f[C]r0\f[R],
\f[C]a\f[R], \f[C]f\f[R], \f[C]vxx\f[R], \f[C]vyy\f[R], \f[C]vzz\f[R],
\f[C]vxy\f[R], \f[C]vxz\f[R], \f[C]vyz\f[R], \f[C]n\f[R],
\f[C]sampling\f[R], \f[C]lmax_calc\f[R], \f[C]extend\f[R],
\f[C]exitstatus\f[R])
.SH Parameters
.TP
\f[B]\f[CB]cilm\f[B]\f[R] : input, real(dp), dimension (2, \f[B]\f[CB]lmax\f[B]\f[R]+1, \f[B]\f[CB]lmax\f[B]\f[R]+1)
The real Schmidt semi-normalized spherical harmonic coefficients of the
magnetic potential.
The coefficients \f[C]c1lm\f[R] and \f[C]c2lm\f[R] refer to the cosine
and sine coefficients, respectively, with \f[C]c1lm=cilm(1,l+1,m+1)\f[R]
and \f[C]c2lm=cilm(2,l+1,m+1)\f[R].
The coefficients are assumed to have units of nT.
.TP
\f[B]\f[CB]lmax\f[B]\f[R] : input, integer
The maximum spherical harmonic degree of the coefficients
\f[C]cilm\f[R].
This determines the number of samples of the output grids,
\f[C]n=2lmax+2\f[R], and the latitudinal sampling interval,
\f[C]90/(lmax+1)\f[R].
.TP
\f[B]\f[CB]r0\f[B]\f[R]: input, real(dp)
The reference radius of the spherical harmonic coefficients.
.TP
\f[B]\f[CB]a\f[B]\f[R] : input, real(dp)
The semi-major axis of the flattened ellipsoid on which the field is
computed.
.TP
\f[B]\f[CB]f\f[B]\f[R] : input, real(dp)
The flattening of the reference ellipsoid:
\f[C]f=(R_equator-R_pole)/R_equator\f[R].
.TP
\f[B]\f[CB]vxx\f[B]\f[R] : output, real(dp), dimension (nlat, nlong)
A 2D map of the \f[C]xx\f[R] component of the magnetic field tensor that
conforms to the sampling theorem of Driscoll and Healy (1994).
If \f[C]sampling\f[R] is 1, the grid is equally sampled and is
dimensioned as (\f[C]n\f[R] by \f[C]n\f[R]), where \f[C]n\f[R] is
\f[C]2lmax+2\f[R].
If sampling is 2, the grid is equally spaced and is dimensioned as
(\f[C]n\f[R] by 2\f[C]n\f[R]).
The first latitudinal band of the grid corresponds to 90 N, the
latitudinal sampling interval is 180/\f[C]n\f[R] degrees, and the
default behavior is to exclude the latitudinal band for 90 S.
The first longitudinal band of the grid is 0 E, by default the
longitudinal band for 360 E is not included, and the longitudinal
sampling interval is 360/\f[C]n\f[R] for an equally sampled and
180/\f[C]n\f[R] for an equally spaced grid, respectively.
If \f[C]extend\f[R] is 1, the longitudinal band for 360 E and the
latitudinal band for 90 S will be included, which increases each of the
dimensions of the grid by 1.
.TP
\f[B]\f[CB]vyy\f[B]\f[R] : output, real(dp), dimension (nlat, nlong)
A 2D equally sampled or equally spaced grid of the \f[C]yy\f[R]
component of the magnetic field tensor.
.TP
\f[B]\f[CB]vzz\f[B]\f[R] : output, real(dp), dimension (nlat, nlong))
A 2D equally sampled or equally spaced grid of the \f[C]zz\f[R]
component of the magnetic field tensor.
.TP
\f[B]\f[CB]vxy\f[B]\f[R] : output, real(dp), dimension (nlat, nlong)
A 2D equally sampled or equally spaced grid of the \f[C]xy\f[R]
component of the magnetic field tensor.
.TP
\f[B]\f[CB]vxz\f[B]\f[R] : output, real(dp), dimension (nlat, nlong)
A 2D equally sampled or equally spaced grid of the \f[C]xz\f[R]
component of the magnetic field tensor.
.TP
\f[B]\f[CB]vyz\f[B]\f[R] : output, real(dp), dimension (nlat, nlong)
A 2D equally sampled or equally spaced grid of the YZ component of the
magnetic field tensor.
.TP
\f[B]\f[CB]n\f[B]\f[R] : output, integer
The number of samples in latitude of the output grids.
This is equal to \f[C]2lmax+2\f[R].
.TP
\f[B]\f[CB]sampling\f[B]\f[R] : optional, input, integer, default = 1
If 1 (default) the output grids are equally sampled (\f[C]n\f[R] by
\f[C]n\f[R]).
If 2, the grids are equally spaced (\f[C]n\f[R] by 2\f[C]n\f[R]).
.TP
\f[B]\f[CB]lmax_calc\f[B]\f[R] : 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[R].
.TP
\f[B]\f[CB]extend\f[B]\f[R] : input, optional, integer, default = 0
If 1, compute the longitudinal band for 360 E and the latitudinal band
for 90 S.
This increases each of the dimensions of \f[C]griddh\f[R] by 1.
.TP
\f[B]\f[CB]exitstatus\f[B]\f[R] : output, optional, integer
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[C]MakeMagGradGridDH\f[R] will create 2-dimensional cylindrical maps
from the spherical harmonic coefficients \f[C]cilm\f[R], equally sampled
(\f[C]n\f[R] by \f[C]n\f[R]) or equally spaced (\f[C]n\f[R] by
2\f[C]n\f[R]) in latitude and longitude, for six components of the
magnetic field tensor (all using geocentric coordinates):
.PP
\f[C](Vxx, Vxy, Vxz)\f[R]
.PD 0
.P
.PD
\f[C](Vyx, Vyy, Vyz)\f[R]
.PD 0
.P
.PD
\f[C](Vzx, Vzy, Vzz)\f[R]
.PP
The reference frame is north-oriented, where \f[C]x\f[R] points north,
\f[C]y\f[R] points west, and \f[C]z\f[R] points upward (all tangent or
perpendicular to a sphere of radius r).
The magnetic potential is defined as
.PP
\f[C]V = r0 Sum_{l=0}\[ha]lmax (r0/r)\[ha](l+1) Sum_{m=-l}\[ha]l C_{lm} Y_{lm}\f[R],
.PP
where \f[C]r0\f[R] is the reference radius of the spherical harmonic
coefficients \f[C]Clm\f[R], and the vector magnetic field is
.PP
\f[C]B = - Grad V\f[R].
.PP
The magnetic field tensor is symmetric, and satisfies
\f[C]Vxx+Vyy+Vzz=0\f[R], though all three diagonal elements are
calculated independently in this routine.
.PP
The components of the magnetic field 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\[cq]s equations are in terms of
latitude and that the \f[C]y\f[R] axis points east):
.PP
\f[C]Vzz = Vrr\f[R]
.PD 0
.P
.PD
\f[C]Vxx = 1/r Vr + 1/r\[ha]2 Vtt\f[R]
.PD 0
.P
.PD
\f[C]Vyy = 1/r Vr + 1/r\[ha]2 /tan(t) Vt + 1/r\[ha]2 /sin(t)\[ha]2 Vpp\f[R]
.PD 0
.P
.PD
\f[C]Vxy = 1/r\[ha]2 /sin(t) Vtp - cos(t)/sin(t)\[ha]2 /r\[ha]2 Vp\f[R]
.PD 0
.P
.PD
\f[C]Vxz = 1/r\[ha]2 Vt - 1/r Vrt\f[R]
.PD 0
.P
.PD
\f[C]Vyz = 1/r\[ha]2 /sin(t) Vp - 1/r /sin(t) Vrp\f[R]
.PP
where \f[C]r\f[R], \f[C]t\f[R], \f[C]p\f[R] stand for radius, theta, and
phi, respectively, and subscripts on \f[C]V\f[R] denote partial
derivatives.
.PP
The output grid are in units of nT / m and are cacluated on a flattened
ellipsoid with semi-major axis \f[C]a\f[R] and flattening \f[C]f\f[R].
The calculated values should be considered exact only when the radii on
the ellipsoid are greater than the maximum radius of the planet (the
potential coefficients are simply downward/upward continued in the
spectral domain).
.PP
The default is to use an input grid that is equally sampled (\f[C]n\f[R]
by \f[C]n\f[R]), but this can be changed to use an equally spaced grid
(\f[C]n\f[R] by 2\f[C]n\f[R]) by the optional argument
\f[C]sampling\f[R].
The redundant longitudinal band for 360 E and the latitudinal band for
90 S are excluded by default, but these can be computed by specifying
the optional argument \f[C]extend\f[R].
.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
makemaggriddh, makegriddh