Quick CIRS RA,Dec to ICRS astrometric place, 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 one of the routines APCI[13], APCG[13], APCO[13] or APCS[13].
Status: support routine.
Only the Sun is taken into account in the light deflection correction.
Iterative techniques are used for the aberration and light deflection corrections so that the routines ATIC13 (or ATICQ) and ATCI13 (or ATCIQ) are accurate inverses; even at the edge of the Sun's disk the discrepancy is only about 1 nanoarcsecond.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | ri | CIRS RA (radians) |
||
real(kind=wp), | intent(in) | :: | di | CIRS 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) | :: | rc | ICRS astrometric RA (radians) |
||
real(kind=wp), | intent(out) | :: | dc | ICRS astrometric Dec (radians) |
subroutine ATICQ ( ri, di, astrom, rc, dc )
implicit none
real(wp),intent(in) :: ri !! CIRS RA (radians)
real(wp),intent(in) :: di !! CIRS 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) :: rc !! ICRS astrometric RA (radians)
real(wp),intent(out) :: dc !! ICRS astrometric Dec (radians)
integer :: j, i
real(wp) :: pi(3), ppr(3), pnat(3), pco(3), w, d(3), &
before(3), r2, r, after(3)
! CIRS RA,Dec to Cartesian.
call S2C ( ri, di, pi )
! Bias-precession-nutation, giving GCRS proper direction.
call TRXP ( astrom(13), pi, ppr )
! Aberration, giving GCRS natural direction.
call ZP ( d )
do j=1,2
r2 = 0.0_wp
do i=1,3
w = ppr(i) - d(i)
before(i) = w
r2 = r2 + w*w
end do
r = sqrt ( r2 )
do i=1,3
before(i) = before(i) / r
end do
call AB ( before, astrom(9), astrom(8), astrom(12), after )
r2 = 0.0_wp
do i=1,3
d(i) = after(i) - before(i)
w = ppr(i) - d(i)
pnat(i) = w
r2 = r2 + w*w
end do
r = sqrt ( r2 )
do i=1,3
pnat(i) = pnat(i) / r
end do
end do
! Light deflection by the Sun, giving BCRS coordinate direction.
call ZP ( d )
do j=1,5
r2 = 0.0_wp
do i=1,3
w = pnat(i) - d(i)
before(i) = w
r2 = r2 + w*w
end do
r = sqrt ( r2 )
do i=1,3
before(i) = before(i) / r
end do
call LDSUN ( before, astrom(5), astrom(8), after )
r2 = 0.0_wp
do i=1,3
d(i) = after(i) - before(i)
w = pnat(i) - d(i)
pco(i) = w
r2 = r2 + w*w
end do
r = sqrt ( r2 )
do i=1,3
pco(i) = pco(i) / r
end do
end do
! ICRS astrometric RA,Dec.
call C2S ( pco, w, dc )
rc = ANP ( w )
end subroutine ATICQ