LTECM Subroutine

public subroutine LTECM(epj, rm)

ICRS equatorial to ecliptic rotation matrix, long-term.

Status: support routine.

Notes

  1. The matrix is in the sense

    E_ep = RM x P_ICRS,
    

    where P_ICRS is a vector with respect to ICRS right ascension and declination axes and E_ep is the same vector with respect to the (inertial) ecliptic and equinox of epoch EPJ.

  2. P_ICRS is a free vector, merely a direction, typically of unit magnitude, and not bound to any particular spatial origin, such as the Earth, Sun or SSB. No assumptions are made about whether it represents starlight and embodies astrometric effects such as parallax or aberration. The transformation is approximately that between mean J2000.0 right ascension and declination and ecliptic longitude and latitude, with only frame bias (always less than 25 mas) to disturb this classical picture.

  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):: rm

ICRS to ecliptic rotation matrix


Calls

proc~~ltecm~~CallsGraph proc~ltecm LTECM proc~ltpequ LTPEQU proc~ltecm->proc~ltpequ proc~ltpecl LTPECL proc~ltecm->proc~ltpecl proc~pxp PXP proc~ltecm->proc~pxp

Called by

proc~~ltecm~~CalledByGraph proc~ltecm LTECM proc~lteceq LTECEQ proc~lteceq->proc~ltecm proc~lteqec LTEQEC proc~lteqec->proc~ltecm

Contents

Source Code


Source Code

    subroutine LTECM ( epj, rm )

    implicit none

    real(wp),intent(in) :: epj !! Julian epoch (TT)
    real(wp),dimension(3,3),intent(out) :: rm !! ICRS to ecliptic rotation matrix

    !  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

    real(wp) :: p(3), z(3), w(3), s, x(3), y(3)

    !  Equator pole.
    call LTPEQU ( epj, p )

    !  Ecliptic pole (bottom row of equatorial to ecliptic matrix).
    call LTPECL ( epj, z )

    !  Equinox (top row of matrix).
    call PXP ( p, z, w )
    call PN ( w, s, x )

    !  Middle row of matrix.
    call PXP ( z, x, y )

    !  Combine with frame bias.
    rm(1,1) =   x(1)    - x(2)*dr + x(3)*dx
    rm(1,2) =   x(1)*dr + x(2)    + x(3)*de
    rm(1,3) = - x(1)*dx - x(2)*de + x(3)
    rm(2,1) =   y(1)    - y(2)*dr + y(3)*dx
    rm(2,2) =   y(1)*dr + y(2)    + y(3)*de
    rm(2,3) = - y(1)*dx - y(2)*de + y(3)
    rm(3,1) =   z(1)    - z(2)*dr + z(3)*dx
    rm(3,2) =   z(1)*dr + z(2)    + z(3)*de
    rm(3,3) = - z(1)*dx - z(2)*de + z(3)

    end subroutine LTECM