FW2M Subroutine

public subroutine FW2M(gamb, phib, psi, eps, r)

Form rotation matrix given the Fukushima-Williams angles.

Status: support routine.

Notes

  1. Naming the following points:

       e = J2000.0 ecliptic pole,
       p = GCRS pole,
       E = ecliptic pole of date,
    

    and P = CIP,

    the four Fukushima-Williams angles are as follows:

    GAMB = gamma = epE
    PHIB = phi = pE
    PSI = psi = pEP
    EPS = epsilon = EP
    
  2. The matrix representing the combined effects of frame bias, precession and nutation is:

    NxPxB = R_1(-EPS).R_3(-PSI).R_1(PHIB).R_3(GAMB)
    
  3. Three different matrices can be constructed, depending on the supplied angles:

    • To obtain the nutation x precession x frame bias matrix, generate the four precession angles, generate the nutation components and add them to the psi_bar and epsilon_A angles, and call the present routine.

    • To obtain the precession x frame bias matrix, generate the four precession angles and call the present routine.

    • To obtain the frame bias matrix, generate the four precession angles for date J2000.0 and call the present routine.

    The nutation-only and precession-only matrices can if necessary be obtained by combining these three appropriately.

Reference

  • Hilton, J. et al., 2006, Celest.Mech.Dyn.Astron. 94, 351

History

  • IAU SOFA revision: 2009 December 15

Arguments

TypeIntentOptionalAttributesName
real(kind=wp), intent(in) :: gamb

F-W angle gamma_bar (radians)

real(kind=wp), intent(in) :: phib

F-W angle phi_bar (radians)

real(kind=wp), intent(in) :: psi

F-W angle psi (radians)

real(kind=wp), intent(in) :: eps

F-W angle epsilon (radians)

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

rotation matrix


Calls

proc~~fw2m~~CallsGraph proc~fw2m FW2M proc~ir IR proc~fw2m->proc~ir proc~rz RZ proc~fw2m->proc~rz proc~rx RX proc~fw2m->proc~rx

Called by

proc~~fw2m~~CalledByGraph proc~fw2m FW2M proc~pnm06a PNM06A proc~pnm06a->proc~fw2m proc~bp06 BP06 proc~bp06->proc~fw2m proc~pmat06 PMAT06 proc~bp06->proc~pmat06 proc~pmat06->proc~fw2m proc~pn06 PN06 proc~pn06->proc~fw2m proc~fw2xy FW2XY proc~fw2xy->proc~fw2m proc~gst06a GST06A proc~gst06a->proc~pnm06a proc~s06a S06A proc~s06a->proc~pnm06a proc~xys06a XYS06A proc~xys06a->proc~pnm06a proc~c2i06a C2I06A proc~c2i06a->proc~pnm06a proc~apco13 APCO13 proc~apco13->proc~pnm06a proc~pn06a PN06A proc~pn06a->proc~pn06 proc~apci13 APCI13 proc~apci13->proc~pnm06a proc~eo06a EO06A proc~eo06a->proc~pnm06a proc~pb06 PB06 proc~pb06->proc~pmat06 proc~ecm06 ECM06 proc~ecm06->proc~pmat06 proc~atco13 ATCO13 proc~atco13->proc~apco13 proc~atci13 ATCI13 proc~atci13->proc~apci13 proc~atic13 ATIC13 proc~atic13->proc~apci13 proc~eceq06 ECEQ06 proc~eceq06->proc~ecm06 proc~ee06a EE06A proc~ee06a->proc~gst06a proc~c2t06a C2T06A proc~c2t06a->proc~c2i06a proc~atoc13 ATOC13 proc~atoc13->proc~apco13 proc~eqec06 EQEC06 proc~eqec06->proc~ecm06

Contents

Source Code


Source Code

    subroutine FW2M ( gamb, phib, psi, eps, r )

    implicit none

    real(wp),intent(in) :: gamb !! F-W angle gamma_bar (radians)
    real(wp),intent(in) :: phib !! F-W angle phi_bar (radians)
    real(wp),intent(in) :: psi !! F-W angle psi (radians)
    real(wp),intent(in) :: eps !! F-W angle epsilon (radians)
    real(wp),dimension(3,3),intent(out) :: r !! rotation matrix

    !  Construct the matrix.
    call IR ( r )
    call RZ ( gamb, r )
    call RX ( phib, r )
    call RZ ( -psi, r )
    call RX ( -eps, r )

    end subroutine FW2M