Skip to content

Commit

Permalink
Merge branch 'framework/i8_interval_division' into develop
Browse files Browse the repository at this point in the history
This merge introduces support at several levels for interval division resulting
in 64-bit quotients. At the lowest level, a new interface for multiplying time
intervals by 64-bit integers has been added to the mpas_timekeeping module.
Building on this, the mpas_interval_division routine now returns quotients as
64-bit integers. These changes make it possible, e.g., to divide time intervals
of many decades to centuries or more by intervals on the order of one second.

* framework/i8_interval_division:
  Add test for 64-bit interval division quotients to testing core
  Enable 64-bit integer quotients in mpas_interval_division
  Add mpas_timekeeping interface for multiplying an interval by an 8-byte integer
  • Loading branch information
mgduda committed Sep 20, 2017
2 parents 6dac06c + 3f09cf9 commit b632938
Show file tree
Hide file tree
Showing 3 changed files with 45 additions and 22 deletions.
29 changes: 19 additions & 10 deletions src/core_test/mpas_test_core_timekeeping_tests.F
Original file line number Diff line number Diff line change
Expand Up @@ -56,56 +56,63 @@ subroutine test_core_test_intervals(domain, threadErrs, err)!{{{
call mpas_log_write(' Performing time interval tests')

call mpas_log_write(' Test 1:')
call test_core_interval_test('0001-01-01_00:00:00', '0000-01-00_10:00:00', '0001_00:00:00', 31, '0000_10:00:00', err_tmp)
call test_core_interval_test('0001-01-01_00:00:00', '0000-01-00_10:00:00', '0001_00:00:00', 31_I8KIND, '0000_10:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if

call mpas_log_write(' Test 2:')
call test_core_interval_test('0001-01-01_00:00:00', '0000-01-00_00:00:00', '0001_00:00:00', 31, '0000_00:00:00', err_tmp)
call test_core_interval_test('0001-01-01_00:00:00', '0000-01-00_00:00:00', '0001_00:00:00', 31_I8KIND, '0000_00:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if

call mpas_log_write(' Test 3:')
call test_core_interval_test('0001-02-01_00:00:00', '0000-01-00_10:00:00', '0001_00:00:00', 28, '0000_10:00:00', err_tmp)
call test_core_interval_test('0001-02-01_00:00:00', '0000-01-00_10:00:00', '0001_00:00:00', 28_I8KIND, '0000_10:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if

call mpas_log_write(' Test 4:')
call test_core_interval_test('0001-02-01_00:00:00', '0000-01-00_00:00:00', '0001_00:00:00', 28, '0000_00:00:00', err_tmp)
call test_core_interval_test('0001-02-01_00:00:00', '0000-01-00_00:00:00', '0001_00:00:00', 28_I8KIND, '0000_00:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if

call mpas_log_write(' Test 5:')
call test_core_interval_test('0001-01-01_00:00:00', '0000-00-00_01:00:00', '0000_00:30:00', 2, '0000_00:00:00', err_tmp)
call test_core_interval_test('0001-01-01_00:00:00', '0000-00-00_01:00:00', '0000_00:30:00', 2_I8KIND, '0000_00:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if

call mpas_log_write(' Test 6:')
call test_core_interval_test('0001-01-01_00:00:00', '0001-01-00_00:00:00', '0001-00-00_00:00:00', 1, '0000-00-31_00:00:00', err_tmp)
call test_core_interval_test('0001-01-01_00:00:00', '0001-01-00_00:00:00', '0001-00-00_00:00:00', 1_I8KIND, '0000-00-31_00:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if


call mpas_log_write(' Test 7:')
call test_core_interval_test('0000-01-01_00:00:00', '1850-00-00_00:00:00', '00:00:01', 58341600000_I8KIND, '0000-00-00_00:00:00', err_tmp)
if ( err_tmp == 0 ) then
call mpas_log_write(' Result: PASSED')
else
call mpas_log_write(' * Result: FAILED', MPAS_LOG_ERR)
end if

call mpas_log_write(' Completed time interval tests')

end if

call mpas_timer_stop('timekeeping tests')
Expand All @@ -114,13 +121,14 @@ end subroutine test_core_test_intervals!}}}

subroutine test_core_interval_test(ref_str, int1_str, int2_str, expected_divs, expected_remainder_str, ierr)!{{{
character (len=*), intent(in) :: ref_str, int1_str, int2_str
integer, intent(in) :: expected_divs
integer (kind=I8KIND), intent(in) :: expected_divs
character (len=*), intent(in) :: expected_remainder_str
integer, intent(out) :: ierr

integer :: divs
integer (kind=I8KIND) :: divs

character (len=StrKIND) :: remainder_str
character (len=StrKIND) :: temp_str

type (mpas_time_type) :: ref_time
type (mpas_timeinterval_type) :: int1, int2, remainder
Expand Down Expand Up @@ -148,7 +156,8 @@ subroutine test_core_interval_test(ref_str, int1_str, int2_str, expected_divs, e
call mpas_get_timeinterval(remainder, startTimeIn=ref_time, timeString=remainder_str)

call mpas_log_write(' Interval Division summary')
call mpas_log_write(' Divisions: $i', intArgs=(/divs/))
write(temp_str,*) ' Divisions: ', divs
call mpas_log_write(trim(temp_str))
call mpas_log_write(' Remainder: ' // trim(remainder_str))
call mpas_log_write('')

Expand Down
3 changes: 2 additions & 1 deletion src/framework/mpas_stream_manager.F
Original file line number Diff line number Diff line change
Expand Up @@ -3892,7 +3892,8 @@ subroutine mpas_build_stream_filename(ref_time, when, filename_interval, filenam
character(len=StrKIND) :: when_string
type (MPAS_Time_type) :: filetime
type (MPAS_TimeInterval_type) :: intv, rem, zeroIntv
integer :: nrecs, nfiles, irec, direction
integer (kind=I8KIND) :: nrecs
integer :: nfiles, irec, direction
logical :: in_future
STREAM_DEBUG_WRITE(' ** Building Filename')
Expand Down
35 changes: 24 additions & 11 deletions src/framework/mpas_timekeeping.F
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,7 @@ module mpas_timekeeping

interface operator (*)
module procedure mul_ti_n
module procedure mul_ti_n8
end interface

interface operator (/)
Expand Down Expand Up @@ -872,7 +873,7 @@ subroutine mpas_reset_clock_alarm(clock, alarmID, interval, ierr)
type (MPAS_Alarm_type), pointer :: alarmPtr

type (MPAS_TimeInterval_type) :: nowInterval, nowRemainder
integer :: nDivs
integer (kind=I8KIND) :: nDivs
integer :: threadNum

threadNum = mpas_threading_get_thread_num()
Expand Down Expand Up @@ -1569,6 +1570,18 @@ type (MPAS_TimeInterval_type) function mul_ti_n(ti, n)
end function mul_ti_n


type (MPAS_TimeInterval_type) function mul_ti_n8(ti, n8)

implicit none

type (MPAS_TimeInterval_type), intent(in) :: ti
integer (kind=I8KIND), intent(in) :: n8

mul_ti_n8 % ti = ti % ti * n8

end function mul_ti_n8


type (MPAS_TimeInterval_type) function div_ti_n(ti, n)

implicit none
Expand Down Expand Up @@ -1603,10 +1616,10 @@ subroutine mpas_interval_division(ref_time, num, den, n, rem)
type (MPAS_Time_type), intent(in) :: ref_time
type (MPAS_TimeInterval_type), intent(in) :: num
type (MPAS_TimeInterval_type), intent(in) :: den
integer, intent(out) :: n
integer (kind=I8KIND), intent(out) :: n
type (MPAS_TimeInterval_type), intent(out) :: rem

integer :: m
integer (kind=I8KIND) :: m
type (MPAS_Time_type) :: target_time
type (MPAS_Time_type) :: updated_time

Expand All @@ -1621,7 +1634,7 @@ subroutine mpas_interval_division(ref_time, num, den, n, rem)
!
if (den == zero) then
call mpas_log_write('mpas_interval_division: Attempting to divide by zero.', MPAS_LOG_WARN)
n = 0
n = 0_I8KIND
rem = zero
return
end if
Expand All @@ -1634,7 +1647,7 @@ subroutine mpas_interval_division(ref_time, num, den, n, rem)
! numerator == denominator
!
if (target_time == updated_time) then
n = 1
n = 1_I8KIND
rem = zero
return
end if
Expand All @@ -1644,37 +1657,37 @@ subroutine mpas_interval_division(ref_time, num, den, n, rem)
! denominator > numerator
!
if ( target_time < updated_time ) then
n = 0
n = 0_I8KIND
rem = num
return
end if

! One interval of den already fits into num
n = 1
n = 1_I8KIND

! Search forward, doubling the interval each time
do while (target_time > updated_time)
n = n * 2
n = n * 2_I8KIND
temp = den * n
updated_time = ref_time + temp
end do

n = n / 2
n = n / 2_I8KIND
m = n
temp = den * n
updated_time = ref_time + temp


! Search backward, halving the interval each time
do while ( m > 0 )
do while ( m > 0_I8KIND )
temp = den * m

if ( updated_time + temp <= target_time ) then
updated_time = updated_time + temp
n = n + m
end if

m = m / 2
m = m / 2_I8KIND
end do

rem = target_time - updated_time
Expand Down

0 comments on commit b632938

Please sign in to comment.