PMSAFE Subroutine

public subroutine PMSAFE(ra1, dec1, pmr1, pmd1, px1, rv1, ep1a, ep1b, ep2a, ep2b, ra2, dec2, pmr2, pmd2, px2, rv2, j)

Star proper motion: update star catalog data for space motion, with special handling to handle the zero parallax case.

Status: support routine.

Notes

  1. The starting and ending TDB epochs EP1A+EP1B and EP2A+EP2B are Julian Dates, apportioned in any convenient way between the two parts (A and B). For example, JD(TDB)=2450123.7 could be expressed in any of these ways, among others:

         EPnA          EPnB
    
     2450123.7D0        0D0        (JD method)
      2451545D0      -1421.3D0     (J2000 method)
     2400000.5D0     50123.2D0     (MJD method)
     2450123.5D0       0.2D0       (date & time method)
    

    The JD method is the most natural and convenient to use in cases where the loss of several decimal digits of resolution is acceptable. The J2000 method is best matched to the way the argument is handled internally and will deliver the optimum resolution. The MJD method and the date & time methods are both good compromises between resolution and convenience.

  2. In accordance with normal star-catalog conventions, the object's right ascension and declination are freed from the effects of secular aberration. The frame, which is aligned to the catalog equator and equinox, is Lorentzian and centered on the SSB.

    The proper motions are the rate of change of the right ascension and declination at the catalog epoch and are in radians per TDB Julian year.

    The parallax and radial velocity are in the same frame.

  3. Care is needed with units. The star coordinates are in radians and the proper motions in radians per Julian year, but the parallax is in arcseconds.

  4. The RA proper motion is in terms of coordinate angle, not true angle. If the catalog uses arcseconds for both RA and Dec proper motions, the RA proper motion will need to be divided by cos(Dec) before use.

  5. Straight-line motion at constant speed, in the inertial frame, is assumed.

  6. An extremely small (or zero or negative) parallax is overridden to ensure that the object is at a finite but very large distance, but not so large that the proper motion is equivalent to a large but safe speed (about 0.1c using the chosen constant). A warning status of 1 is added to the status if this action has been taken.

  7. If the space velocity is a significant fraction of c (see the constant VMAX in the routine STARPV), it is arbitrarily set to zero. When this action occurs, 2 is added to the status.

  8. The relativistic adjustment carried out in the STARPV routine involves an iterative calculation. If the process fails to converge within a set number of iterations, 4 is added to the status.

History

  • IAU SOFA revision: 2013 June 6

Arguments

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

right ascension (radians), before

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

declination (radians), before

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

RA proper motion (radians/year), before

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

Dec proper motion (radians/year), before

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

parallax (arcseconds), before

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

radial velocity (km/s, +ve = receding), before

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

"before" epoch, part A (Note 1)

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

"before" epoch, part B (Note 1)

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

"after" epoch, part A (Note 1)

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

"after" epoch, part B (Note 1)

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

right ascension (radians), after

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

declination (radians), after

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

RA proper motion (radians/year), after

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

Dec proper motion (radians/year), after

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

parallax (arcseconds), after

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

radial velocity (km/s, +ve = receding), after

integer, intent(out) :: j
  • -1 = system error (should not occur)
  • 0 = no warnings or errors
  • 1 = distance overridden (Note 6)
  • 2 = excessive velocity (Note 7)
  • 4 = solution didn't converge (Note 8)
  • else = binary logical OR of the above warnings

Calls

proc~~pmsafe~~CallsGraph proc~pmsafe PMSAFE proc~seps SEPS proc~pmsafe->proc~seps proc~starpm STARPM proc~pmsafe->proc~starpm proc~s2c S2C proc~seps->proc~s2c proc~sepp SEPP proc~seps->proc~sepp proc~starpv STARPV proc~starpm->proc~starpv proc~pdp PDP proc~starpm->proc~pdp proc~pvu PVU proc~starpm->proc~pvu proc~pvstar PVSTAR proc~starpm->proc~pvstar proc~starpv->proc~pdp proc~zp ZP proc~starpv->proc~zp proc~sxp SXP proc~starpv->proc~sxp proc~ppp PPP proc~starpv->proc~ppp proc~pmp PMP proc~starpv->proc~pmp proc~s2pv S2PV proc~starpv->proc~s2pv proc~ppsp PPSP proc~pvu->proc~ppsp proc~sepp->proc~pdp proc~pxp PXP proc~sepp->proc~pxp proc~pvstar->proc~pdp proc~pvstar->proc~sxp proc~pvstar->proc~ppp proc~pv2s PV2S proc~pvstar->proc~pv2s proc~pvstar->proc~pmp proc~anp ANP proc~pvstar->proc~anp

Contents

Source Code


Source Code

    subroutine PMSAFE ( ra1, dec1, pmr1, pmd1, px1, rv1, &
                        ep1a, ep1b, ep2a, ep2b, &
                        ra2, dec2, pmr2, pmd2, px2, rv2, j )

    implicit none

    real(wp),intent(in) :: ra1 !! right ascension (radians), before
    real(wp),intent(in) :: dec1 !! declination (radians), before
    real(wp),intent(in) :: pmr1 !! RA proper motion (radians/year), before
    real(wp),intent(in) :: pmd1 !! Dec proper motion (radians/year), before
    real(wp),intent(in) :: px1 !! parallax (arcseconds), before
    real(wp),intent(in) :: rv1 !! radial velocity (km/s, +ve = receding), before
    real(wp),intent(in) :: ep1a !! "before" epoch, part A (Note 1)
    real(wp),intent(in) :: ep1b !! "before" epoch, part B (Note 1)
    real(wp),intent(in) :: ep2a !! "after" epoch, part A (Note 1)
    real(wp),intent(in) :: ep2b !! "after" epoch, part B (Note 1)
    real(wp),intent(out) :: ra2 !! right ascension (radians), after
    real(wp),intent(out) :: dec2 !! declination (radians), after
    real(wp),intent(out) :: pmr2 !! RA proper motion (radians/year), after
    real(wp),intent(out) :: pmd2 !! Dec proper motion (radians/year), after
    real(wp),intent(out) :: px2 !! parallax (arcseconds), after
    real(wp),intent(out) :: rv2 !! radial velocity (km/s, +ve = receding), after
    integer,intent(out) :: j !! status:
                             !!  * -1 = system error (should not occur)
                             !!  *  0 = no warnings or errors
                             !!  *  1 = distance overridden (Note 6)
                             !!  *  2 = excessive velocity (Note 7)
                             !!  *  4 = solution didn't converge (Note 8)
                             !!  *  else = binary logical OR of the above warnings

    !  Minimum allowed parallax (arcsec)
    real(wp),parameter :: pxmin = 5.0e-7_wp

    !  Factor giving maximum allowed transverse speed of about 1% c
    real(wp),parameter :: f = 326.0_wp

    integer :: jpx
    real(wp) :: pm, px1a

    !  Proper motion in one year (radians).
    call SEPS ( ra1, dec1, ra1+pmr1, dec1+pmd1, pm )

    !  Override the parallax to reduce the chances of a warning status.
    jpx = 0
    px1a = px1
    pm = pm * f
    if ( px1a < pm ) then
       jpx = 1
       px1a = pm
    end if
    if ( px1a < pxmin ) then
       jpx = 1
       px1a = pxmin
    end if

    !  Carry out the transformation using the modified parallax.
    call STARPM ( ra1, dec1, pmr1, pmd1, px1a, rv1, &
                  ep1a, ep1b, ep2a, ep2b, &
                  ra2, dec2, pmr2, pmd2, px2, rv2, j )

    !  Revise the status.
    if ( mod(j,2) == 0 ) j = j + jpx

    end subroutine PMSAFE