Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Feature/wdboggs/mapl3 freq aspect ts reftime 3358 #3404

Merged
Merged
Show file tree
Hide file tree
Changes from 5 commits
Commits
Show all changes
19 commits
Select commit Hold shift + click to select a range
e7dbf6a
Make module private default.
darianboggs Feb 7, 2025
63516f3
Initial impl of freq aspect compatibility
darianboggs Feb 7, 2025
539127d
Successfully compiles
darianboggs Feb 7, 2025
cc2a5bb
Fix failing FrequencyAspect tests
darianboggs Feb 10, 2025
0018085
Change FrequencyAspect test to different reference times
darianboggs Feb 10, 2025
c9de6c1
Facilitate merge with release/MAPLv3
darianboggs Feb 11, 2025
cf6adc6
Merge branch 'release/MAPL-v3' into feature/wdboggs/mapl3_freq_aspect…
darianboggs Feb 11, 2025
8754b1e
Use setters in constructor; update set_reference_time()
darianboggs Feb 11, 2025
7f4b271
Update test to use updated constructor
darianboggs Feb 11, 2025
7465456
Update fail test
darianboggs Feb 11, 2025
063e97d
Merge branch 'release/MAPL-v3' into feature/wdboggs/mapl3_freq_aspect…
darianboggs Feb 11, 2025
d026f41
Remove 2 setters in constructor
darianboggs Feb 11, 2025
89cf54e
Move contents of check_compatibility into supports_conversion_specific
darianboggs Feb 11, 2025
d200d8b
Clean up argument names
darianboggs Feb 11, 2025
4246c45
Correct check on accumulation type: should be src
darianboggs Feb 12, 2025
bf844f6
Hand merge release/MAPL-v3 changes in
darianboggs Feb 12, 2025
64c1ef9
Merge branch 'release/MAPL-v3' into feature/wdboggs/mapl3_freq_aspect…
darianboggs Feb 12, 2025
9f3c3e1
Reverse src and dst in times_and_intervals_are_compatible
darianboggs Feb 13, 2025
ef99af4
Merge branch 'release/MAPL-v3' into feature/wdboggs/mapl3_freq_aspect…
darianboggs Feb 13, 2025
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion esmf_utils/ESMF_Time_Utilities.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,7 @@ module mapl3g_ESMF_Time_Utilities
use esmf, I4 => ESMF_KIND_I4
use mapl_ErrorHandling
implicit none (type, external)
! private !wdb fixme deleteme should this be private
private

public :: zero_time_interval
public :: intervals_are_compatible
Expand Down
93 changes: 49 additions & 44 deletions generic3g/specs/FrequencyAspect.F90
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@ module mapl3g_FrequencyAspect
use mapl3g_AspectId
use mapl3g_StateItemAspect
use mapl3g_AccumulatorActionInterface
use mapl3g_ESMF_Time_Utilities, only: times_and_intervals_are_compatible, zero_time_interval
use esmf
implicit none
private
Expand All @@ -14,6 +15,7 @@ module mapl3g_FrequencyAspect
type, extends(StateItemAspect) :: FrequencyAspect
private
type(ESMF_TimeInterval) :: timestep_
type(ESMF_Time) :: reference_time_
character(len=:), allocatable :: accumulation_type_
contains
! These are implementations of extended derived type.
Expand All @@ -28,22 +30,18 @@ module mapl3g_FrequencyAspect
procedure :: set_timestep
procedure :: get_accumulation_type
procedure :: set_accumulation_type
procedure :: get_reference_time
procedure :: set_reference_time
procedure, private :: zero_timestep
end type FrequencyAspect

interface FrequencyAspect
module procedure :: construct_frequency_aspect
end interface FrequencyAspect

interface operator(.divides.)
module procedure :: aspect_divides
end interface operator(.divides.)

! This value should not be accessed directly. Use get_zero() instead.
! There is no constructor for ESMF_TimeInterval, so the value cannot be initialized
! at construction. The get_zero() function initializes the value the first time
! and returns a pointer to the value.
type(ESMF_TimeInterval), target :: ZERO_TI
interface check_compatibility
module procedure :: check_freq_aspect_compatibility
end interface check_compatibility

contains

Expand Down Expand Up @@ -77,6 +75,22 @@ subroutine set_timestep(this, timestep)

end subroutine set_timestep

function get_reference_time(this) result(rf)
type(ESMF_Time) :: rf
class(FrequencyAspect), intent(in) :: this

rf = this%reference_time_

end function get_reference_time

subroutine set_reference_time(this, reference_time)
class(FrequencyAspect), intent(inout) :: this
type(ESMF_Time), intent(in) :: reference_time

this%reference_time_ = reference_time

end subroutine set_reference_time

subroutine zero_timestep(this)
class(FrequencyAspect), intent(inout) :: this

Expand Down Expand Up @@ -110,15 +124,16 @@ logical function matches(src, dst) result(does_match)
type(ESMF_TimeInterval), pointer :: zero

does_match = .TRUE.
zero => get_zero()
zero => zero_time_interval()
src_timestep = src%get_timestep()
if(src_timestep == zero) return
select type(dst)
class is (FrequencyAspect)
dst_timestep = dst%get_timestep()
if(dst_timestep == zero) return
if(.not. accumulation_type_is_valid(dst%get_accumulation_type())) return
does_match = dst_timestep == src_timestep
does_match = dst_timestep == src_timestep .and. &
& src%get_reference_time() == dst%get_reference_time()
end select

end function matches
Expand Down Expand Up @@ -167,50 +182,40 @@ end function supports_conversion_general
logical function supports_conversion_specific(src, dst) result(supports)
class(FrequencyAspect), intent(in) :: src
class(StateItemAspect), intent(in) :: dst
integer :: status
logical :: compatible

select type(dst)
class is (FrequencyAspect)
supports = src .divides. dst
call check_compatibility(dst, src, compatible, rc=status)
supports = compatible .and. status == _SUCCESS
end select

end function supports_conversion_specific

logical function aspect_divides(factor, base)
class(FrequencyAspect), intent(in) :: factor
class(FrequencyAspect), intent(in) :: base

aspect_divides = interval_divides(factor%get_timestep(), base%get_timestep())

end function aspect_divides

logical function interval_divides(factor, base) result(lval)
type(ESMF_TimeInterval), intent(in) :: factor
type(ESMF_TimeInterval), intent(in) :: base
type(ESMF_TimeInterval), pointer :: zero

lval = .FALSE.
zero => get_zero()
if(factor == zero) return
lval = mod(base, factor) == zero

end function interval_divides

function get_zero() result(zero)
type(ESMF_TimeInterval), pointer :: zero
logical, save :: zero_is_uninitialized = .TRUE.

if(zero_is_uninitialized) then
call ESMF_TimeIntervalSet(ZERO_TI, ns=0)
zero_is_uninitialized = .FALSE.
end if
zero => ZERO_TI

end function get_zero

function get_aspect_id() result(aspect_id)
type(AspectId) :: aspect_id
aspect_id = FREQUENCY_ASPECT_ID
end function get_aspect_id

subroutine check_freq_aspect_compatibility(child, parent, compatible, rc)
class(FrequencyAspect), intent(in) :: child
class(FrequencyAspect), intent(in) :: parent
logical, intent(out) :: compatible
integer, optional, intent(out) :: rc
integer :: status
type(ESMF_TimeInterval) :: child_step, parent_step
type(ESMF_Time) :: child_reference, parent_reference

child_step = child%get_timestep()
child_reference = child%get_reference_time()
parent_step = parent%get_timestep()
parent_reference = parent%get_reference_time()

call times_and_intervals_are_compatible(child_step, child_reference, &
& parent_step, parent_reference, compatible, _RC)
_RETURN(_SUCCESS)

end subroutine check_freq_aspect_compatibility

end module mapl3g_FrequencyAspect
8 changes: 8 additions & 0 deletions generic3g/tests/Test_Aspects.pf
Original file line number Diff line number Diff line change
Expand Up @@ -285,12 +285,20 @@ contains
type(FrequencyAspect) :: import, export

type(ESMF_TimeInterval) :: dt1, dt2
type(ESMF_Time) :: ref1, ref2

integer :: status
call ESMF_TimeIntervalSet(dt1, s=4)
call ESMF_TimeIntervalSet(dt2, s=2) ! commensurate

import = FrequencyAspect(dt1, 'mean')
call ESMF_TimeSet(ref1, s=0, rc=status)
@assertEqual(0, status, 'Nonzero status')
call import%set_reference_time(ref1)
export = FrequencyAspect(dt2)
call ESMF_TimeSet(ref2, m=0, rc=status)
@assertEqual(0, status, 'Nonzero status')
call export%set_reference_time(ref2)

@assert_that(export%can_connect_to(import), is(true()))

Expand Down