H2FK5 Subroutine

public subroutine H2FK5(rh, dh, drh, ddh, pxh, rvh, r5, d5, dr5, dd5, px5, rv5)

Transform Hipparcos star data into the FK5 (J2000.0) system.

Status: support routine.

Notes

  1. This routine transforms Hipparcos star positions and proper motions into FK5 J2000.0.

  2. The proper motions in RA are dRA/dt rather than cos(Dec)*dRA/dt, and are per year rather than per century.

  3. The FK5 to Hipparcos transformation is modeled as a pure rotation and spin; zonal errors in the FK5 catalog are not taken into account.

  4. See also FK52H, FK5HZ, HFK5Z.

Reference

  • F. Mignard & M. Froeschle, Astron.Astrophys., 354, 732-739 (2000).

History

  • IAU SOFA revision: 2017 October 12

Arguments

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

RA (radians) [Hipparcos, epoch J2000.0]

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

Dec (radians) [Hipparcos, epoch J2000.0]

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

proper motion in RA (dRA/dt, rad/Jyear) [Hipparcos, epoch J2000.0]

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

proper motion in Dec (dDec/dt, rad/Jyear) [Hipparcos, epoch J2000.0]

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

parallax (arcsec) [Hipparcos, epoch J2000.0]

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

radial velocity (km/s, positive = receding) [Hipparcos, epoch J2000.0]

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

RA (radians) [FK5, equinox J2000.0, epoch J2000.0]

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

Dec (radians) [FK5, equinox J2000.0, epoch J2000.0]

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

proper motion in RA (dRA/dt, rad/Jyear) [FK5, equinox J2000.0, epoch J2000.0]

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

proper motion in Dec (dDec/dt, rad/Jyear) [FK5, equinox J2000.0, epoch J2000.0]

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

parallax (arcsec) [FK5, equinox J2000.0, epoch J2000.0]

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

radial velocity (km/s, positive = receding) [FK5, equinox J2000.0, epoch J2000.0]


Calls

proc~~h2fk5~~CallsGraph proc~h2fk5 H2FK5 proc~starpv STARPV proc~h2fk5->proc~starpv proc~rxp RXP proc~h2fk5->proc~rxp proc~fk5hip FK5HIP proc~h2fk5->proc~fk5hip proc~trxp TRXP proc~h2fk5->proc~trxp proc~pxp PXP proc~h2fk5->proc~pxp proc~pmp PMP proc~h2fk5->proc~pmp proc~pvstar PVSTAR proc~h2fk5->proc~pvstar proc~starpv->proc~pmp proc~pdp PDP proc~starpv->proc~pdp proc~sxp SXP proc~starpv->proc~sxp proc~zp ZP proc~starpv->proc~zp proc~ppp PPP proc~starpv->proc~ppp proc~s2pv S2PV proc~starpv->proc~s2pv proc~rv2m RV2M proc~fk5hip->proc~rv2m proc~trxp->proc~rxp proc~tr TR proc~trxp->proc~tr proc~pvstar->proc~pmp proc~pvstar->proc~pdp proc~pvstar->proc~sxp proc~pvstar->proc~ppp proc~pv2s PV2S proc~pvstar->proc~pv2s proc~anp ANP proc~pvstar->proc~anp

Contents

Source Code


Source Code

    subroutine H2FK5 ( rh, dh, drh, ddh, pxh, rvh, &
                       r5, d5, dr5, dd5, px5, rv5 )

    implicit none

    real(wp),intent(in) :: rh !! RA (radians) [Hipparcos, epoch J2000.0]
    real(wp),intent(in) :: dh !! Dec (radians) [Hipparcos, epoch J2000.0]
    real(wp),intent(in) :: drh !! proper motion in RA (dRA/dt, rad/Jyear) [Hipparcos, epoch J2000.0]
    real(wp),intent(in) :: ddh !! proper motion in Dec (dDec/dt, rad/Jyear) [Hipparcos, epoch J2000.0]
    real(wp),intent(in) :: pxh !! parallax (arcsec) [Hipparcos, epoch J2000.0]
    real(wp),intent(in) :: rvh !! radial velocity (km/s, positive = receding) [Hipparcos, epoch J2000.0]
    real(wp),intent(out) :: r5 !! RA (radians) [FK5, equinox J2000.0, epoch J2000.0]
    real(wp),intent(out) :: d5 !! Dec (radians) [FK5, equinox J2000.0, epoch J2000.0]
    real(wp),intent(out) :: dr5 !! proper motion in RA (dRA/dt, rad/Jyear) [FK5, equinox J2000.0, epoch J2000.0]
    real(wp),intent(out) :: dd5 !! proper motion in Dec (dDec/dt, rad/Jyear) [FK5, equinox J2000.0, epoch J2000.0]
    real(wp),intent(out) :: px5 !! parallax (arcsec) [FK5, equinox J2000.0, epoch J2000.0]
    real(wp),intent(out) :: rv5 !! radial velocity (km/s, positive = receding) [FK5, equinox J2000.0, epoch J2000.0]

    real(wp) :: pvh(3,2), r5h(3,3), s5h(3), sh(3), wxp(3), &
                vv(3), pv5(3,2)
    integer :: j, i

    !  Hipparcos barycentric position/velocity pv-vector (normalized).
    call STARPV ( rh, dh, drh, ddh, pxh, rvh, pvh, j )

    !  FK5 to Hipparcos orientation matrix and spin vector.
    call FK5HIP ( r5h, s5h )

    !  Make spin units per day instead of per year.
    do i=1,3
       s5h(i) = s5h(i) / 365.25_wp
    end do

    !  Orient the spin into the Hipparcos system.
    call RXP ( r5h, s5h, sh )

    !  De-orient the Hipparcos position into the FK5 system.
    call TRXP ( r5h, pvh(1,1), pv5(1,1) )

    !  Apply spin to the position giving an extra space motion component.
    call PXP ( pvh(1,1), sh, wxp )

    !  Subtract this component from the Hipparcos space motion.
    call PMP ( pvh(1,2), wxp, vv )

    !  De-orient the Hipparcos space motion into the FK5 system.
    call TRXP ( r5h, vv, pv5(1,2) )

    !  FK5 pv-vector to spherical.
    call PVSTAR ( pv5, r5, d5, dr5, dd5, px5, rv5, j )

    end subroutine H2FK5