RV2M Subroutine

public subroutine RV2M(w, r)

Form the r-matrix corresponding to a given r-vector.

Status: vector/matrix support routine.

Notes

  1. A rotation matrix describes a rotation through some angle about some arbitrary axis called the Euler axis. The "rotation vector" supplied to this routine has the same direction as the Euler axis, and its magnitude is the angle in radians.

  2. If W is null, the unit matrix is returned.

  3. The reference frame rotates clockwise as seen looking along the rotation vector from the origin.

History

  • IAU SOFA revision: 2015 January 30

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in), dimension(3):: w

rotation vector (Note 1)

real(kind=wp), intent(out), dimension(3,3):: r

rotation matrix


Called by

proc~~rv2m~~CalledByGraph proc~rv2m RV2M proc~fk5hip FK5HIP proc~fk5hip->proc~rv2m proc~hfk5z HFK5Z proc~hfk5z->proc~rv2m proc~hfk5z->proc~fk5hip proc~fk5hz FK5HZ proc~fk5hz->proc~rv2m proc~fk5hz->proc~fk5hip proc~fk52h FK52H proc~fk52h->proc~fk5hip proc~h2fk5 H2FK5 proc~h2fk5->proc~fk5hip

Contents

Source Code


Source Code

    subroutine RV2M ( w, r )

    implicit none

    real(wp),dimension(3),intent(in) :: w !! rotation vector (Note 1)
    real(wp),dimension(3,3),intent(out) :: r !! rotation matrix

    real(wp) :: x, y, z, phi, s, c, f

    !  Euler angle (magnitude of rotation vector) and functions.
    x = w(1)
    y = w(2)
    z = w(3)
    phi = sqrt(x*x + y*y + z*z)
    s = sin(phi)
    c = cos(phi)
    f = 1.0_wp - c

    !  Euler axis (direction of rotation vector), perhaps null.
    if ( phi > 0.0_wp ) then
       x = x / phi
       y = y / phi
       z = z / phi
    end if

    !  Form the rotation matrix.
    r(1,1) = x*x*f + c
    r(1,2) = x*y*f + z*s
    r(1,3) = x*z*f - y*s
    r(2,1) = y*x*f - z*s
    r(2,2) = y*y*f + c
    r(2,3) = y*z*f + x*s
    r(3,1) = z*x*f + y*s
    r(3,2) = z*y*f - x*s
    r(3,3) = z*z*f + c

    end subroutine RV2M