LTPB Subroutine

public subroutine LTPB(epj, rpb)

Long-term precession matrix, including ICRS frame bias.

Status: support routine.

Notes

  1. The matrix is in the sense

    P_date = RPB x P_ICRS,
    

    where P_J2000 is a vector in the International Celestial Reference System, and P_date is the vector with respect to the Celestial Intermediate Reference System at that date but with nutation neglected.

  2. A first order frame bias formulation is used, of sub- microarcsecond accuracy compared with a full 3D rotation.

  3. 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):: rpb

precession-bias matrix, J2000.0 to date


Calls

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

Contents

Source Code


Source Code

    subroutine LTPB ( epj, rpb )

    implicit none

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

    !  Frame bias (IERS Conventions 2010, Eqs. 5.21 and 5.33)
    real(wp),parameter :: dx = -0.016617_wp * das2r
    real(wp),parameter :: de = -0.0068192_wp * das2r
    real(wp),parameter :: dr = -0.0146_wp * das2r

    integer :: i
    real(wp) :: rp(3,3)

    !  Precession matrix.
    call LTP  ( epj, rp )

    !  Apply the bias.
    do i=1,3
       rpb(i,1) =   rp(i,1)    - rp(i,2)*dr + rp(i,3)*dx
       rpb(i,2) =   rp(i,1)*dr + rp(i,2)    + rp(i,3)*de
       rpb(i,3) = - rp(i,1)*dx - rp(i,2)*de + rp(i,3)
    end do

    end subroutine LTPB