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

Force a read of 'forcing' streams on init #10

Merged
merged 2 commits into from
Jul 15, 2021
Merged
Changes from all commits
Commits
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
55 changes: 54 additions & 1 deletion components/mpas-albany-landice/src/mode_forward/mpas_li_core.F
Original file line number Diff line number Diff line change
Expand Up @@ -112,7 +112,9 @@ function li_core_init(domain, startTimeStamp) result(err)
real (kind=RKIND), pointer :: daysSinceStart
integer, dimension(:), pointer :: cellProcID
type (MPAS_Time_type) :: xtime_timeType, simulationStartTime_timeType

type (MPAS_stream_list_type), pointer :: stream_cursor
type (MPAS_Alarm_type), pointer :: alarm_cursor
character (len=StrKIND) :: actualWhen

err = 0
err_tmp = 0
Expand All @@ -127,9 +129,11 @@ function li_core_init(domain, startTimeStamp) result(err)
err = ior(err, err_tmp)

if (config_do_restart) then
call mpas_log_write("This is a restart: read stream 'restart'.")
call mpas_stream_mgr_read(domain % streamManager, streamID='restart', ierr=err_tmp)
err = ior(err, err_tmp)
else
call mpas_log_write("This is not a restart: read stream 'input'.")
call mpas_stream_mgr_read(domain % streamManager, streamID='input', ierr=err_tmp)
err = ior(err, err_tmp)
end if
Expand All @@ -143,6 +147,7 @@ function li_core_init(domain, startTimeStamp) result(err)
!
! Read the remaining input streams
!
call mpas_log_write("Looking for other input streams set for 'initial_only'.")
call mpas_timer_start('io_read', .false.)
call mpas_stream_mgr_read(domain % streamManager, ierr=err_tmp)
err = ior(err, err_tmp)
Expand All @@ -151,6 +156,54 @@ function li_core_init(domain, startTimeStamp) result(err)
call mpas_stream_mgr_reset_alarms(domain % streamManager, direction=MPAS_STREAM_INPUT, ierr=err_tmp)
err = ior(err, err_tmp)
call mpas_timer_stop('reset_io_alarms')
call mpas_log_write("Finished reading 'initial_only' input streams.")

!
! Read time-varying inputs, if present (i.e., forcing)
!
! Note: The stream manager will not read streams with input_interval set to a time interval on init
! But in general, we want that behavior - we do not want to have to wait a whole time step (or input_interval)
! before having forcing applied. Therefore, we must force those streams to be read here.
call mpas_timer_start('io_read', .false.)
call mpas_log_write("Looking for recurring input streams (forcing) that should be forced to be read at the initial time.")
stream_cursor => domain % streamManager % streams % head
do while (associated(stream_cursor))
if (stream_cursor % direction == MPAS_STREAM_INPUT) then
! Only consider streams of direction input (don't consider input/output streams)
if ((trim(stream_cursor % name) == 'input') .or. &
(trim(stream_cursor % name) == 'restart')) then
! Don't attempt to interact with these two special streams
! (even though restart is not type input, including to be safe)
stream_cursor => stream_cursor % next
cycle
endif
! Only force a read of streams that are recurring. To find out if a stream is recurring, we need to
! get its input alarm. To find the alarm, loop over the alarms in the stream manager's clock until we find it.
alarm_cursor => domain % streamManager % streamClock % alarmListHead
do while (associated(alarm_cursor))
if (trim(alarm_cursor % alarmID) == (trim(stream_cursor % name)//'_input')) then
! We found the stream we are looking for
if (alarm_cursor % isRecurring) then
! Force a read of this stream, and use the latest before time available in the file
call mpas_stream_mgr_read(domain % streamManager, streamID = stream_cursor % name, rightNow = .true., &
whence = MPAS_STREAM_LATEST_BEFORE, actualWhen=actualWhen, ierr=err_tmp)
call mpas_log_write(" * Forced an initial read of input stream '" // trim(stream_cursor%name) // &
"' from time: " // trim(actualWhen))
err = ior(err, err_tmp)
endif
exit ! skip the rest of this loop - we processed the alarm we were looking for
endif
alarm_cursor => alarm_cursor % next
end do
end if ! if stream direction is INPUT
stream_cursor => stream_cursor % next
end do
call mpas_timer_stop('io_read')
call mpas_timer_start('reset_io_alarms', .false.)
call mpas_stream_mgr_reset_alarms(domain % streamManager, direction=MPAS_STREAM_INPUT, ierr=err_tmp)
err = ior(err, err_tmp)
call mpas_timer_stop('reset_io_alarms')
call mpas_log_write("Finished processing recurring input streams at initial time.")

! ===
! Initialize some time stuff on each block
Expand Down