/
pymakegravgriddh.html
114 lines (103 loc) · 9.01 KB
/
pymakegravgriddh.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
106
107
108
109
110
111
112
113
114
<!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="../../python-routines.html" class="dir">Python</a> > <a href="../../pygravmag.html" class="dir">Gravity and Magnetics</a></p>
<h1 id="makegravgriddh">MakeGravGridDH</h1>
<p>Create 2D cylindrical maps on a flattened and rotating ellipsoid of all three components of the gravity field, the gravity disturbance, and the gravitational potential.</p>
<h1 id="usage">Usage</h1>
<p><code>rad</code>, <code>theta</code>, <code>phi</code>, <code>total</code> = pyshtools.MakeGravGridDH (<code>cilm</code>, <code>gm</code>, <code>r0</code>, [<code>a</code>, <code>f</code>, <code>lmax</code>, <code>sampling</code>, <code>lmax_calc</code>, <code>omega</code>, <code>normal_gravity</code>])</p>
<h1 id="returns">Returns</h1>
<dl>
<dt><code>rad</code> : float, dimension (2*<code>lmax</code>+2, <code>sampling</code>*(2*<code>lmax</code>+2))</dt>
<dd>A 2D equally sampled (<code>n</code> by <code>n</code>) or equally spaced (<code>n</code> by 2<code>n</code>) grid of the radial component of the gravity field corresponding to the input spherical harmonic coefficients <code>cilm</code>. The first latitudinal band corresponds to 90 N, the latitudinal band for 90 S is not included, and the latitudinal sampling interval is 180/<code>n</code> degrees. The first longitudinal band is 0 E, the longitudinal band for 360 E is not included, and the longitudinal sampling interval is 360/<code>n</code> for an equally sampled and 180/<code>n</code> for an equally spaced grid, respectively.
</dd>
<dt><code>theta</code> : float, dimension (2*<code>lmax</code>+2, <code>sampling</code>*(2*<code>lmax</code>+2))</dt>
<dd>A 2D equally sampled or equally spaced grid of the theta component of the gravity field.
</dd>
<dt><code>phi</code> : float, dimension (2*<code>lmax</code>+2, <code>sampling</code>*(2*<code>lmax</code>+2))</dt>
<dd>A 2D equally sampled or equally spaced grid of the phi component of the gravity field.
</dd>
<dt><code>total</code> : float, dimension (2*<code>lmax</code>+2, <code>sampling</code>*(2*<code>lmax</code>+2))</dt>
<dd>A 2D equally sampled or equally spaced grid of the magnitude of the gravity acceleration.
</dd>
</dl>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>cilm</code> : float, dimension (2, <code>lmaxin</code>+1, <code>lmaxin</code>+1)</dt>
<dd>The real gravitational potential spherical harmonic coefficients. The coefficients c1lm and c2lm refer to the cosine and sine coefficients, respectively, with <code>c1lm=cilm[0,l,m]</code> and <code>c2lm=cilm[1,l,m]</code>.
</dd>
<dt><code>gm</code> : float</dt>
<dd>The gravitational constant multiplied by the mass of the planet.
</dd>
<dt><code>r0</code>: float</dt>
<dd>The reference radius of the spherical harmonic coefficients.
</dd>
<dt><code>a</code> : optional, float, default = <code>r0</code></dt>
<dd>The semi-major axis of the flattened ellipsoid on which the field is computed.
</dd>
<dt><code>f</code> : optional, float, default = 0</dt>
<dd>The flattening of the reference ellipsoid: <code>f=(R_equator-R_pole)/R_equator</code>.
</dd>
<dt><code>lmax</code> : optional, integer, default = <code>lmaxin</code></dt>
<dd>The maximum spherical harmonic degree of the coefficients <code>cilm</code>. This determines the number of samples of the output grids, <code>n=2lmax+2</code>, and the latitudinal sampling interval, <code>90/(lmax+1)</code>.
</dd>
<dt><code>sampling</code> : optional, integer, default = 2</dt>
<dd>If 1 (default) the output grids are equally sampled (<code>n</code> by <code>n</code>). If 2, the grids are equally spaced (<code>n</code> by 2<code>n</code>).
</dd>
<dt><code>lmax_calc</code> : optional, integer, default = <code>lmax</code></dt>
<dd>The maximum spherical harmonic degree used in evaluating the functions. This must be less than or equal to <code>lmax</code>.
</dd>
<dt><code>omega</code> : optional, float, default = 0</dt>
<dd>The angular rotation rate of the planet.
</dd>
<dt><code>normal_gravity</code> : optional, integer, default = 1</dt>
<dd>If 1, the normal gravity (the gravitational acceleration on the ellipsoid) will be subtracted from the total gravity, yielding the "gravity disturbance." This is done using Somigliana's formula (after converting geocentric to geodetic coordinates).
</dd>
</dl>
<h1 id="description">Description</h1>
<p><code>MakeGravGridDH</code> will create 2-dimensional cylindrical maps from the spherical harmonic coefficients <code>cilm</code>, equally sampled (<code>n</code> by <code>n</code>) or equally spaced (<code>n</code> by 2<code>n</code>) in latitude and longitude, for the three vector components of the gravity field, the magnitude of the gravity field, and the potential (all using geocentric coordinates). The gravitational potential is given by</p>
<p><code>V = GM/r Sum_{l=0}^lmax (r0/r)^l Sum_{m=-l}^l C_{lm} Y_{lm}</code>,</p>
<p>and the gravitational acceleration is</p>
<p><code>B = Grad V</code>.</p>
<p>The coefficients are referenced to a radius <code>r0</code>, and the function is computed on a flattened ellipsoid with semi-major axis <code>a</code> (i.e., the mean equatorial radius) and flattening <code>f</code>. All grids are output in SI units, and the sign of the radial components is positive when directed upwards. If the optional angular rotation rate <code>omega</code> is specified, the potential and radial gravitational acceleration will be calculated in a body-fixed rotating reference frame.</p>
<p>To remove the "normal gravity" (the total gravitational acceleration on the ellipsoid) from the magnitude of the total gravity field (to obtain the "gravity disturbance"), set <code>normal_gravity</code> to 1. To convert m/s^2 to mGals, multiply the gravity grids by 10^5.</p>
<p>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). The components of the gravity vector are calculated using spherical coordinates whose origin is at the center of the planet, and the components are thus not normal to the reference ellipsoid.</p>
<p>The default is to calculate grids for use in the Driscoll and Healy routines that are equally spaced (<code>n</code> by <code>2n</code>), but this can be changed to calculate equally sampled grids (<code>n</code> by <code>n</code>) by setting the optional argument <code>sampling</code> to 1. The input value of <code>lmax</code> determines the number of samples, <code>n=2lmax+2</code>, and the latitudinal sampling interval, 90/(<code>lmax</code>+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/<code>n</code> degrees. The first longitudinal band is 0 E, the longitudinal band for 360 E is not calculated, and the longitudinal sampling interval is 360/<code>n</code> for equally sampled and 180/<code>n</code> for equally spaced grids, respectively.</p>
<h1 id="references">References</h1>
<p>Driscoll, J.R. and D.M. Healy, Computing Fourier transforms and convolutions on the 2-sphere, Adv. Appl. Math., 15, 202-250, 1994.</p>
<h1 id="see-also">See also</h1>
<p><a href="pymakegeoidgrid.html">makegeoidgrid</a>, <a href="pymakegravgradgriddh.html">makegravgradgriddh</a>, <a href="pynormalgravity.html">normalgravity</a>, <a href="pymakegriddh.html">makegriddh</a></p>
<p class="dir">
> <a href="../../../index.html" class="dir">Home</a> > <a href="../../documentation.html" class="dir">Documentation</a> > <a href="../../python-routines.html" class="dir">Python</a> > <a href="../../pygravmag.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>