aboutsummaryrefslogtreecommitdiff
path: root/src/slalib/dav2m.f
blob: dc3f46caf90e196bb05ea5284c929bc087fe8609 (plain) (blame)
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
SUBROUTINE sla_DAV2M (AXVEC, RMAT)
*+
*     - - - - - -
*      D A V 2 M
*     - - - - - -
*
*  Form the rotation matrix corresponding to a given axial vector.
*  (double precision)
*
*  A rotation matrix describes a rotation about some arbitrary axis.
*  The axis is called the Euler axis, and the angle through which the
*  reference frame rotates is called the Euler angle.  The axial
*  vector supplied to this routine has the same direction as the
*  Euler axis, and its magnitude is the Euler angle in radians.
*
*  Given:
*    AXVEC  d(3)     axial vector (radians)
*
*  Returned:
*    RMAT   d(3,3)   rotation matrix
*
*  If AXVEC is null, the unit matrix is returned.
*
*  The reference frame rotates clockwise as seen looking along
*  the axial vector from the origin.
*
*  P.T.Wallace   Starlink   June 1989
*
*  Copyright (C) 1995 Rutherford Appleton Laboratory
*-

      IMPLICIT NONE

      DOUBLE PRECISION AXVEC(3),RMAT(3,3)

      DOUBLE PRECISION X,Y,Z,PHI,S,C,W



*  Euler angle - magnitude of axial vector - and functions
      X = AXVEC(1)
      Y = AXVEC(2)
      Z = AXVEC(3)
      PHI = SQRT(X*X+Y*Y+Z*Z)
      S = SIN(PHI)
      C = COS(PHI)
      W = 1D0-C

*  Euler axis - direction of axial vector (perhaps null)
      IF (PHI.NE.0D0) THEN
         X = X/PHI
         Y = Y/PHI
         Z = Z/PHI
      END IF

*  Compute the rotation matrix
      RMAT(1,1) = X*X*W+C
      RMAT(1,2) = X*Y*W+Z*S
      RMAT(1,3) = X*Z*W-Y*S
      RMAT(2,1) = X*Y*W-Z*S
      RMAT(2,2) = Y*Y*W+C
      RMAT(2,3) = Y*Z*W+X*S
      RMAT(3,1) = X*Z*W+Y*S
      RMAT(3,2) = Y*Z*W-X*S
      RMAT(3,3) = Z*Z*W+C

      END