MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
releaseschedulemodule Module Reference

Particle release scheduling.

Data Types

type  releasescheduletype
 Particle release scheduling utility. More...
 

Functions/Subroutines

type(releasescheduletype) function, pointer, public create_release_schedule (tol)
 Create a new release schedule object. More...
 
subroutine deallocate (this)
 Deallocate the release schedule. More...
 
subroutine log (this, iout)
 Write the release schedule to the given output unit. More...
 
subroutine schedule (this, trelease)
 Add a release time to the schedule. More...
 
subroutine advance (this, lines)
 Refresh the schedule for the current time step. More...
 
logical function any (this)
 Check if any releases are scheduled. More...
 
integer function count (this)
 Return the number of releases scheduled. More...
 

Function/Subroutine Documentation

◆ advance()

subroutine releaseschedulemodule::advance ( class(releasescheduletype), intent(inout)  this,
character(len=linelength), dimension(:), intent(in), optional  lines 
)
private

This involves several tasks: first, advance the time selection. Then, if period-block release setting lines are provided, reinitialize the time step selection for the given period. Finally, refresh the schedule array, deduplicating any times closer than the set tolerance.

This routine is idempotent.

Definition at line 110 of file ReleaseSchedule.f90.

111  use tdismodule, only: totimc, kstp, endofperiod
112  class(ReleaseScheduleType), intent(inout) :: this
113  character(len=LINELENGTH), intent(in), optional :: lines(:)
114  integer(I4B) :: it, i
115  real(DP) :: tprevious
116  real(DP) :: trelease
117 
118  ! Advance the time selection.
119  call this%time_select%advance()
120 
121  ! Reinitialize the time step selection if new
122  ! period-block release settings are provided.
123  if (present(lines)) then
124  call this%step_select%init()
125  do i = 1, size(lines)
126  call this%step_select%read(lines(i))
127  end do
128  end if
129 
130  ! Reinitialize the release time schedule.
131  if (allocated(this%times)) deallocate (this%times)
132  allocate (this%times(0))
133 
134  tprevious = -done
135  trelease = -done
136 
137  ! Add a release time configured by period-block
138  ! settings, if one is scheduled this time step.
139  if (this%step_select%is_selected(kstp, endofperiod=endofperiod)) then
140  trelease = totimc
141  call this%schedule(trelease)
142  tprevious = trelease
143  end if
144 
145  ! Add explicitly configured release times, up
146  ! to the configured tolerance of coincidence.
147  if (this%time_select%any()) then
148  do it = this%time_select%selection(1), this%time_select%selection(2)
149  trelease = this%time_select%times(it)
150  ! Skip the release time if it's too close
151  ! to the previously scheduled release time.
152  if (tprevious >= dzero .and. is_close( &
153  tprevious, &
154  trelease, &
155  atol=this%tolerance)) cycle
156  call this%schedule(trelease)
157  tprevious = trelease
158  end do
159  end if
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
Here is the call graph for this function:

◆ any()

logical function releaseschedulemodule::any ( class(releasescheduletype this)

Note: be sure to call advance() before calling this function, or the result may still be associated with a prior time step.

Definition at line 167 of file ReleaseSchedule.f90.

168  class(ReleaseScheduleType) :: this
169  a = this%count() > 0

◆ count()

integer function releaseschedulemodule::count ( class(releasescheduletype this)
private

Note: be sure to call advance() before calling this function, or the result may still be associated with a prior time step.

Definition at line 177 of file ReleaseSchedule.f90.

178  class(ReleaseScheduleType) :: this
179  n = size(this%times)

◆ create_release_schedule()

type(releasescheduletype) function, pointer, public releaseschedulemodule::create_release_schedule ( real(dp), intent(in)  tol)
Parameters
[in]tolcoincident release time tolerance
Returns
schedule pointer

Definition at line 46 of file ReleaseSchedule.f90.

47  real(DP), intent(in) :: tol !< coincident release time tolerance
48  type(ReleaseScheduleType), pointer :: sched !< schedule pointer
49 
50  allocate (sched)
51  allocate (sched%times(0))
52  allocate (sched%time_select)
53  allocate (sched%step_select)
54  call sched%time_select%init()
55  call sched%step_select%init()
56  sched%tolerance = tol
57 
Here is the caller graph for this function:

◆ deallocate()

subroutine releaseschedulemodule::deallocate ( class(releasescheduletype), intent(inout)  this)
private
Parameters
[in,out]thisthis instance

Definition at line 61 of file ReleaseSchedule.f90.

62  class(ReleaseScheduleType), intent(inout) :: this !< this instance
63 
64  deallocate (this%times)
65  call this%time_select%deallocate()
66  call this%step_select%deallocate()
67  deallocate (this%time_select)
68  deallocate (this%step_select)
69 

◆ log()

subroutine releaseschedulemodule::log ( class(releasescheduletype), intent(inout)  this,
integer(i4b), intent(in)  iout 
)
private
Parameters
[in,out]thisthis instance
[in]ioutoutput unit

Definition at line 73 of file ReleaseSchedule.f90.

74  class(ReleaseScheduleType), intent(inout) :: this !< this instance
75  integer(I4B), intent(in) :: iout !< output unit
76  character(len=*), parameter :: fmt = &
77  &"(6x,A,': ',50(G0,' '))"
78 
79  if (this%any()) then
80  write (iout, fmt) 'RELEASE SCHEDULE', this%times
81  else
82  write (iout, "(1x,a,1x,a)") 'NO RELEASES SCHEDULED'
83  end if

◆ schedule()

subroutine releaseschedulemodule::schedule ( class(releasescheduletype), intent(inout)  this,
real(dp), intent(in)  trelease 
)
private

To schedule multiple release times at once, expand and populate the time selection object by hand. DO NOT attempt to manipulate the times array; this is a read-only property which the schedule maintains.

Definition at line 93 of file ReleaseSchedule.f90.

94  class(ReleaseScheduleType), intent(inout) :: this
95  real(DP), intent(in) :: trelease
96  call expandarray(this%times)
97  this%times(size(this%times)) = trelease