FW2XY Subroutine

public subroutine FW2XY(gamb, phib, psi, eps, x, y)

CIP X,Y given Fukushima-Williams bias-precession-nutation 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(-EPSA).R_3(-PSI).R_1(PHIB).R_3(GAMB)
    

    The returned values x,y are elements (3,1) and (3,2) of the matrix. Near J2000.0, they are essentially angles in radians.

Reference

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

History

  • IAU SOFA revision: 2013 September 2

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) :: x

CIP unit vector X,Y

real(kind=wp), intent(out) :: y

CIP unit vector X,Y


Calls

proc~~fw2xy~~CallsGraph proc~fw2xy FW2XY proc~fw2m FW2M proc~fw2xy->proc~fw2m proc~bpn2xy BPN2XY proc~fw2xy->proc~bpn2xy proc~ir IR proc~fw2m->proc~ir proc~rz RZ proc~fw2m->proc~rz proc~rx RX proc~fw2m->proc~rx

Contents

Source Code


Source Code

    subroutine FW2XY ( gamb, phib, psi, eps, x, y )

    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),intent(out) :: x !! CIP unit vector X,Y
    real(wp),intent(out) :: y !! CIP unit vector X,Y

    real(wp) :: r(3,3)

    !  Form NxPxB matrix.
    call FW2M ( gamb, phib, psi, eps, r )

    !  Extract CIP X,Y.
    call BPN2XY ( r, x, y )

    end subroutine FW2XY