FK52H Subroutine

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

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

Status: support routine.

Notes

  1. This routine transforms FK5 star positions and proper motions into the system of the Hipparcos catalog.

  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 H2FK5, 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) :: r5

RA (radians)

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

Dec (radians)

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

proper motion in RA (dRA/dt, rad/Jyear)

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

proper motion in Dec (dDec/dt, rad/Jyear)

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

parallax (arcsec)

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

radial velocity (km/s, positive = receding)

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

RA (radians)

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

Dec (radians)

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

proper motion in RA (dRA/dt, rad/Jyear)

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

proper motion in Dec (dDec/dt, rad/Jyear)

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

parallax (arcsec)

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

radial velocity (km/s, positive = receding)


Calls

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

Contents

Source Code


Source Code

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

    implicit none

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

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

    !  FK5 barycentric position/velocity pv-vector (normalized).
    call STARPV ( r5, d5, dr5, dd5, px5, rv5, pv5, 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 FK5 position into the Hipparcos system.
    call RXP ( r5h, pv5(1,1), pvh(1,1) )

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

    !  Add this component to the FK5 space motion.
    call PPP ( wxp, pv5(1,2), vv )

    !  Orient the FK5 space motion into the Hipparcos system.
    call RXP ( r5h, vv, pvh(1,2) )

    !  Hipparcos pv-vector to spherical.
    call PVSTAR ( pvh, rh, dh, drh, ddh, pxh, rvh, j )

    end subroutine FK52H