A2AF Subroutine

public subroutine A2AF(ndp, angle, sign, idmsf)

Decompose radians into degrees, arcminutes, arcseconds, fraction.

Status: vector/matrix support routine.

Notes

  1. NDP is interpreted as follows:
     NDP         resolution
      :      ...0000 00 00
     -7         1000 00 00
     -6          100 00 00
     -5           10 00 00
     -4            1 00 00
     -3            0 10 00
     -2            0 01 00
     -1            0 00 10
      0            0 00 01
      1            0 00 00.1
      2            0 00 00.01
      3            0 00 00.001
      :            0 00 00.000...
  1. The largest positive useful value for NDP is determined by the size of ANGLE, the format of REAL(WP) floating-point numbers on the target platform, and the risk of overflowing IDMSF(4). On a typical platform, for ANGLE up to 2pi, the available floating-point precision might correspond to NDP=12. However, the practical limit is typically NDP=9, set by the capacity of a 32-bit IDMSF(4).

  2. The absolute value of ANGLE may exceed 2pi. In cases where it does not, it is up to the caller to test for and handle the case where ANGLE is very nearly 2pi and rounds up to 360 degrees, by testing for IDMSF(1)=360 and setting IDMSF(1-4) to zero.

History

  • IAU SOFA revision: 2007 December 3

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: ndp

resolution (Note 1)

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

angle in radians

character(len=*), intent(out) :: sign

'+' or '-'

integer, intent(out), dimension(4):: idmsf

degrees, arcminutes, arcseconds, fraction


Calls

proc~~a2af~~CallsGraph proc~a2af A2AF proc~d2tf D2TF proc~a2af->proc~d2tf

Contents

Source Code


Source Code

    subroutine A2AF ( ndp, angle, sign, idmsf )

    implicit none

    integer,intent(in) :: ndp !! resolution (Note 1)
    real(wp),intent(in) :: angle !! angle in radians
    character(len=*),intent(out) :: sign !! '+' or '-'
    integer,dimension(4),intent(out) :: idmsf !! degrees, arcminutes, arcseconds, fraction

    !  Hours to degrees * radians to turns
    real(wp),parameter :: f = 15.0_wp/d2pi

    !  Scale then use days to h,m,s routine.
    call D2TF ( ndp, angle*f, sign, idmsf )

    end subroutine A2AF