NUTM80 Subroutine

public subroutine NUTM80(date1, date2, rmatn)

Form the matrix of nutation for a given date, IAU 1980 model.

Status: support routine.

Notes

  1. The TDB date DATE1+DATE2 is a Julian Date, apportioned in any convenient way between the two arguments. For example, JD(TDB)=2450123.7 could be expressed in any of these ways, among others:

        DATE1          DATE2
    
     2450123.7D0        0D0        (JD method)
      2451545D0      -1421.3D0     (J2000 method)
     2400000.5D0     50123.2D0     (MJD method)
     2450123.5D0       0.2D0       (date & time method)
    

    The JD method is the most natural and convenient to use in cases where the loss of several decimal digits of resolution is acceptable. The J2000 method is best matched to the way the argument is handled internally and will deliver the optimum resolution. The MJD method and the date & time methods are both good compromises between resolution and convenience.

  2. The matrix operates in the sense V(true) = RMATN * V(mean), where the p-vector V(true) is with respect to the true equatorial triad of date and the p-vector V(mean) is with respect to the mean equatorial triad of date.

History

  • IAU SOFA revision: 2012 September 5

Arguments

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

TDB date (Note 1)

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

TDB date (Note 1)

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

nutation matrix


Calls

proc~~nutm80~~CallsGraph proc~nutm80 NUTM80 proc~nut80 NUT80 proc~nutm80->proc~nut80 proc~obl80 OBL80 proc~nutm80->proc~obl80 proc~numat NUMAT proc~nutm80->proc~numat proc~anpm ANPM proc~nut80->proc~anpm proc~ir IR proc~numat->proc~ir proc~rz RZ proc~numat->proc~rz proc~rx RX proc~numat->proc~rx

Called by

proc~~nutm80~~CalledByGraph proc~nutm80 NUTM80 proc~pnm80 PNM80 proc~pnm80->proc~nutm80

Contents

Source Code


Source Code

    subroutine NUTM80 ( date1, date2, rmatn )

    implicit none

    real(wp),intent(in) :: date1 !! TDB date (Note 1)
    real(wp),intent(in) :: date2 !! TDB date (Note 1)
    real(wp),dimension(3,3),intent(out) :: rmatn !! nutation matrix

    real(wp) :: dpsi, deps, epsa

    !  Nutation components and mean obliquity.
    call NUT80 ( date1, date2, dpsi, deps )
    epsa = OBL80 ( date1, date2 )

    !  Build the rotation matrix.
    call NUMAT ( epsa, dpsi, deps, rmatn )

    end subroutine NUTM80