/
shrotaterealcoef.1
118 lines (118 loc) · 4.01 KB
/
shrotaterealcoef.1
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
115
116
117
118
.\" Automatically generated by Pandoc 2.0.3
.\"
.TH "shrotaterealcoef" "1" "2016\-12\-15" "Fortran 95" "SHTOOLS 4.1"
.hy
.SH SHRotateRealCoef
.PP
Determine the spherical harmonic coefficients of a real function rotated
by three Euler angles.
.SH Usage
.PP
call SHRotateRealCoef (\f[C]cilmrot\f[], \f[C]cilm\f[], \f[C]lmax\f[],
\f[C]x\f[], \f[C]dj\f[], \f[C]exitstatus\f[])
.SH Parameters
.TP
.B \f[C]cilmrot\f[] : output, real*8, dimension (2, \f[C]lmax\f[]+1, \f[C]lmax\f[]+1)
The spherical harmonic coefficients of the rotated function, normalized
for use with the geodesy 4\-pi spherical harmonics.
.RS
.RE
.TP
.B \f[C]cilm\f[] : input, real*8, dimension (2, \f[C]lmax\f[]+1, \f[C]lmax\f[]+1)
The input real spherical harmonic coefficients.
The coefficients must correspond to geodesy 4\-pi normalized spherical
harmonics that do not possess the Condon\-Shortley phase convention.
.RS
.RE
.TP
.B \f[C]x\f[] : input, real*8, dimension(3)
The three Euler angles, alpha, beta, and gamma, in radians.
.RS
.RE
.TP
.B \f[C]dj\f[] : input, real*8, dimension (\f[C]lmax\f[]+1, \f[C]lmax\f[]+1, \f[C]lmax\f[]+1)
The rotation matrix \f[C]dj(pi/2)\f[], obtained from a call to
\f[C]djpi2\f[].
.RS
.RE
.TP
.B \f[C]lmax\f[] : input, integer
The maximum spherical harmonic degree of the input and output
coefficients.
.RS
.RE
.TP
.B \f[C]exitstatus\f[] : output, optional, integer
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.
.RS
.RE
.SH Description
.PP
\f[C]SHRotateRealCoef\f[] will take the real spherical harmonic
coefficients of a function, rotate it according to the three Euler
anlges in \f[C]x\f[], and output the spherical harmonic coefficients of
the rotated function.
The input and output coefficients must correspond to geodesy 4\-pi
normalized spherical harmonics that do not possess the Condon\-Shortley
phase convention.
The input rotation matrix \f[C]dj\f[] is computed by a call to
\f[C]djpi2\f[].
.PP
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.
.PP
\f[C]Scheme\ A:\f[]
.IP " (I)" 6
Rotation about the z axis by alpha.
.IP " (II)" 6
Rotation about the new y axis by beta.
.IP "(III)" 6
Rotation about the new z axis by gamma.
.PP
\f[C]Scheme\ B:\f[]
.IP " (I)" 6
Rotation about the z axis by gamma.
.IP " (II)" 6
Rotation about the initial y axis by beta.
.IP "(III)" 6
Rotation about the initial z axis by alpha.
.PP
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
.PP
\f[C]x(alpha,\ beta,\ gamma)\f[].
.PP
For a rotation of the physical body without rotation of the coordinate
system, use
.PP
\f[C]x(\-gamma,\ \-beta,\ \-alpha)\f[].
.PP
To perform the inverse transform of \f[C]x(alpha,\ beta,\ gamma)\f[],
use \f[C]x(\-gamma,\ \-beta,\ \-alpha)\f[].
.PP
Note that this routine uses the \[lq]y convention\[rq], where the second
rotation is with respect to the new y axis.
If alpha, beta, and gamma were orginally defined in terms of the \[lq]x
convention\[rq], where the second rotation was with respect to the newx
axis, the Euler angles according to the y convention would be
\f[C]alpha_y=alpha_x\-pi/2\f[], \f[C]beta_x=beta_y\f[], and
\f[C]gamma_y=gamma_x+pi/2\f[].
.PP
This routine first converts the real coefficients to complex form using
\f[C]SHrtoc\f[].
Then the coefficients are converted to indexed form using
\f[C]SHCilmToCindex\f[], these are sent to \f[C]SHRotateCoef\f[], the
result if converted back to \f[C]cilm\f[] complex form using
\f[C]SHCindexToCilm\f[], and these are finally converted back to real
form using \f[C]SHctor\f[].
.SH See also
.PP
djpi2, shrotatecoef, shctor, shrtoc, shcilmtocindex, shcindextocilm