MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
particlereleaseschedulemodule Module Reference

Particle release scheduling.

Data Types

type  particlereleasescheduletype
 Particle release scheduling utility. More...
 

Functions/Subroutines

type(particlereleasescheduletype) function, pointer, public create_release_schedule (tolerance)
 Create a new release schedule. More...
 
subroutine destroy (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 particlereleaseschedulemodule::advance ( class(particlereleasescheduletype), 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 108 of file ParticleReleaseSchedule.f90.

109  use tdismodule, only: totimc, kstp, endofperiod
110  class(ParticleReleaseScheduleType), intent(inout) :: this
111  character(len=LINELENGTH), intent(in), optional :: lines(:)
112  integer(I4B) :: it, i
113  real(DP) :: tprevious
114  real(DP) :: trelease
115 
116  ! Advance the time selection.
117  call this%time_select%advance()
118 
119  ! Reinitialize the time step selection if new
120  ! period-block release settings are provided.
121  if (present(lines)) then
122  call this%step_select%init()
123  do i = 1, size(lines)
124  call this%step_select%read(lines(i))
125  end do
126  end if
127 
128  ! Reinitialize the release time schedule.
129  if (allocated(this%times)) deallocate (this%times)
130  allocate (this%times(0))
131 
132  tprevious = -done
133  trelease = -done
134 
135  ! Add a release time configured by period-block
136  ! settings, if one is scheduled this time step.
137  if (this%step_select%is_selected(kstp, endofperiod=endofperiod)) then
138  trelease = totimc
139  call this%schedule(trelease)
140  tprevious = trelease
141  end if
142 
143  ! Schedule explicitly specified release times, up
144  ! to the configured tolerance of coincidence
145  if (this%time_select%any()) then
146  do it = this%time_select%selection(1), this%time_select%selection(2)
147  trelease = this%time_select%times(it)
148  if (tprevious >= dzero .and. is_close( &
149  tprevious, &
150  trelease, &
151  atol=this%tolerance)) cycle
152 
153  call this%schedule(trelease)
154  tprevious = trelease
155  end do
156  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 particlereleaseschedulemodule::any ( class(particlereleasescheduletype 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 164 of file ParticleReleaseSchedule.f90.

165  class(ParticleReleaseScheduleType) :: this
166  a = this%count() > 0

◆ count()

integer function particlereleaseschedulemodule::count ( class(particlereleasescheduletype 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 174 of file ParticleReleaseSchedule.f90.

175  class(ParticleReleaseScheduleType) :: this
176  n = size(this%times)

◆ create_release_schedule()

type(particlereleasescheduletype) function, pointer, public particlereleaseschedulemodule::create_release_schedule ( real(dp), intent(in)  tolerance)
Parameters
[in]tolerancecoincident release time tolerance
Returns
schedule pointer

Definition at line 46 of file ParticleReleaseSchedule.f90.

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

◆ destroy()

subroutine particlereleaseschedulemodule::destroy ( class(particlereleasescheduletype), intent(inout)  this)
private
Parameters
[in,out]thisthis instance

Definition at line 60 of file ParticleReleaseSchedule.f90.

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

◆ log()

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

Definition at line 71 of file ParticleReleaseSchedule.f90.

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

◆ schedule()

subroutine particlereleaseschedulemodule::schedule ( class(particlereleasescheduletype), 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 91 of file ParticleReleaseSchedule.f90.

92  class(ParticleReleaseScheduleType), intent(inout) :: this
93  real(DP), intent(in) :: trelease
94  call expandarray(this%times)
95  this%times(size(this%times)) = trelease