-
Notifications
You must be signed in to change notification settings - Fork 103
/
shmultiply.html
88 lines (77 loc) · 5.14 KB
/
shmultiply.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="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../transformations.html" class="dir">Spherical Harmonic Transformations</a></p>
<h1 id="shmultiply">SHMultiply</h1>
<p>Multiply two functions and determine the spherical harmonic coefficients of the resulting function.</p>
<h1 id="usage">Usage</h1>
<p>call SHMultiply (<code>shout</code>, <code>sh1</code>, <code>lmax1</code>, <code>sh2</code>, <code>lmax2</code>, <code>precomp</code>, <code>norm</code>, <code>csphase</code>, <code>exitstatus</code>)</p>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>shout</code> : output, real*8, dimension (2, <code>lmax1</code>+<code>lmax2</code>+1, <code>lmax1</code>+<code>lmax2</code>+1)</dt>
<dd>The real spherical harmonic coefficients corresponding to the multiplication of <code>sh1</code> and <code>sh2</code> in the space domain.
</dd>
<dt><code>sh1</code> : input, real*8, dimension (2, <code>lmax1</code>+1, <code>lmax1</code>+1)</dt>
<dd>The spherical harmonic coefficients of the first function.
</dd>
<dt><code>lmax1</code> : input, integer</dt>
<dd>The maximum spherical harmonic degree used in evaluting <code>sh1</code>.
</dd>
<dt><code>sh2</code> : input, real*8, dimension (2, <code>lmax2</code>+1, <code>lmax2</code>+1)</dt>
<dd>The spherical harmonic coefficients of the second function.
</dd>
<dt><code>lmax2</code> : input, integer</dt>
<dd>The maximum spherical harmonic degree used in evaluting <code>sh2</code>.
</dd>
<dt><code>precomp</code> : input, optional, integer, default = 0</dt>
<dd>If 1, the array of Legendre functions <code>plx</code> will be precomputed on the Gauss-Legendre quadrature nodes.
</dd>
<dt><code>norm</code> : input, 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> : 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>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>SHMultiply</code> will take two sets of spherical harmonic coefficients, multiply the functions in the space domain, and expand the resulting field in spherical harmonics using <code>SHExpandGLQ</code>. The spherical harmonic bandwidth of the resulting field is <code>lmax1+lmax2</code>, where <code>lmax1</code> and <code>lmax2</code> are the bandwidths of the input fields. If the optional parameter <code>precomp</code> is set, then the array of Legendre functions <code>plx</code> will be precomputed on the Gauss-Legendre quadrature nodes.</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.</p>
<h1 id="see-also">See also</h1>
<p><a href="shexpandglq.html">shexpandglq</a>, <a href="makegridglq.html">makegridglq</a>, <a href="shglq.html">shglq</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">© 2017 <a href="https://github.com/SHTOOLS">SHTOOLS</a></td>
</tr>
</tbody>
</table>
</div>
</body>
</html>