Quick observed place to CIRS, given the star-independent astrometry parameters.
Use of this routine is appropriate when efficiency is important and where many star positions are all to be transformed for one date. The star-independent astrometry parameters can be obtained by calling APIO[13] or APCO[13].
Status: support routine.
"Observed" Az,El means the position that would be seen by a perfect geodetically aligned theodolite. 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. By removing from the observed place the effects of atmospheric refraction and diurnal aberration, the CIRS RA,Dec is obtained.
Only the first character of the type argument is significant. 'R' or 'r' indicates that OB1 and OB2 are the observed right ascension and declination; 'H' or 'h' indicates that they are hour angle (west +ve) and declination; anything else ('A' or 'a' is recommended) indicates that OB1 and OB2 are azimuth (north zero, east 90 deg) and zenith distance. (Zenith distance is used rather than altitude in order to reflect the fact that no allowance is made for depression of the horizon.)
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 0D05 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 ATIOQ and ATOIQ 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.
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.
The star-independent astrometry parameters in ASTROM may be computed with APIO13 (or APIO). If nothing has changed significantly except the time, APER13 (or APER) may be used to perform the requisite adjustment to the ASTROM array.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=*), | intent(in) | :: | type | type of coordinates: 'R', 'H' or 'A' (Note 2) |
||
real(kind=wp), | intent(in) | :: | ob1 | observed Az, HA or RA (radians; Az is N=0,E=90) |
||
real(kind=wp), | intent(in) | :: | ob2 | observed ZD or Dec (radians) |
||
real(kind=wp), | intent(in), | dimension(30) | :: | astrom | star-independent astrometry parameters: (1) PM time interval (SSB, Julian years) (2-4) SSB to observer (vector, au) (5-7) Sun to observer (unit vector) (8) distance from Sun to observer (au) (9-11) v: barycentric observer velocity (vector, c) (12) sqrt(1-|v|^2): reciprocal of Lorenz factor (13-21) bias-precession-nutation matrix (22) longitude + s' (radians) (23) polar motion xp wrt local meridian (radians) (24) polar motion yp wrt local meridian (radians) (25) sine of geodetic latitude (26) cosine of geodetic latitude (27) magnitude of diurnal aberration vector (28) "local" Earth rotation angle (radians) (29) refraction constant A (radians) (30) refraction constant B (radians) |
|
real(kind=wp), | intent(out) | :: | ri | CIRS right ascension (CIO-based, radians) |
||
real(kind=wp), | intent(out) | :: | di | CIRS declination (radians) |
subroutine ATOIQ ( type, ob1, ob2, astrom, ri, di )
implicit none
character(len=*),intent(in) :: type !! type of coordinates: 'R', 'H' or 'A' (Note 2)
real(wp),intent(in) :: ob1 !! observed Az, HA or RA (radians; Az is N=0,E=90)
real(wp),intent(in) :: ob2 !! observed ZD or Dec (radians)
real(wp),dimension(30),intent(in) :: astrom !! star-independent astrometry parameters:
!! (1) PM time interval (SSB, Julian years)
!! (2-4) SSB to observer (vector, au)
!! (5-7) Sun to observer (unit vector)
!! (8) distance from Sun to observer (au)
!! (9-11) v: barycentric observer velocity (vector, c)
!! (12) sqrt(1-|v|^2): reciprocal of Lorenz factor
!! (13-21) bias-precession-nutation matrix
!! (22) longitude + s' (radians)
!! (23) polar motion xp wrt local meridian (radians)
!! (24) polar motion yp wrt local meridian (radians)
!! (25) sine of geodetic latitude
!! (26) cosine of geodetic latitude
!! (27) magnitude of diurnal aberration vector
!! (28) "local" Earth rotation angle (radians)
!! (29) refraction constant A (radians)
!! (30) refraction constant B (radians)
real(wp),intent(out) :: ri !! CIRS right ascension (CIO-based, radians)
real(wp),intent(out) :: di !! CIRS declination (radians)
character(len=1) :: c
real(wp) :: c1, c2, sphi, cphi, ce, xaeo, yaeo, zaeo, v(3), &
xmhdo, ymhdo, zmhdo, az, sz, zdo, refa, refb, &
tz, dref, zdt, xaet, yaet, zaet, &
xmhda, ymhda, zmhda, f, xhd, yhd, zhd, &
xpl, ypl, w, hma
! Coordinate type.
c = type(:1)
! Coordinates.
c1 = ob1
c2 = ob2
! Sin, cos of latitude.
sphi = astrom(25)
cphi = astrom(26)
! Standardize coordinate type.
if ( c=='r' .or. c=='R' ) then
c = 'R'
else if ( c=='h' .or. c=='H' ) then
c = 'H'
else
c = 'A'
end if
! If Az,ZD, convert to Cartesian (S=0,E=90).
if ( c=='A' ) then
ce = sin(c2)
xaeo = - cos(c1) * ce
yaeo = sin(c1) * ce
zaeo = cos(c2)
else
! If RA,Dec, convert to HA,Dec.
if ( c=='R' ) c1 = astrom(28) - c1
! To Cartesian -HA,DeC.
call S2C ( -c1, c2, v )
xmhdo = v(1)
ymhdo = v(2)
zmhdo = v(3)
! To Cartesian Az,El (S=0,E=90).
xaeo = sphi*xmhdo - cphi*zmhdo
yaeo = ymhdo
zaeo = cphi*xmhdo + sphi*zmhdo
end if
! Azimuth (S=0,E=90).
if ( xaeo/=0.0_wp .or. yaeo/=0.0_wp ) then
az = atan2(yaeo,xaeo)
else
az = 0.0_wp
end if
! Sine of observed ZD, and observed ZD.
sz = sqrt ( xaeo*xaeo + yaeo*yaeo )
zdo = atan2 ( sz, zaeo )
!
! Refraction
! ----------
! Fast algorithm using two constant model.
refa = astrom(29)
refb = astrom(30)
tz = sz / zaeo
dref = ( refa + refb*tz*tz ) * tz
zdt = zdo + dref
! To Cartesian Az,ZD.
ce = sin(zdt)
xaet = cos(az) * ce
yaet = sin(az) * ce
zaet = cos(zdt)
! Cartesian Az,ZD to Cartesian -HA,DeC.
xmhda = sphi*xaet + cphi*zaet
ymhda = yaet
zmhda = - cphi*xaet + sphi*zaet
! Diurnal aberration.
f = ( 1.0_wp + astrom(27)*ymhda )
xhd = f * xmhda
yhd = f * ( ymhda - astrom(27) )
zhd = f * zmhda
! Polar motion.
xpl = astrom(23)
ypl = astrom(24)
w = xpl*xhd - ypl*yhd + zhd
v(1) = xhd - xpl*w
v(2) = yhd + ypl*w
v(3) = w - ( xpl*xpl + ypl*ypl ) * zhd
! To spherical -HA,DeC.
call C2S ( v, hma, di )
! Right ascension.
ri = ANP ( astrom(28) + hma )
end subroutine ATOIQ