The Elk Code
roteuler.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: roteuler
8 ! !INTERFACE:
9 pure subroutine roteuler(rot,ang)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! rot : proper rotation matrix (in,real(3,3))
12 ! ang : Euler angles (alpha, beta, gamma) (out,real(3))
13 ! !DESCRIPTION:
14 ! Given a proper $3\times 3$ rotation matrix
15 ! \begin{align*}
16 ! &R(\alpha,\beta,\gamma)=\\
17 ! &\left(\begin{matrix}
18 ! \cos\gamma\cos\beta\cos\alpha-\sin\gamma\sin\alpha &
19 ! \cos\gamma\cos\beta\sin\alpha+\sin\gamma\cos\alpha &
20 ! -\cos\gamma\sin\beta \\
21 ! -\sin\gamma\cos\beta\cos\alpha-\cos\gamma\sin\alpha &
22 ! -\sin\gamma\cos\beta\sin\alpha+\cos\gamma\cos\alpha &
23 ! \sin\gamma\sin\beta \\
24 ! \sin\beta\cos\alpha &
25 ! \sin\beta\sin\alpha &
26 ! \cos\beta
27 ! \end{matrix}\right),
28 ! \end{align*}
29 ! this routine determines the Euler angles, $(\alpha,\beta,\gamma)$. This
30 ! corresponds to the so-called `$y$-convention', which involves the following
31 ! successive rotations of the coordinate system:
32 ! \begin{itemize}
33 ! \item[1.]{The $x_1$-, $x_2$-, $x_3$-axes are rotated anticlockwise through
34 ! an angle $\alpha$ about the $x_3$ axis}
35 ! \item[2.]{The $x_1'$-, $x_2'$-, $x_3'$-axes are rotated anticlockwise
36 ! through an angle $\beta$ about the $x_2'$ axis}
37 ! \item[3.]{The $x_1''$-, $x_2''$-, $x_3''$-axes are rotated anticlockwise
38 ! through an angle $\gamma$ about the $x_3''$ axis}
39 ! \end{itemize}
40 ! Note that the Euler angles are not necessarily unique for a given rotation
41 ! matrix. It is the responsibility of the calling routine to provide a valid
42 ! proper rotation matrix, i.e. $R$ is orthogonal and $|R|=1$, since no
43 ! checking is performed by this routine.
44 !
45 ! !REVISION HISTORY:
46 ! Created May 2003 (JKD)
47 ! Fixed problem thanks to Frank Wagner, June 2013 (JKD)
48 ! Removed check for invalid rotation matrices, March 2026 (JKD)
49 !EOP
50 !BOC
51 implicit none
52 ! arguments
53 real(8), intent(in) :: rot(3,3)
54 real(8), intent(out) :: ang(3)
55 ! local variables
56 real(8), parameter :: eps=1.d-8
57 real(8), parameter :: pi=3.1415926535897932385d0
58 if ((abs(rot(3,1)) > eps).or.(abs(rot(3,2)) > eps)) then
59  ang(1)=atan2(rot(3,2),rot(3,1))
60  if (abs(rot(3,1)) > abs(rot(3,2))) then
61  ang(2)=atan2(rot(3,1)/cos(ang(1)),rot(3,3))
62  else
63  ang(2)=atan2(rot(3,2)/sin(ang(1)),rot(3,3))
64  end if
65  ang(3)=atan2(rot(2,3),-rot(1,3))
66 else
67  ang(1)=atan2(rot(1,2),rot(1,1))
68  if (rot(3,3) > 0.d0) then
69  ang(2)=0.d0
70  ang(3)=0.d0
71  else
72  ang(2)=pi
73  ang(3)=pi
74  end if
75 end if
76 end subroutine
77 !EOC
78 
pure subroutine roteuler(rot, ang)
Definition: roteuler.f90:10
real(8), parameter pi
Definition: modmain.f90:1224