ATCO13 Subroutine

public subroutine ATCO13(rc, dc, pr, pd, px, rv, utc1, utc2, dut1, elong, phi, hm, xp, yp, phpa, tc, rh, wl, aob, zob, hob, dob, rob, eo, j)

ICRS RA,Dec to observed place. The caller supplies UTC, site coordinates, ambient air conditions and observing wavelength.

Status: support routine.

Notes

  1. Star data for an epoch other than J2000.0 (for example from the Hipparcos catalog, which has an epoch of J1991.25) will require a preliminary call to PMSAFE before use.

  2. The proper motion in RA is dRA/dt rather than cos(Dec)*dRA/dt.

  3. UTC1+UTC2 is quasi Julian Date (see Note 2), apportioned in any convenient way between the two arguments, for example where UTC1 is the Julian Day Number and UTC2 is the fraction of a day.

    However, JD cannot unambiguously represent UTC during a leap second unless special measures are taken. The convention in the present routine is that the JD day represents UTC days whether the length is 86399, 86400 or 86401 SI seconds.

    Applications should use the routine DTF2D to convert from calendar date and time of day into 2-part quasi Julian Date, as it implements the leap-second-ambiguity convention just described.

  4. The warning status "dubious year" flags UTCs that predate the introduction of the time scale or that are too far in the future to be trusted. See DAT for further details.

  5. UT1-UTC is tabulated in IERS bulletins. It increases by exactly one second at the end of each positive UTC leap second, introduced in order to keep UT1-UTC within +/- 0.9s. n.b. This practice is under review, and in the future UT1-UTC may grow essentially without limit.

  6. The geographical coordinates are with respect to the WGS84 reference ellipsoid. TAKE CARE WITH THE LONGITUDE SIGN: the longitude required by the present routine is east-positive (i.e. right-handed), in accordance with geographical convention.

  7. The polar motion XP,YP can be obtained from IERS bulletins. The values are the coordinates (in radians) of the Celestial Intermediate Pole with respect to the International Terrestrial Reference System (see IERS Conventions 2003), measured along the meridians 0 and 90 deg west respectively. For many applications, XP and YP can be set to zero.

  8. If hm, the height above the ellipsoid of the observing station in meters, is not known but phpa, the pressure in hPa (=mB), is available, an adequate estimate of hm can be obtained from the expression

        hm = -29.3 * tsl * log ( phpa / 1013.25 );
    

    where tsl is the approximate sea-level air temperature in K (See Astrophysical Quantities, C.W.Allen, 3rd edition, section 52). Similarly, if the pressure phpa is not known, it can be estimated from the height of the observing station, hm, as follows:

        phpa = 1013.25 * exp ( -hm / ( 29.3 * tsl ) );
    

    Note, however, that the refraction is nearly proportional to the pressure and that an accurate phpa value is important for precise work.

  9. The argument WL specifies the observing wavelength in micrometers. The transition from optical to radio is assumed to occur at 100 micrometers (about 3000 GHz).

  10. The accuracy of the result is limited by the corrections for refraction, which use a simple Atan(z) + Btan^3(z) model. Providing the meteorological parameters are known accurately and there are no gross local effects, the predicted observed coordinates should be within 0.05 arcsec (optical) or 1 arcsec (radio) for a zenith distance of less than 70 degrees, better than 30 arcsec (optical or radio) at 85 degrees and better than 20 arcmin (optical) or 30 arcmin (radio) at the horizon.

    Without refraction, the complementary routines ATCO13 and ATOC13 are self-consistent to better than 1 microarcsecond all over the celestial sphere. With refraction included, consistency falls off at high zenith distances, but is still better than 0.05 arcsec at 85 degrees.

  11. "Observed" Az,ZD means the position that would be seen by a perfect geodetically aligned theodolite. (Zenith distance is used rather than altitude in order to reflect the fact that no allowance is made for depression of the horizon.) This is related to the observed HA,Dec via the standard rotation, using the geodetic latitude (corrected for polar motion), while the observed HA and RA are related simply through the Earth rotation angle and the site longitude. "Observed" RA,Dec or HA,Dec thus means the position that would be seen by a perfect equatorial with its polar axis aligned to the Earth's axis of rotation.

  12. It is advisable to take great care with units, as even unlikely values of the input parameters are accepted and processed in accordance with the models used.

History

  • IAU SOFA revision: 2016 February 2

Arguments

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

ICRS right ascension at J2000.0 (radians, Note 1)

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

ICRS declination at J2000.0 (radians, Note 1)

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

RA proper motion (radians/year; Note 2)

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

Dec proper motion (radians/year)

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

parallax (arcsec)

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

radial velocity (km/s, +ve if receding)

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

UTC as a 2-part...

real(kind=wp), intent(in) :: utc2
real(kind=wp), intent(in) :: dut1

UT1-UTC (seconds, Note 5)

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

longitude (radians, east +ve, Note 6)

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

latitude (geodetic, radians, Note 6)

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

height above ellipsoid (m, geodetic, Notes 6,8)

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

polar motion coordinate (radians, Note 7)

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

polar motion coordinate (radians, Note 7)

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

pressure at the observer (hPa = mB, Note 8)

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

ambient temperature at the observer (deg C)

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

relative humidity at the observer (range 0-1)

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

wavelength (micrometers, Note 9)

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

observed azimuth (radians: N=0,E=90)

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

observed zenith distance (radians)

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

observed hour angle (radians)

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

observed declination (radians)

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

observed right ascension (CIO-based, radians)

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

equation of the origins (ERA-GST)

integer, intent(out) :: j
  • +1 = dubious year (Note 4)
  • 0 = OK
  • -1 = unacceptable date

Calls

proc~~atco13~~CallsGraph proc~atco13 ATCO13 proc~apco13 APCO13 proc~atco13->proc~apco13 proc~atioq ATIOQ proc~atco13->proc~atioq proc~atciq ATCIQ proc~atco13->proc~atciq proc~era00 ERA00 proc~apco13->proc~era00 proc~apco APCO proc~apco13->proc~apco proc~epv00 EPV00 proc~apco13->proc~epv00 proc~taitt TAITT proc~apco13->proc~taitt proc~bpn2xy BPN2XY proc~apco13->proc~bpn2xy proc~sp00 SP00 proc~apco13->proc~sp00 proc~refco REFCO proc~apco13->proc~refco proc~utcut1 UTCUT1 proc~apco13->proc~utcut1 proc~pnm06a PNM06A proc~apco13->proc~pnm06a proc~eors EORS proc~apco13->proc~eors proc~utctai UTCTAI proc~apco13->proc~utctai proc~s06 S06 proc~apco13->proc~s06 proc~s2c S2C proc~atioq->proc~s2c proc~c2s C2S proc~atioq->proc~c2s proc~anp ANP proc~atioq->proc~anp proc~ldsun LDSUN proc~atciq->proc~ldsun proc~ab AB proc~atciq->proc~ab proc~atciq->proc~c2s proc~pmpx PMPX proc~atciq->proc~pmpx proc~rxp RXP proc~atciq->proc~rxp proc~atciq->proc~anp proc~ld LD proc~ldsun->proc~ld proc~era00->proc~anp proc~cr CR proc~apco->proc~cr proc~pvtob PVTOB proc~apco->proc~pvtob proc~apcs APCS proc~apco->proc~apcs proc~trxpv TRXPV proc~apco->proc~trxpv proc~c2ixys C2IXYS proc~apco->proc~c2ixys proc~aper APER proc~apco->proc~aper proc~pdp PDP proc~ab->proc~pdp proc~utcut1->proc~utctai proc~dat DAT proc~utcut1->proc~dat proc~jd2cal JD2CAL proc~utcut1->proc~jd2cal proc~taiut1 TAIUT1 proc~utcut1->proc~taiut1 proc~pfw06 PFW06 proc~pnm06a->proc~pfw06 proc~fw2m FW2M proc~pnm06a->proc~fw2m proc~nut06a NUT06A proc~pnm06a->proc~nut06a proc~utctai->proc~dat proc~utctai->proc~jd2cal proc~cal2jd CAL2JD proc~utctai->proc~cal2jd proc~fae03 FAE03 proc~s06->proc~fae03 proc~faf03 FAF03 proc~s06->proc~faf03 proc~faom03 FAOM03 proc~s06->proc~faom03 proc~fal03 FAL03 proc~s06->proc~fal03 proc~fave03 FAVE03 proc~s06->proc~fave03 proc~falp03 FALP03 proc~s06->proc~falp03 proc~fad03 FAD03 proc~s06->proc~fad03 proc~fapa03 FAPA03 proc~s06->proc~fapa03 proc~pmpx->proc~pdp proc~dat->proc~cal2jd proc~gd2gc GD2GC proc~pvtob->proc~gd2gc proc~pom00 POM00 proc~pvtob->proc~pom00 proc~trxp TRXP proc~pvtob->proc~trxp proc~obl06 OBL06 proc~pfw06->proc~obl06 proc~pn PN proc~apcs->proc~pn proc~ir IR proc~apcs->proc~ir proc~rxpv RXPV proc~trxpv->proc~rxpv proc~tr TR proc~trxpv->proc~tr proc~rx RX proc~fw2m->proc~rx proc~fw2m->proc~ir proc~rz RZ proc~fw2m->proc~rz proc~ry RY proc~c2ixys->proc~ry proc~c2ixys->proc~ir proc~c2ixys->proc~rz proc~ld->proc~pdp proc~pxp PXP proc~ld->proc~pxp proc~nut00a NUT00A proc~nut06a->proc~nut00a proc~gd2gce GD2GCE proc~gd2gc->proc~gd2gce proc~zp ZP proc~gd2gc->proc~zp proc~eform EFORM proc~gd2gc->proc~eform proc~sxp SXP proc~pn->proc~sxp proc~pn->proc~zp proc~rxpv->proc~rxp proc~nut00a->proc~fae03 proc~nut00a->proc~faf03 proc~nut00a->proc~faom03 proc~nut00a->proc~fal03 proc~nut00a->proc~fave03 proc~nut00a->proc~fapa03 proc~fame03 FAME03 proc~nut00a->proc~fame03 proc~faju03 FAJU03 proc~nut00a->proc~faju03 proc~fasa03 FASA03 proc~nut00a->proc~fasa03 proc~fama03 FAMA03 proc~nut00a->proc~fama03 proc~faur03 FAUR03 proc~nut00a->proc~faur03 proc~pom00->proc~rx proc~pom00->proc~ry proc~pom00->proc~ir proc~pom00->proc~rz proc~trxp->proc~rxp proc~trxp->proc~tr

Contents

Source Code


Source Code

    subroutine ATCO13 ( rc, dc, pr, pd, px, rv, &
                        utc1, utc2, dut1, elong, phi, hm, xp, yp, &
                        phpa, tc, rh, wl, &
                        aob, zob, hob, dob, rob, eo, j )

    implicit none

    real(wp),intent(in) :: rc !! ICRS right ascension at J2000.0 (radians, Note 1)
    real(wp),intent(in) :: dc !! ICRS declination at J2000.0 (radians, Note 1)
    real(wp),intent(in) :: pr !! RA proper motion (radians/year; Note 2)
    real(wp),intent(in) :: pd !! Dec proper motion (radians/year)
    real(wp),intent(in) :: px !! parallax (arcsec)
    real(wp),intent(in) :: rv !! radial velocity (km/s, +ve if receding)
    real(wp),intent(in) :: utc1 !! UTC as a 2-part...
    real(wp),intent(in) :: utc2 !! ...quasi Julian Date (Notes 3-4)
    real(wp),intent(in) :: dut1 !! UT1-UTC (seconds, Note 5)
    real(wp),intent(in) :: elong !! longitude (radians, east +ve, Note 6)
    real(wp),intent(in) :: phi !! latitude (geodetic, radians, Note 6)
    real(wp),intent(in) :: hm !! height above ellipsoid (m, geodetic, Notes 6,8)
    real(wp),intent(in) :: xp !! polar motion coordinate (radians, Note 7)
    real(wp),intent(in) :: yp !! polar motion coordinate (radians, Note 7)
    real(wp),intent(in) :: phpa !! pressure at the observer (hPa = mB, Note 8)
    real(wp),intent(in) :: tc !! ambient temperature at the observer (deg C)
    real(wp),intent(in) :: rh !! relative humidity at the observer (range 0-1)
    real(wp),intent(in) :: wl !! wavelength (micrometers, Note 9)
    real(wp),intent(out) :: aob !! observed azimuth (radians: N=0,E=90)
    real(wp),intent(out) :: zob !! observed zenith distance (radians)
    real(wp),intent(out) :: hob !! observed hour angle (radians)
    real(wp),intent(out) :: dob !! observed declination (radians)
    real(wp),intent(out) :: rob !! observed right ascension (CIO-based, radians)
    real(wp),intent(out) :: eo !! equation of the origins (ERA-GST)
    integer,intent(out) :: j !! status:
                             !! * +1 = dubious year (Note 4)
                             !! *  0 = OK
                             !! * -1 = unacceptable date

    integer :: js
    real(wp) :: astrom(30)
    real(wp) :: ri, di

    !  Star-independent astrometry parameters.
    call APCO13 ( utc1, utc2, dut1, elong, phi, hm, xp, yp, &
                  phpa, tc, rh, wl, astrom, eo, js )
    if ( js>=0 ) then

      !  Transform ICRS to CIRS.
      call ATCIQ ( rc, dc, pr, pd, px, rv, astrom, ri, di )

      !  Transform CIRS to observed.
      call ATIOQ ( ri, di, astrom, aob, zob, hob, dob, rob )

    end if

    !  Return OK/warning status.
    j = js

    end subroutine ATCO13