-
Notifications
You must be signed in to change notification settings - Fork 106
/
shexpandglqc.html
94 lines (83 loc) · 7.68 KB
/
shexpandglqc.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
<!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="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../transformations.html" class="dir">Spherical Harmonic Transformations</a></p>
<h1 id="shexpandglqc">SHExpandGLQC</h1>
<p>Expand a 2D grid sampled on the Gauss-Legendre quadrature nodes into spherical harmonics.</p>
<h1 id="usage">Usage</h1>
<p>call SHExpandGLQC (<code>cilm</code>, <code>lmax</code>, <code>gridglq</code>, <code>w</code>, <code>plx</code>, <code>zero</code>, <code>norm</code>, <code>csphase</code>, <code>lmax_calc</code>, <code>exitstatus</code>)</p>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>cilm</code> : output, complex*16, dimension (2, <code>lmax</code>+1, <code>lmax</code>+1) or (2, <code>lmax_calc</code>+1, <code>lmax_calc</code>+1)</dt>
<dd>The complex spherical harmonic coefficients of the complex function. The first index specifies the coefficient corresponding to the positive and negative order of <code>m</code>, respectively, with <code>Clm=cilm(1,l+1,m+1)</code> and <code>Cl,-m =cilm(2,l+1,m+1)</code>.
</dd>
<dt><code>lmax</code> : input, integer</dt>
<dd>The spherical harmonic bandwidth of the grid. If <code>lmax_calc</code> is not specified, this also corresponds to the maximum spherical harmonic degree of the coefficients <code>cilm</code>.
</dd>
<dt><code>gridglq</code> : input, complex*16, dimension(<code>lmax</code>+1, 2*<code>lmax</code>+1)</dt>
<dd>A 2D grid of complex data sampled on the Gauss-Legendre quadrature nodes. 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. See also <code>GLQGridCoord></code>.
</dd>
<dt><code>w</code> : input, real*8, dimension (<code>lmax</code>+1)</dt>
<dd>The Gauss-Legendre quadrature weights used in the integration over latitude. These are obtained from a call to <code>SHGLQ</code>.
</dd>
<dt><code>plx</code> : input, optional, 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 Gauss-Legendre quadrature nodes. These are determined from a call to <code>SHGLQ</code> with the option <code>cnorm=1</code>. Either <code>plx</code> or <code>zero</code> must be present, but not both.
</dd>
<dt><code>zero</code> : input, optional, real*8, 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>. Either <code>plx</code> or <code>zero</code> must be present, but not both.
</dd>
<dt><code>norm</code> : input, optional, integer, default = 1</dt>
<dd>1 (default) = 4-pi (geodesy) normalized harmonics; 2 = Schmidt semi-normalized harmonics; 3 = unnormalized harmonics; 4 = orthonormal harmonics.
</dd>
<dt><code>csphase</code> : input, 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> : input, optional, integer, default = <code>lmax</code></dt>
<dd>The maximum spherical harmonic degree calculated in the spherical harmonic expansion.
</dd>
<dt><code>exitstatus</code> : output, optional, integer</dt>
<dd>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.
</dd>
</dl>
<h1 id="description">Description</h1>
<p><code>SHExpandGLQC</code> will expand a 2-dimensional grid of complex data sampled on the Gauss-Legendre quadrature nodes into complex spherical harmonics. This is the inverse of the routine <code>MakeGridGLQC</code>. The latitudinal nodes of the input grid 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. It is implicitly assumed that the function is bandlimited to degree <code>lmax</code>. If the optional parameter <code>lmax_calc</code> is specified, the spherical harmonic coefficients will be calculated up to this degree, instead of <code>lmax</code>.</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 in this routine 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 spherical harmonic transformation may be speeded up by precomputing the Legendre functions on the Gauss-Legendre quadrature nodes in the routine <code>SHGLQ</code> with the optional parameter <code>cnorm</code> set to 1. 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="makegridglqc.html"><code>makegridglqc</code></a>, <a href="shexpandglq.html"><code>shexpandglq</code></a>, <a href="makegridglq.html"><code>makegridglq</code></a>, <a href="shexpanddh.html"><code>shexpanddh</code></a>, <a href="makegriddh.html"><code>makegriddh</code></a>, <a href="shexpanddhc.html"><code>shexpanddhc</code></a>, <a href="makegriddhc.html"><code>makegriddhc</code></a>, <a href="shexpandlsq.html"><code>shexpandlsq</code></a>, <a href="glqgridcoord.html"><code>glqgridcoord</code></a>, <a href="shglq.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="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../transformations.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>