LTEQEC Subroutine

public subroutine LTEQEC(epj, dr, dd, dl, db)

Transformation from ICRS equatorial coordinates to ecliptic coordinates (mean equinox and ecliptic of date), using a long-term precession model.

Status: support routine.

  1. No assumptions are made about whether the coordinates represent starlight and embody astrometric effects such as parallax or aberration.

  2. The transformation is approximately that from mean J2000.0 right ascension and declination to ecliptic longitude and latitude (mean equinox and ecliptic of date), 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: 2016 February 9

Arguments

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

Julian epoch (TT)

real(kind=wp), intent(in) :: dr

ICRS right ascension (radians)

real(kind=wp), intent(in) :: dd

ICRS right declination (radians)

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

ecliptic longitude (radians)

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

ecliptic latitude (radians)


Calls

proc~~lteqec~~CallsGraph proc~lteqec LTEQEC proc~c2s C2S proc~lteqec->proc~c2s proc~anpm ANPM proc~lteqec->proc~anpm proc~rxp RXP proc~lteqec->proc~rxp proc~ltecm LTECM proc~lteqec->proc~ltecm proc~s2c S2C proc~lteqec->proc~s2c proc~anp ANP proc~lteqec->proc~anp proc~ltpequ LTPEQU proc~ltecm->proc~ltpequ proc~ltpecl LTPECL proc~ltecm->proc~ltpecl proc~pxp PXP proc~ltecm->proc~pxp

Contents

Source Code


Source Code

    subroutine LTEQEC ( epj, dr, dd, dl, db )

    implicit none

    real(wp),intent(in) :: epj !! Julian epoch (TT)
    real(wp),intent(in) :: dr !! ICRS right ascension (radians)
    real(wp),intent(in) :: dd !! ICRS right declination (radians)
    real(wp),intent(out) :: dl !! ecliptic longitude (radians)
    real(wp),intent(out) :: db !! ecliptic latitude (radians)

    real(wp) :: rm(3,3), v1(3), v2(3), a, b

    !  Spherical to Cartesian.
    call S2C ( dr, dd, v1 )

    !  Rotation matrix, ICRS equatorial to ecliptic.
    call LTECM ( epj, rm )

    !  The transformation from ICRS to ecliptic.
    call RXP ( rm, v1, v2 )

    !  Cartesian to spherical.
    call C2S ( v2, a, b )

    !  Express in conventional ranges.
    dl = ANP ( a )
    db = ANPM ( b )

    end subroutine LTEQEC