LTP Subroutine

public subroutine LTP(epj, rp)

Long-term precession matrix.

Status: support routine.

Notes

  1. The matrix is in the sense

    P_date = RP x P_J2000,
    

    where P_J2000 is a vector with respect to the J2000.0 mean equator and equinox and P_date is the same vector with respect to the equator and equinox of epoch EPJ.

  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: 2015 December 6

Arguments

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

Julian epoch (TT)

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

precession matrix, J2000.0 to date


Calls

proc~~ltp~~CallsGraph proc~ltp LTP proc~ltpequ LTPEQU proc~ltp->proc~ltpequ proc~ltpecl LTPECL proc~ltp->proc~ltpecl proc~pxp PXP proc~ltp->proc~pxp

Called by

proc~~ltp~~CalledByGraph proc~ltp LTP proc~ltpb LTPB proc~ltpb->proc~ltp

Contents

Source Code

LTP

Source Code

    subroutine LTP ( epj, rp )

    implicit none

    real(wp),intent(in) :: epj !! Julian epoch (TT)
    real(wp),dimension(3,3),intent(out) :: rp !! precession matrix, J2000.0 to date

    integer :: i
    real(wp) :: peqr(3), pecl(3), v(3), w, eqx(3)

    !  Equator pole (bottom row of matrix).
    call LTPEQU ( epj, peqr )

    !  Ecliptic pole.
    call LTPECL ( epj, pecl )

    !  Equinox (top row of matrix).
    call PXP ( peqr, pecl, v )
    call PN ( v, w, eqx )

    !  Middle row of matrix.
    call PXP ( peqr, eqx, v )

    !  Assemble the matrix.
    do i=1,3
       rp(1,i) = eqx(i)
       rp(2,i) = v(i)
       rp(3,i) = peqr(i)
    end do

    end subroutine LTP