LTPECL Subroutine

public subroutine LTPECL(epj, vec)

Long-term precession of the ecliptic.

Status: support routine.

Notes

  1. The returned vector is with respect to the J2000.0 mean equator and equinox.

  2. The Vondrak et al. (2011, 2012) 400 millennia precession model agrees with the IAU 2006 precession at J2000.0 and stays within 100 microarcseconds during the 20th and 21st centuries. It is accurate to a few arcseconds throughout the historical period, worsening to a few tenths of a degree at the end of the +/- 200,000 year time span.

References

  • Vondrak, J., Capitaine, N. and Wallace, P., 2011, New precession expressions, valid for long time intervals, Astron.Astrophys. 534, A22

  • Vondrak, J., Capitaine, N. and Wallace, P., 2012, New precession expressions, valid for long time intervals (Corrigendum), Astron.Astrophys. 541, C1

History

  • IAU SOFA revision: 2016 February 9

Arguments

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

Julian epoch (TT)

real(kind=wp), intent(out), dimension(3):: vec

ecliptic pole unit vector


Called by

proc~~ltpecl~~CalledByGraph proc~ltpecl LTPECL proc~ltecm LTECM proc~ltecm->proc~ltpecl proc~ltp LTP proc~ltp->proc~ltpecl proc~lteceq LTECEQ proc~lteceq->proc~ltecm proc~ltpb LTPB proc~ltpb->proc~ltp proc~lteqec LTEQEC proc~lteqec->proc~ltecm

Contents

Source Code


Source Code

    subroutine LTPECL ( epj, vec )

    implicit none

    real(wp),intent(in) :: epj !! Julian epoch (TT)
    real(wp),dimension(3),intent(out) :: vec !! ecliptic pole unit vector

    !  Obliquity at J2000.0 (radians).
    real(wp),parameter :: eps0 = 84381.406_wp * das2r

    !  Number of polynomial terms
    integer,parameter :: npol = 4

    !  Number of periodic terms
    integer,parameter :: nper = 8

    !  Miscellaneous
    integer :: i, j
    real(wp) :: t, p, q, w, a, s, c

    !  Polynomial and periodic coefficients

    !  Polynomials
    real(wp),dimension(npol,2),parameter :: pqpol = reshape([&
           +5851.607687_wp, &
              -0.1189000_wp, &
              -0.00028913_wp, &
              +0.000000101_wp, &
           -1600.886300_wp, &
              +1.1689818_wp, &
              -0.00000020_wp, &
              -0.000000437_wp ], [npol,2])

    !  Periodics
    real(wp),dimension(5,nper),parameter :: pqper = reshape([&
       708.15_wp,  -5486.751211_wp,   -684.661560_wp, &
                     667.666730_wp,  -5523.863691_wp, &
      2309.00_wp,    -17.127623_wp,   2446.283880_wp, &
                   -2354.886252_wp,   -549.747450_wp, &
      1620.00_wp,   -617.517403_wp,    399.671049_wp, &
                    -428.152441_wp,   -310.998056_wp, &
       492.20_wp,    413.442940_wp,   -356.652376_wp, &
                     376.202861_wp,    421.535876_wp, &
      1183.00_wp,     78.614193_wp,   -186.387003_wp, &
                     184.778874_wp,    -36.776172_wp, &
       622.00_wp,   -180.732815_wp,   -316.800070_wp, &
                     335.321713_wp,   -145.278396_wp, &
       882.00_wp,    -87.676083_wp,    198.296701_wp, &
                    -185.138669_wp,    -34.744450_wp, &
       547.00_wp,     46.140315_wp,    101.135679_wp, &
                    -120.972830_wp,     22.885731_wp ], [5,nper])

    !  Centuries since J2000.
    t = (epj-2000.0_wp)/100.0_wp

    !  Initialize P_A and Q_A accumulators.
    p = 0.0_wp
    q = 0.0_wp

    !  Periodic terms.
    w = d2pi*t
    do i=1,nper
       a = w/pqper(1,i)
       s = sin(a)
       c = cos(a)
       p = p + c*pqper(2,i) + s*pqper(4,i)
       q = q + c*pqper(3,i) + s*pqper(5,i)
    end do

    !  Polynomial terms.
    w = 1.0_wp
    do i=1,npol
       p = p + pqpol(i,1)*w
       q = q + pqpol(i,2)*w
       w = w*t
    end do

    !  P_A and Q_A (radians).
    p = p*das2r
    q = q*das2r

    !  Form the ecliptic pole vector.
    w = sqrt(max(1.0_wp-p*p-q*q,0.0_wp))
    s = sin(eps0)
    c = cos(eps0)
    vec(1) = p
    vec(2) = - q*c - w*s
    vec(3) = - q*s + w*c

    end subroutine LTPECL