-
Notifications
You must be signed in to change notification settings - Fork 106
/
shrotatecoef.html
97 lines (86 loc) · 6.36 KB
/
shrotatecoef.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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01//EN"
"http://www.w3.org/TR/html4/strict.dtd">
<html>
<head>
<title>SHTOOLS - Spherical harmonic rotations</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="../../rotations.html" class="dir">Spherical Harmonic Rotations</a></p>
<h1 id="shrotatecoef">SHRotateCoef</h1>
<p>Determine the spherical harmonic coefficients of a complex function rotated by three Euler angles.</p>
<h1 id="usage">Usage</h1>
<p>call SHRotateCoef (<code>x</code>, <code>coef</code>, <code>rcoef</code>, <code>dj</code>, <code>lmax</code>, <code>exitstatus</code>)</p>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>x</code> : input, real*8, dimension(3)</dt>
<dd>The three Euler angles, alpha, beta, and gamma, in radians.
</dd>
<dt><code>coef</code> : input, real*8, dimension (2, (<code>lmax</code>+1)*(<code>lmax</code>+2)/2)</dt>
<dd>The input complex spherical harmonic coefficients. This is an indexed array where the real and complex components are given by <code>coef(1,:)</code> and <code>coef(2,:)</code>, respectively. The functions <code>SHCilmToCindex</code> and <code>SHCindexToCilm</code> are used to convert to and from indexed and <code>cilm(2,:,:)</code> form. The coefficients must correspond to unit-normalized spherical harmonics that possess the Condon-Shortley phase convention.
</dd>
<dt><code>rcoef</code> : output, real*8, dimension (2, (<code>lmax</code>+1)*(<code>lmax</code>+2)/2)</dt>
<dd>The spherical harmonic coefficients of the rotated function in indexed form.
</dd>
<dt><code>dj</code> : input, real*8, dimension (<code>lmax</code>+1, <code>lmax</code>+1, <code>lmax</code>+1)</dt>
<dd>The rotation matrix <code>dj(pi/2)</code> obtained from a call to <code>djpi2</code>.
</dd>
<dt><code>lmax</code> : input, integer</dt>
<dd>The maximum spherical harmonic degree of the input and output coefficients.
</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>SHRotateCoef</code> will take the complex spherical harmonic coefficients of a function, rotate it according to the three Euler anlges in <code>x</code>, and output the spherical harmonic coefficients of the rotated function. The input and output coefficients are in an indexed form that can be converted to and from <code>cilm(2,:,:)</code> form by using the functions <code>SHCilmToCindex</code> and <code>SHCindexToCilm</code>. The coefficients must correspond to unit-normalized spherical harmonics that possess the Condon-Shortley phase convention. Real spherical harmonics can be converted to and from complex form using <code>SHrtoc</code> and <code>SHctor</code>. The input rotation matrix <code>dj</code> is computed by a call to <code>djpi2</code>.</p>
<p>The rotation of a coordinate system or body can be viewed in two complementary ways involving three successive rotations. Both methods have the same initial and final configurations, and the angles listed in both schemes are the same.</p>
<p><code>Scheme A:</code></p>
<ol style="list-style-type: upper-roman">
<li>Rotation about the z axis by alpha.</li>
<li>Rotation about the new y axis by beta.</li>
<li>Rotation about the new z axis by gamma.</li>
</ol>
<p><code>Scheme B:</code></p>
<ol style="list-style-type: upper-roman">
<li>Rotation about the z axis by gamma.</li>
<li>Rotation about the initial y axis by beta.</li>
<li>Rotation about the initial z axis by alpha.</li>
</ol>
<p>The rotations can further be viewed either as a rotation of the coordinate system or the physical body. For a rotation of the coordinate system without rotation of the physical body, use</p>
<p><code>x(alpha, beta, gamma)</code>.</p>
<p>For a rotation of the physical body without rotation of the coordinate system, use</p>
<p><code>x(-gamma, -beta, -alpha)</code>.</p>
<p>To perform the inverse transform of <code>x(alpha, beta, gamma)</code>, use <code>x(-gamma, -beta, -alpha)</code>.</p>
<p>Note that this routine uses the "y convention", where the second rotation is with respect to the new y axis. If alpha, beta, and gamma were orginally defined in terms of the "x convention", where the second rotation was with respect to the newx axis, the Euler angles according to the y convention would be <code>alpha_y=alpha_x-pi/2</code>, <code>beta_x=beta_y</code>, and <code>gamma_y=gamma_x+pi/2</code>.</p>
<h1 id="see-also">See also</h1>
<p><a href="djpi2.html">djpi2</a>, <a href="shrotaterealcoef.html">shrotaterealcoef</a>, <a href="shctor.html">shctor</a>, <a href="shrtoc.html">shrtoc</a>, <a href="shcilmtocindex.html">shcilmtocindex</a>, <a href="shcindextocilm.html">shcindextocilm</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="../../rotations.html" class="dir">Spherical Harmonic Rotations</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>