-
Notifications
You must be signed in to change notification settings - Fork 106
/
batohilmrhoh.1
169 lines (169 loc) · 7.86 KB
/
batohilmrhoh.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
.\" Automatically generated by Pandoc 2.7.3
.\"
.TH "batohilmrhoh" "1" "2019-09-17" "Fortran 95" "SHTOOLS 4.5"
.hy
.SH BAtoHilmRhoH
.PP
Calculate iteratively the relief along an interface with lateral density
variations that corresponds to a given Bouguer anomaly using the
algorithm of Wieczorek and Phillips (1998).
.SH Usage
.PP
call BAtoHilmRhoH (\f[C]cilm\f[R], \f[C]ba\f[R], \f[C]grid\f[R],
\f[C]lmax\f[R], \f[C]nmax\f[R], \f[C]mass\f[R], \f[C]r0\f[R],
\f[C]rho\f[R], \f[C]gridtype\f[R], \f[C]w\f[R], \f[C]plx\f[R],
\f[C]zero\f[R], \f[C]filtertype\f[R], \f[C]filterdeg\f[R],
\f[C]lmaxcalc\f[R], \f[C]exitstatus\f[R])
.SH Parameters
.TP
.B \f[C]cilm\f[R] : output, real(dp), dimension (2, \f[C]lmaxcalc\f[R]+1, \f[C]lmaxcalc\f[R]+1)
An estimate of the real spherical harmonic coefficients (geodesy
normalized) of relief along an interface with lateraly varying density
contrast \f[C]rho\f[R] that satisfies the Bouguer anomaly \f[C]ba\f[R].
The degree-zero term corresponds to the mean radius of the relief.
.TP
.B \f[C]ba\f[R] : input, real(dp), dimension (2, \f[C]lmaxcalc\f[R]+1, \f[C]lmaxcalc\f[R]+1)
The real spherical harmonic coefficients of the Bouguer anomaly
referenced to a spherical interface \f[C]r0\f[R].
.TP
.B \f[C]grid\f[R] : input, real(dp), dimension (\f[C]lmax\f[R]+1, 2*\f[C]lmax\f[R]+1) for \f[C]gridtype\f[R] 1, (2*\f[C]lmax\f[R]+2, 2*\f[C]lmax\f[R]+2) for \f[C]gridtype\f[R] 2, (2*\f[C]lmax\f[R]+2, 4*\f[C]lmax\f[R]+4) for \f[C]gridtype\f[R] 3
The initial estimate for the radii of the interface evaluated on a grid
corresponding to a function of maximum spherical harmonic degree
\f[C]lmax\f[R].
This is calculated by a call to either \f[C]MakeGridGLQ\f[R] or
\f[C]MakeGridDH\f[R].
This grid must contain the degree-0 average radius of the interface.
.TP
.B \f[C]lmax\f[R] : input, integer
The spherical harmonic bandwidth of the input relief \f[C]grid\f[R],
which determines the dimensions of \f[C]grid\f[R].
If \f[C]lmaxcalc\f[R] is not set, this determines also the maximum
spherical harmonic degree of the output spherical harmonic coefficients
of the relief and the input spherical harmonics of the Bouguer anomaly.
.TP
.B \f[C]nmax\f[R] : input, integer
The maximum order used in the Taylor-series expansion used in
calculating the potential coefficients.
.TP
.B \f[C]mass\f[R] : input, real(dp)
The mass of the planet in kg.
.TP
.B \f[C]r0\f[R] : input, real(dp)
The reference radius of the Bouguer anomaly \f[C]ba\f[R].
.TP
.B \f[C]rho\f[R] : input, real(dp), dimension (\f[C]lmax\f[R]+1, 2*\f[C]lmax\f[R]+1) for \f[C]gridtype\f[R] 1, (2*\f[C]lmax\f[R]+2, 2*\f[C]lmax\f[R]+2) for \f[C]gridtype\f[R] 2, (2*\f[C]lmax\f[R]+2, 4*\f[C]lmax\f[R]+4) for \f[C]gridtype\f[R] 3
The density contrast of the relief in kg/m\[ha]3, with the same
dimensions as \f[C]grid\f[R].
.TP
.B \f[C]gridtype\f[R] : input, integer
1 = Gauss-Legendre grids, calculated using \f[C]SHGLQ\f[R] and
\f[C]MakeGridGLQ\f[R].
2 = Equally sampled Driscoll-Healy grids, \f[C]n\f[R] by \f[C]n\f[R],
calculated using \f[C]MakeGridDH\f[R].
3 = Equally spaced Driscoll-Healy grids, \f[C]n\f[R] by 2\f[C]n\f[R],
calculated using \f[C]MakeGridDH\f[R].
.TP
.B \f[C]w\f[R] : optional, input, real(dp), dimension (\f[C]lmax\f[R]+1)
The weights used in the Gauss-Legendre quadrature.
These are calculated from a call to \f[C]SHGLQ\f[R].
If present, one of \f[C]plx\f[R] or \f[C]zero\f[R] must also be present.
.TP
.B \f[C]plx\f[R] : optional, input, real(dp), dimension (\f[C]lmax\f[R]+1, (\f[C]lmax\f[R]+1)*(\f[C]lmax\f[R]+2)/2)
An array of the associated Legendre functions calculated at the nodes
used in the Gauss-Legendre quadrature.
These are determined from a call to \f[C]SHGLQ\f[R].
.TP
.B \f[C]zero\f[R] : optional, input, real(dp), dimension (\f[C]lmax\f[R]+1)
The nodes used in the Gauss-Legendre quadrature over latitude,
calculated by a call to \f[C]SHGLQ\f[R].
.TP
.B \f[C]filtertype\f[R] : optional, input, integer, default = 0
Apply a filter when calculating the relief in order to minimize the
destabilizing effects of downward continuation which amplify
uncertainties in the Bouguer anomaly.
If 0, no filtering is applied.
If 1, use the minimum amplitude filter \f[C]DownContFilterMA\f[R].
If 2, use the minimum curvature filter \f[C]DownContFilterMC\f[R].
.TP
.B \f[C]filterdeg\f[R] : optional, input, integer
The spherical harmonic degree for which the filter is 0.5.
.TP
.B \f[C]lmaxcalc\f[R] : optional, input, integer, default = \f[C]lmax\f[R]
The maximum degree that will be calculated in the spherical harmonic
expansions.
.TP
.B \f[C]exitstatus\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]BAtoHilmRhoH\f[R] is used to solve iteratively for the relief along
an interface with lateral variations in density that corresponds to a
given Bouguer anomaly.
This is equation 18 of Wieczorek and Phillips (1998), modified to
account for density variations as in equation 30 of Wieczorek 2007,
which implicitly takes into consideration the finite-amplitude
correction.
Each iteration takes as input a guess for the relief and outputs the
iteratively improved spherical harmonic coefficients of this relief.
These coefficients can then be re-expanded and re-input into this
routine as the next guess.
For the initial guess, it is often sufficient to use the relief
predicted using the first-order \[lq]mass sheet\[rq] approximation.
The input relief \f[C]grid\f[R] can be of one of three type specified by
\f[C]gridtype\f[R]: 1 for Gauss-Legendre grids, 2 for equally sampled
Driscoll-Healy grids (\f[C]n\f[R] by \f[C]n\f[R]), and 3 for equally
spaced Driscoll-Healy grids (\f[C]n\f[R] by 2\f[C]n\f[R]).
.PP
If the algorithm does not converge, one might want to try damping the
initial estimate.
Alternatively, iterations of the following form have proven successfulin
in damping oscilations between successive iterations:
.PP
\f[C]h3 = (h2+h1)/2\f[R] \f[C]h4 = f(h3)\f[R]
.PP
It is important to understand that as an intermediate step, this routine
calculates the spherical harmonic coefficients of the density multiplied
by the relief raised to the nth power.
As such, if the input function is bandlimited to degree \f[C]L\f[R], the
resulting function will thus be bandlimited to degree \f[C]L*nmax\f[R].
This subroutine implicitly assumes that \f[C]lmax\f[R] is greater than
or equal to \f[C]L*nmax\f[R].
If this is not the case, then aliasing will occur.
In practice, for accurate results, it is found that \f[C]lmax\f[R] needs
only to be about twice the size of \f[C]L\f[R], though this should be
verified for each application.
Thus, if the input function is considered to be bandlimited to degree
\f[C]L\f[R], the function should be evaluated on a grid corresponding to
a maximum degree of about \f[C]2L\f[R].
.PP
If the input grid is evaluated on the Gauss-Legendre points, it is
necessary to specify the optional parameters \f[C]w\f[R]
and\f[C]zero\f[R], or \f[C]w\f[R] and \f[C]plx\f[R], which are
calculated by a call to \f[C]SHGLQ\f[R].
If memory is not an issue, the algorithm can be speeded up by inputing
the optional array \f[C]plx\f[R] of precomputed associated Legendre
functions on the Gauss-Legendre nodes.
If \f[C]plx\f[R] is not specified, then it is necessary to input the
optional array \f[C]zero\f[R] that contains the latitudinal
Gauss-Legendre quadrature nodes.
.PP
This routine uses geodesy 4-pi normalized spherical harmonics that
exclude the Condon-Shortley phase; This can not be modified.
.SH References
.PP
Wieczorek, M.
A.
and R.
J.
Phillips, Potential anomalies on a sphere: applications to the thickness
of the lunar crust, J.
Geophys.
Res., 103, 1715-1724, 1998.
.SH See also
.PP
batohilm, cilmplus, cilmplusrhoh, shexpandglq, makegridglq, shglq,
shexpanddh, makegriddh, glqgridcoord, downcontfilterma, downcontfiltermc