/
pymakegridglq.html
88 lines (77 loc) · 6.41 KB
/
pymakegridglq.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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01//EN"
"http://www.w3.org/TR/html4/strict.dtd">
<html>
<head>
<title>SHTOOLS - Spherical harmonic transformations</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="../../../index.html">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="../../pytransformations.html" class="dir">Spherical Harmonic Transformations</a></p>
<h1 id="makegridglq">MakeGridGLQ</h1>
<p>Create a 2D map from a set of spherical harmonic coefficients sampled on the Gauss-Legendre quadrature nodes.</p>
<h1 id="usage">Usage</h1>
<p><code>gridglq</code> = MakeGridGLQ (<code>cilm</code>, <code>zero</code>, [<code>lmax</code>, <code>norm</code>, <code>csphase</code>, <code>lmax_calc</code>])</p>
<h1 id="returns">Returns</h1>
<dl>
<dt><code>gridglq</code> : float, dimension (<code>lmax</code>+1, 2*<code>lmax</code>+1)</dt>
<dd>A 2D map of the function sampled on the Gauss-Legendre quadrature nodes.
</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 spherical harmonic coefficients of the function. When evaluating the function, the maximum spherical harmonic degree considered is the minimum of <code>lmax</code>, <code>lmaxin</code> or <code>lmax_calc</code> (if specified). The coefficients <code>C1lm</code> and <code>C2lm</code> refer to the "cosine" (<code>Clm</code>) and "sine" (<code>Slm</code>) coefficients, respectively, with <code>Clm=cilm[0,l,m]</code> and <code>Slm=cilm[1,l,m]</code>.
</dd>
<dt><code>zero</code> : float, dimension (<code>lmax</code>+1)</dt>
<dd>The nodes used in the Gauss-Legendre quadrature over latitude, calculated by a call to <code>SHGLQ</code>.
</dd>
<dt><code>lmax</code> : optional, integer, default = <code>lamxin</code></dt>
<dd>The maximum spherical harmonic bandwidth of the function. This determines the sampling nodes of the dimensions of the output grid.
</dd>
<dt><code>norm</code> : optional, integer, default = 1</dt>
<dd>1 (default) = Geodesy 4-pi normalized harmonics; 2 = Schmidt semi-normalized harmonics; 3 = unnormalized harmonics; 4 = orthonormal harmonics.
</dd>
<dt><code>csphase</code> : optional, integer, default = 1</dt>
<dd>1 (default) = do not apply the Condon-Shortley phase factor to the associated Legendre functions; -1 = append the Condon-Shortley phase factor of (-1)^m to the associated Legendre functions.
</dd>
<dt><code>lmax_calc</code> : optional, integer, default = <code>lmax</code></dt>
<dd>The maximum spherical harmonic degree used in evaluating the function. This must be less than or equal to <code>lmax</code>.
</dd>
</dl>
<h1 id="description">Description</h1>
<p><code>MakeGridGLQ</code> will create a 2-dimensional map from a set of input spherical harmonic coefficients sampled on the Gauss-Legendre quadrature nodes. This is the inverse of the routine <code>SHExpandGLQ</code>. The latitudinal nodes correspond to the zeros of the Legendre polynomial of degree <code>lmax+1</code>, and the longitudinal nodes are equally spaced with an interval of <code>360/(2*lmax+1)</code> degrees. When evaluating the function, the maximum spherical harmonic degree that is considered is the minimum of <code>lmax</code>, the size of <code>cilm-1</code>, or <code>lmax_calc</code> (if specified).</p>
<p>The employed spherical harmonic normalization and Condon-Shortley phase convention can be set by the optional arguments <code>norm</code> and <code>csphase</code>; if not set, the default is to use geodesy 4-pi normalized harmonics that exclude the Condon-Shortley phase of (-1)^m. The normalized legendre functions are calculated using the scaling algorithm of Holmes and Featherstone (2002), which are accurate to about degree 2800. The unnormalized functions are accurate only to about degree 15.</p>
<p>The reconstruction of the spherical harmonic function may be speeded up by precomputing the Legendre functions on the Gauss-Legendre quadrature nodes in the routine <code>SHGLQ</code>. However, given that this array contains on the order of <code>lmax</code>**3 entries, this is only feasible for moderate values of <code>lmax</code>.</p>
<h1 id="references">References</h1>
<p>Holmes, S. A., and W. E. Featherstone, A unified approach to the Clenshaw summation and the recursive computation of very high degree and order normalised associated Legendre functions, J. Geodesy, 76, 279- 299, 2002.</p>
<h1 id="see-also">See also</h1>
<p><a href="pyshexpandglq.html"><code>shexpandglq</code></a>, <a href="pyshexpandglqc.html"><code>shexpandglqc</code></a>, <a href="pymakegridglqc.html"><code>makegridglqc</code></a>, <a href="pyshexpanddh.html"><code>shexpanddh</code></a>, <a href="pymakegriddh.html"><code>makegriddh</code></a>, <a href="pyshexpanddhc.html"><code>shexpanddhc</code></a>, <a href="pymakegriddhc.html"><code>makegriddhc</code></a>, <a href="pyshexpandlsq.html"><code>shexpandlsq</code></a>, <a href="pyglqgridcoord.html"><code>glqgridcoord</code></a>, <a href="pyshglq.html"><code>shglq</code></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="../../pytransformations.html" class="dir">Spherical Harmonic Transformations</a></p>
<table class="footer2">
<tbody>
<tr>
<td class="c1"><a href="https://lagrange.oca.eu/?lang=en">Laboratoire Lagrange</a></td>
<td class="c2"><a href="https://www.oca.eu/?lang=en">Observatoire de la Côte d'Azur</a></td>
<td class="c3">© 2016 <a href="https://github.com/SHTOOLS">SHTOOLS</a></td>
</tr>
</tbody>
</table>
</div>
</body>
</html>