For a terrestrial observer, prepare star-independent astrometry parameters for transformations between CIRS and observed coordinates. The caller supplies the Earth orientation information and the refraction constants as well as the site coordinates.
Status: support routine.
SP, the TIO locator s', is a tiny quantity needed only by the most precise applications. It can either be set to zero or predicted using the SOFA routine SP00.
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.
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.
Internally, the polar motion is stored in a form rotated onto the local meridian.
The refraction constants REFA and REFB are for use in a dZ = Atan(Z)+Btan^3(Z) model, where Z is the observed (i.e. refracted) zenith distance and dZ is the amount of refraction.
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.
In cases where the caller does not wish to provide the Earth rotation information and refraction constants, the routine APIO13 can be used instead of the present routine. This starts from UTC and weather readings etc. and computes suitable values using other SOFA routines.
This is one of several routines that inserts into the ASTROM array star-independent parameters needed for the chain of astrometric transformations ICRS <-> GCRS <-> CIRS <-> observed.
The various routines support different classes of observer and portions of the transformation chain:
routines observer transformation
APCG APCG13 geocentric ICRS <-> GCRS
APCI APCI13 terrestrial ICRS <-> CIRS
APCO APCO13 terrestrial ICRS <-> observed
APCS APCS13 space ICRS <-> GCRS
APER APER13 terrestrial update Earth rotation
APIO APIO13 terrestrial CIRS <-> observed
Those with names ending in "13" use contemporary SOFA models to compute the various ephemerides. The others accept ephemerides supplied by the caller.
The transformation from ICRS to GCRS covers space motion, parallax, light deflection, and aberration. From GCRS to CIRS comprises frame bias and precession-nutation. From CIRS to observed takes account of Earth rotation, polar motion, diurnal aberration and parallax (unless subsumed into the ICRS <-> GCRS transformation), and atmospheric refraction.
The context array ASTROM produced by this routine is used by ATIOQ and ATOIQ.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | sp | the TIO locator s' (radians, Note 1) |
||
real(kind=wp), | intent(in) | :: | theta | Earth rotation angle (radians) |
||
real(kind=wp), | intent(in) | :: | elong | longitude (radians, east +ve, Note 2) |
||
real(kind=wp), | intent(in) | :: | phi | geodetic latitude (radians, Note 2) |
||
real(kind=wp), | intent(in) | :: | hm | height above ellipsoid (m, geodetic Note 2) |
||
real(kind=wp), | intent(in) | :: | xp | polar motion coordinate (radians, Note 3) |
||
real(kind=wp), | intent(in) | :: | yp | polar motion coordinate (radians, Note 3) |
||
real(kind=wp), | intent(in) | :: | refa | refraction constant A (radians, Note 4) |
||
real(kind=wp), | intent(in) | :: | refb | refraction constant B (radians, Note 4) |
||
real(kind=wp), | intent(inout), | dimension(30) | :: | astrom | star-independent astrometry parameters: (1) unchanged (2-4) unchanged (5-7) unchanged (8) unchanged (9-11) unchanged (12) unchanged (13-21) unchanged (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) |
subroutine APIO ( sp, theta, elong, phi, hm, xp, yp, &
refa, refb, astrom )
implicit none
real(wp),intent(in) :: sp !! the TIO locator s' (radians, Note 1)
real(wp),intent(in) :: theta !! Earth rotation angle (radians)
real(wp),intent(in) :: elong !! longitude (radians, east +ve, Note 2)
real(wp),intent(in) :: phi !! geodetic latitude (radians, Note 2)
real(wp),intent(in) :: hm !! height above ellipsoid (m, geodetic Note 2)
real(wp),intent(in) :: xp !! polar motion coordinate (radians, Note 3)
real(wp),intent(in) :: yp !! polar motion coordinate (radians, Note 3)
real(wp),intent(in) :: refa !! refraction constant A (radians, Note 4)
real(wp),intent(in) :: refb !! refraction constant B (radians, Note 4)
real(wp),dimension(30),intent(inout) :: astrom !! star-independent astrometry parameters:
!!
!! (1) unchanged
!! (2-4) unchanged
!! (5-7) unchanged
!! (8) unchanged
!! (9-11) unchanged
!! (12) unchanged
!! (13-21) unchanged
!! (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) :: sl, cl, pv(3,2)
! Longitude with adjustment for TIO locator s'.
astrom(22) = elong + sp
! Polar motion, rotated onto the local meridian.
sl = sin(astrom(22))
cl = cos(astrom(22))
astrom(23) = xp*cl - yp*sl
astrom(24) = xp*sl + yp*cl
! Functions of latitude.
astrom(25) = sin(phi)
astrom(26) = cos(phi)
! Observer's geocentric position and velocity (m, m/s, CIRS).
call PVTOB ( elong, phi, hm, xp, yp, sp, theta, pv )
! Magnitude of diurnal aberration vector.
astrom(27) = sqrt ( pv(1,2)*pv(1,2) + pv(2,2)*pv(2,2) ) / cmps
! Refraction constants.
astrom(29) = refa
astrom(30) = refb
! Local Earth rotation angle.
call APER ( theta, astrom )
end subroutine APIO