/
cilmminusrhoh.html
105 lines (94 loc) · 9.93 KB
/
cilmminusrhoh.html
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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01//EN"
"http://www.w3.org/TR/html4/strict.dtd">
<html>
<head>
<title>SHTOOLS - Gravity and magnetics routines</title>
<meta http-equiv="Content-Type" content="text/html; charset=UTF-8">
<link rel="stylesheet" type="text/css" href="../../css/sh.css">
<link rel="icon" type="image/vnd.microsoft.icon" href="../../images/favicon.ico">
</head>
<body>
<div class="main">
<p class="centeredimage"><img src="../../images/logo.jpg" width=894 height=135 alt="SHTOOLS --- Tools for working with spherical harmonics"></p>
<table class="menu">
<tbody>
<tr>
<td><a href="http://shtools.ipgp.fr/">HOME</a></td>
<td><a href="https://github.com/SHTOOLS/SHTOOLS/releases">DOWNLOAD</a></td>
<td class="selected"><a href="../../documentation.html">DOCUMENTATION</a></td>
<td><a href="../../faq.html">FAQ</a> </td>
</tr>
</tbody>
</table>
<p class="dir">
> <a href="../../../index.html" class="dir">Home</a> > <a href="../../documentation.html" class="dir">Documentation</a> > <a href="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../gravmag.html" class="dir">Gravity and Magnetics</a></p>
<h1 id="cilmminusrhoh">CilmMinusRhoH</h1>
<p>Calculate the gravitational potential interior to relief referenced to a spherical interface with laterally varying density using the finite amplitude algorithm of Wieczorek (2007).</p>
<h1 id="usage">Usage</h1>
<p>call CilmMinusRhoH (<code>cilm</code>, <code>gridin</code>, <code>lmax</code>, <code>nmax</code>, <code>mass</code>, <code>d</code>, <code>rho</code>, <code>gridtype</code>, <code>w</code>, <code>zero</code>, <code>plx</code>, <code>n</code>, <code>dref</code>)</p>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>cilm</code> : output, real*8, dimension (2, <code>lmax</code>+1, <code>lmax</code>+1)</dt>
<dd>The real spherical harmonic coefficients (geodesy normalized) of the gravitational potential corresponding to constant density relief referenced to a spherical interface of radius <code>d</code>.
</dd>
<dt><code>gridin</code> : input, real*8, dimension (<code>lmax</code>+1, 2*<code>lmax</code>+1) for gridtype 1, (<code>n</code>, <code>n</code>) for gridtype 2, (<code>n</code>, 2*<code>n</code>) for gridtype 3</dt>
<dd>The radii of the interface evaluated on a grid corresponding to a function of maximum spherical harmonic degree <code>lmax</code>. This is calculated by a call to either <code>MakeGridGLQ</code> or <code>MakeGridDH</code>.
</dd>
<dt><code>lmax</code> : input, integer</dt>
<dd>The maximum spherical harmonic degree of the output spherical harmonic coefficients. This degree also determines the dimension of the input relief <code>gridin</code> for <code>gridtype</code> 1. For Driscoll-Healy grids, <code>lmax</code> must be less than or equal to <code>n/2-1</code>.
</dd>
<dt><code>nmax</code> : input, integer</dt>
<dd>The maximum order used in the Taylor-series expansion used in calculating the potential coefficients.
</dd>
<dt><code>mass</code> : input, real*8</dt>
<dd>The mass of the planet in kg.
</dd>
<dt><code>d</code> : output, real*8</dt>
<dd>The mean radius of the relief in meters.
</dd>
<dt><code>rho</code> : input, real*8, dimension (<code>lmax</code>+1, 2*<code>lmax</code>+1) for gridtype 1, (<code>n</code>, <code>n</code>) for gridtype 2, (<code>n</code>, 2*<code>n</code>) for gridtype 3</dt>
<dd>The density contrast of the relief in kg/m^3 expressed on a grid with the same dimensions as <code>gridin</code>.
</dd>
<dt><code>gridtype</code> : input, integer</dt>
<dd>1 = Gauss-Legendre grids, calculated using <code>SHGLQ</code> and <code>MakeGridGLQ></code>. 2 = Equally sampled Driscoll-Healy grids, <code>n</code> by <code>n</code>, calculated using <code>MakeGridDH</code>. 3 = Equally spaced Driscoll-Healy grids, <code>n</code> by 2<code>n</code>, calculated using <code>MakeGridDH</code>.
</dd>
<dt><code>w</code> : optional, input, real*8, dimension (<code>lmax</code>+1)</dt>
<dd>The weights used in the Gauss-Legendre quadrature, which are required for <code>gridtype</code> = 1. These are calculated from a call to <code>SHGLQ</code>.
</dd>
<dt><code>zero</code> : optional, input, real*8, dimension (<code>lmax</code>+1)</dt>
<dd>The nodes used in the Gauss-Legendre quadrature over latitude for <code>gridtype</code> 1, calculated by a call to <code>SHGLQ</code>. One of <code>plx</code> or <code>zero</code> must be present when <code>gridtype=1</code>, but not both.
</dd>
<dt><code>plx</code> : optional, input, real*8, dimension (<code>lmax</code>+1, (<code>lmax</code>+1)*(<code>lmax</code>+2)/2)</dt>
<dd>An array of the associated Legendre functions calculated at the nodes used in the Gauss-Legendre quadrature for <code>gridtype</code> 1. These are determined from a call to <code>SHGLQ</code>. One of <code>plx</code> or <code>zero</code> must be present when <code>gridtype=1</code>, but not both.
</dd>
<dt><code>n</code> : optional, input, integer</dt>
<dd>The number of samples in latitude when using Driscoll-Healy grids. For a function bandlimited to <code>lmax</code>, <code>n=2(lmax+1)</code>. This is required for gridtypes 2 and 3.
</dd>
<dt><code>dref</code> : optional, input, real*8</dt>
<dd>The reference radius to be used when calculating both the relief and spherical harmonic coefficients. If this is not specified, this parameter will be set equal to the mean radius of <code>gridin</code>.
</dd>
</dl>
<h1 id="description">Description</h1>
<p><code>CilmMinusRhoH</code> will calculate the spherical harmonic coefficients of the gravitational potential interior to relief referenced to a spherical interface with laterally varying density. This is a generalization of equation 30 of Wieczorek (2007), and equation 12 of Wieczorek and Phillips (1998). The potential is strictly valid only when the coefficients are evaluated at a radius greater than the maximum radius of the relief. The relief and laterally varying density are input as a grid, whose type is specified by <code>gridtype</code> (1 for Gauss-Legendre quadrature grids, 2 for <code>n</code> by <code>n</code> Driscoll and Healy sampled grids, and 3 for <code>n</code> by 2<code>n</code> Driscoll and Healy sampled grids). The input relief <code>gridin</code> must correspond to absolute radii. The parameter <code>nmax</code> is the order of the Taylor series used in the algorithm to approximate the potential coefficients. By default, the relief and spherical harmonic coefficients will be referenced to the mean radius of <code>gridin</code>. However, if the optional parameter <code>dref</code> is specified, this will be used instead as the reference radius.</p>
<p>It is important to understand that as an intermediate step, this routine calculates the spherical harmonic coefficients of the density multiplied by the relief (referenced to the mean radius of <code>gridin</code> or <code>dref</code>) raised to the nth power. As such, if the input function is bandlimited to degree <code>L</code>, the resulting function will be bandlimited to degree <code>L*nmax</code>. This subroutine implicitly assumes that <code>gridin</code> and <code>rho</code> have an effective spherical harmonic bandwidth greater or equal to this value. (The effective bandwidth is equal to <code>lmax</code> for <code>gridtype</code> 1, and is <code>n/2-1</code> for <code>gridtype</code> 2 or 3.) If this is not the case, aliasing will occur. In practice, for accurate results, it is found that the effective bandwidth needs only to be about three times the size of <code>L</code>, though this should be verified for each application. Thus, if the input function is considered to be bandlimited to degree <code>L</code>, the function should be evaluated on a grid corresponding to a maximum degree of about <code>3*</code>L. Aliasing effects can be partially mitigated by using Driscoll and Healy <code>n</code> by 2<code>n</code> grids.</p>
<p>If the input grid is evaluated on the Gauss-Legendre points, it is necessary to specify the optional parameters <code>w</code> and <code>zero</code>, or <code>w</code> and <code>plx</code>, which are calculated by a call to <code>SHGLQ</code>. In contast, if Driscoll-Healy grids are used, it is necessary to specify the optional parameter <code>n</code>. If memory is not an issue, the algorithm can be speeded up when using Gauss-Lengendre grids by inputing the optional array <code>plx</code> (along with <code>w</code>) of precomputed associated Legendre functions on the Gauss-Legendre nodes. Both of these variables are computed by a call to <code>SHGLQ</code>.</p>
<p>This routine uses geodesy 4-pi normalized spherical harmonics that exclude the Condon-Shortley phase.</p>
<h1 id="references">References</h1>
<p>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.</p>
<p>Wieczorek, M. A., Gravity and topography of the terrestrial planets, Treatise on Geophysics, 10, 165-206, doi:10.1016/B978-044452748-6/00156-5, 2007.</p>
<h1 id="see-also">See also</h1>
<p><a href="cilmminus.html">cilmminus</a>, <a href="cilmplus.html">cilmplus</a>, <a href="cilmplusrhoh.html">cilmplusrhoh</a>, <a href="makegridglq.html">makegridglq</a>, <a href="makegriddh.html">makegriddh</a>, <a href="shglq.html">shglq</a>, <a href="glqgridcoord.html">glqgridcoord</a></p>
<p class="dir">
> <a href="../../../index.html" class="dir">Home</a> > <a href="../../documentation.html" class="dir">Documentation</a> > <a href="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../gravmag.html" class="dir">Gravity and Magnetics</a></p>
<table class="footer2" summary = "SHTOOLS; Fortran and Python spherical harmonic transform software package">
<tbody>
<tr>
<td class="c1"><a href="http://www.ipgp.fr/">Institut de Physique du Globe de Paris</a></td>
<td class="c2"><a href="http://www.sorbonne-paris-cite.fr/index.php/en">University of Sorbonne Paris Cité</a></td>
<td class="c3">© 2016 <a href="https://github.com/SHTOOLS">SHTOOLS</a></td>
</tr>
</tbody>
</table>
</div>
</body>
</html>