Decompose days to hours, minutes, seconds, fraction.
Status: vector/matrix support routine.
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...
The largest positive useful value for NDP is determined by the size of DAYS, the format of REAL(WP) floating-point numbers on the target platform, and the risk of overflowing IHMSF(4). On a typical platform, for DAYS up to 1D0, 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 IHMSF(4).
The absolute value of DAYS may exceed 1D0. In cases where it does not, it is up to the caller to test for and handle the case where DAYS is very nearly 1D0 and rounds up to 24 hours, by testing for IHMSF(1)=24 and setting IHMSF(1-4) to zero.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ndp | resolution (Note 1) |
||
real(kind=wp), | intent(in) | :: | days | interval in days |
||
character(len=*), | intent(out) | :: | sign | '+' or '-' |
||
integer, | intent(out), | dimension(4) | :: | ihmsf | hours, minutes, seconds, fraction |
subroutine D2TF ( ndp, days, sign, ihmsf )
implicit none
integer,intent(in) :: ndp !! resolution (Note 1)
real(wp),intent(in) :: days !! interval in days
character(len=*),intent(out) :: sign !! '+' or '-'
integer,dimension(4),intent(out) :: ihmsf !! hours, minutes, seconds, fraction
integer :: nrs, n
real(wp) :: rs, rm, rh, a, ah, am, as, af
! Handle sign.
if ( days >= 0.0_wp ) then
sign = '+'
else
sign = '-'
end if
! Interval in seconds.
a = d2s * abs(days)
! Pre-round if resolution coarser than 1 second (then pretend NDP=1).
if ( ndp < 0 ) then
nrs = 1
do n=1,-ndp
if ( n==2 .or. n==4 ) then
nrs = nrs * 6
else
nrs = nrs * 10
end if
end do
rs = real(nrs, wp)
a = rs * anint(a/rs)
end if
! Express the unit of each field in resolution units.
nrs = 1
do n=1,ndp
nrs = nrs * 10
end do
rs = real(nrs, wp)
rm = rs * 60.0_wp
rh = rm * 60.0_wp
! Round the interval and express in resolution units.
a = anint(rs*a)
! Break into fields.
ah = aint(a/rh)
a = a - ah*rh
am = aint(a/rm)
a = a - am*rm
as = aint(a/rs)
af = a - as*rs
! Return results.
ihmsf(1) = nint(ah)
ihmsf(2) = nint(am)
ihmsf(3) = nint(as)
ihmsf(4) = nint(af)
end subroutine D2TF