ATICQN Subroutine

public subroutine ATICQN(ri, di, astrom, n, b, rc, dc)

Quick CIRS to ICRS astrometric place transformation, given the star-independent astrometry parameters plus a list of light- deflecting bodies.

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].

If the only light-deflecting body to be taken into account is the Sun, the ATICQ routine can be used instead.

Status: support routine.

Notes

  1. Iterative techniques are used for the aberration and light deflection corrections so that the routines ATICQN and ATCIQN are accurate inverses; even at the edge of the Sun's disk the discrepancy is only about 1 nanoarcsecond.

  2. If the only light-deflecting body to be taken into account is the Sun, the ATICQ routine can be used instead.

  3. The array B contains N entries, one for each body to be considered. If N = 0, no gravitational light deflection will be applied, not even for the Sun.

  4. The array B should include an entry for the Sun as well as for any planet or other body to be taken into account. The entries should be in the order in which the light passes the body.

  5. In the entry in the B array for body I, the mass parameter B(1,I) can, as required, be adjusted in order to allow for such effects as quadrupole field.

  6. The deflection limiter parameter B(2,I) is phi^2/2, where phi is the angular separation (in radians) between star and body at which limiting is applied. As phi shrinks below the chosen threshold, the deflection is artificially reduced, reaching zero for phi = 0. Example values suitable for a terrestrial observer, together with masses, are as follows:

    body I     B(1,I)         B(2,I)
    
    Sun        1D0            6D-6
    Jupiter    0.00095435D0   3D-9
    Saturn     0.00028574D0   3D-10
    
  7. For efficiency, validation of the contents of the B array is omitted. The supplied masses must be greater than zero, the position and velocity vectors must be right, and the deflection limiter greater than zero.

History

  • IAU SOFA revision: 2013 September 30

Arguments

TypeIntentOptionalAttributesName
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)

integer, intent(in) :: n

number of bodies (Note 3)

real(kind=wp), intent(in), dimension(8,n):: b

data for each of the NB bodies (Notes 3,4): (1,I) mass of the body (solar masses, Note 5) (2,I) deflection limiter (Note 6) (3-5,I) barycentric position of the body (au) (6-8,I) barycentric velocity of the body (au/day)

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

ICRS astrometric RA (radians)

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

ICRS astrometric Dec (radians)


Calls

proc~~aticqn~~CallsGraph proc~aticqn ATICQN proc~ldn LDN proc~aticqn->proc~ldn proc~zp ZP proc~aticqn->proc~zp proc~c2s C2S proc~aticqn->proc~c2s proc~ab AB proc~aticqn->proc~ab proc~s2c S2C proc~aticqn->proc~s2c proc~trxp TRXP proc~aticqn->proc~trxp proc~anp ANP proc~aticqn->proc~anp proc~pdp PDP proc~ldn->proc~pdp proc~ppsp PPSP proc~ldn->proc~ppsp proc~pmp PMP proc~ldn->proc~pmp proc~ld LD proc~ldn->proc~ld proc~ab->proc~pdp proc~rxp RXP proc~trxp->proc~rxp proc~tr TR proc~trxp->proc~tr proc~ld->proc~pdp proc~pxp PXP proc~ld->proc~pxp

Contents

Source Code


Source Code

    subroutine ATICQN ( ri, di, astrom, n, b, 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)
    integer,intent(in) :: n !! number of bodies (Note 3)
    real(wp),dimension(8,n),intent(in) :: b !! data for each of the NB bodies (Notes 3,4):
                                            !! (1,I)    mass of the body (solar masses, Note 5)
                                            !! (2,I)    deflection limiter (Note 6)
                                            !! (3-5,I)  barycentric position of the body (au)
                                            !! (6-8,I)  barycentric velocity of the body (au/day)
    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, 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 LDN ( n, b, astrom(2), before, 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 ATICQN