MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
timeselectmodule Module Reference

Specify times for some event to occur.

Data Types

type  timeselecttype
 Represents a series of instants at which some event should occur. More...
 

Functions/Subroutines

subroutine deallocate (this)
 Deallocate the time selection object. More...
 
subroutine expand (this, increment)
 Expand capacity by the given amount. Resets the current slice. More...
 
subroutine init (this)
 Initialize or clear the time selection object. More...
 
logical(lgp) function increasing (this)
 Determine if times strictly increase. More...
 
subroutine log (this, iout, verb)
 Show the current time selection, if any. More...
 
subroutine select (this, t0, t1, changed)
 Select times in the interval (t0, t1] (exclusive lower, inclusive upper). More...
 
subroutine advance (this)
 Update the selection to the current time step. More...
 
logical(lgp) function any (this)
 Check if any times are currently selected. More...
 
integer(i4b) function count (this)
 Return the number of times currently selected. More...
 
subroutine sort (this)
 Sort the time selection in increasing order. More...
 
subroutine extend (this, a)
 Extend the time selection with the given array. More...
 

Function/Subroutine Documentation

◆ advance()

subroutine timeselectmodule::advance ( class(timeselecttype this)

Definition at line 204 of file TimeSelect.f90.

205  ! modules
206  use tdismodule, only: kper, kstp, nper, nstp, totimc, delt
207  ! dummy
208  class(TimeSelectType) :: this
209  ! local
210  real(DP) :: l, u
211 
212  if (kper == 1 .and. kstp == 1) then
213  ! For first time step, use a small negative lower bound
214  ! capture times at t=0.0 despite exclusive lower bound
215  l = -epsilon(dzero)
216  else
217  l = totimc
218  end if
219  if (kper == nper .and. kstp == nstp(kper)) then
220  ! For last time step, use a large upper bound to
221  ! capture times beyond the end of the simulation
222  u = huge(done)
223  else
224  u = totimc + delt
225  end if
226  call this%select(l, u)
integer(i4b), dimension(:), pointer, public, contiguous nstp
number of time steps in each stress period
Definition: tdis.f90:39
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
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21

◆ any()

logical(lgp) function timeselectmodule::any ( class(timeselecttype this)

Indicates whether any times are selected for the current time step.

Note that this routine does NOT indicate whether the times array has nonzero size; use the size intrinsic for that.

Definition at line 238 of file TimeSelect.f90.

239  class(TimeSelectType) :: this
240  logical(LGP) :: a
241 
242  a = all(this%selection > 0)

◆ count()

integer(i4b) function timeselectmodule::count ( class(timeselecttype this)

Returns the number of times selected for the current time step.

Note that this routine does NOT return the total size of the times array; use the size intrinsic for that.

Definition at line 253 of file TimeSelect.f90.

254  class(TimeSelectType) :: this
255  integer(I4B) :: n
256 
257  if (this%any()) then
258  n = this%selection(2) - this%selection(1)
259  else
260  n = 0
261  end if

◆ deallocate()

subroutine timeselectmodule::deallocate ( class(timeselecttype this)

Definition at line 54 of file TimeSelect.f90.

55  class(TimeSelectType) :: this
56  deallocate (this%times)

◆ expand()

subroutine timeselectmodule::expand ( class(timeselecttype this,
integer(i4b), intent(in), optional  increment 
)

Definition at line 60 of file TimeSelect.f90.

61  class(TimeSelectType) :: this
62  integer(I4B), optional, intent(in) :: increment
63 
64  call expandarray(this%times, increment=increment)
65  this%selection = (/1, size(this%times)/)

◆ extend()

subroutine timeselectmodule::extend ( class(timeselecttype this,
real(dp), dimension(:)  a 
)

This routine sorts the selection after appending the elements of the given array, but users should likely still call increasing() to check for duplicate times.

Definition at line 284 of file TimeSelect.f90.

285  class(TimeSelectType) :: this
286  real(DP) :: a(:)
287 
288  this%times = [this%times, a]
289  call this%sort()

◆ increasing()

logical(lgp) function timeselectmodule::increasing ( class(timeselecttype this)

Returns true if the times array strictly increases, as well as if the times array is empty, or not yet allocated. Note that this function operates on the entire times array, not the current selection. Note also that this function conducts exact comparisons; deduplication with tolerance must be done manually.

Definition at line 86 of file TimeSelect.f90.

87  class(TimeSelectType) :: this
88  logical(LGP) :: inc
89  integer(I4B) :: i
90  real(DP) :: l, t
91 
92  inc = .true.
93  if (.not. allocated(this%times)) return
94  do i = 1, size(this%times)
95  t = this%times(i)
96  if (i /= 1) then
97  if (l >= t) then
98  inc = .false.
99  return
100  end if
101  end if
102  l = t
103  end do

◆ init()

subroutine timeselectmodule::init ( class(timeselecttype this)

Definition at line 69 of file TimeSelect.f90.

70  class(TimeSelectType) :: this
71 
72  if (allocated(this%times)) deallocate (this%times)
73  allocate (this%times(0))
74  this%selection = (/0, 0/)

◆ log()

subroutine timeselectmodule::log ( class(timeselecttype this,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  verb 
)
Parameters
thisthis instance
[in]ioutoutput unit
[in]verbselection name

Definition at line 107 of file TimeSelect.f90.

108  ! dummy
109  class(TimeSelectType) :: this !< this instance
110  integer(I4B), intent(in) :: iout !< output unit
111  character(len=*), intent(in) :: verb !< selection name
112  ! formats
113  character(len=*), parameter :: fmt = &
114  &"(6x,'THE FOLLOWING TIMES WILL BE ',A,': ',50(G0,' '))"
115 
116  if (this%any()) then
117  write (iout, fmt) verb, this%times(this%selection(1):this%selection(2))
118  else
119  write (iout, "(a,1x,a)") 'NO TIMES WILL BE', verb
120  end if

◆ select()

subroutine timeselectmodule::select ( class(timeselecttype this,
real(dp), intent(in)  t0,
real(dp), intent(in)  t1,
logical(lgp), intent(inout), optional  changed 
)

Finds and stores the index of the first time strictly after the start time, and of the last time at or before the end time. This interval convention ensures that times at exact time step boundaries are captured by the later (earlier in simulation time) step without duplication.

Allows filtering the times for e.g. a particular stress period and time step. Array indices are assumed to start at 1. If no times are found to fall within the selection (i.e. the interval falls entirely between two consecutive times or beyond the time range), indices are set to [-1, -1].

The given start and end times are first checked against currently stored indices to avoid recalculating them if possible, allowing multiple consuming components (e.g., subdomain particle tracking solutions) to share the object efficiently, provided all proceed through stress periods and time steps in lockstep, i.e. they all solve any given period/step before any will proceed to the next.

Definition at line 142 of file TimeSelect.f90.

143  ! dummy
144  class(TimeSelectType) :: this
145  real(DP), intent(in) :: t0, t1
146  logical(LGP), intent(inout), optional :: changed
147  ! local
148  integer(I4B) :: i, i0, i1
149  integer(I4B) :: l, u, lp, up
150  real(DP) :: t
151 
152  ! by default, need to iterate over all times
153  i0 = 1
154  i1 = size(this%times)
155 
156  ! if no times fall within the slice, set to [-1, -1]
157  l = -1
158  u = -1
159 
160  ! previous bounding indices
161  lp = this%selection(1)
162  up = this%selection(2)
163 
164  ! Check if we can reuse either the lower or upper bound.
165  ! The lower doesn't need to change if it indexes the first
166  ! time strictly after the slice's start (i.e., the previous
167  ! time is at or before t0, and this time is after t0).
168  ! The upper doesn't need to change if it indexes the last
169  ! time at or before the slice's end.
170  if (lp > 0 .and. up > 0) then
171  if (lp > 1) then
172  if (this%times(lp - 1) <= t0 .and. &
173  this%times(lp) > t0) then
174  l = lp
175  i0 = l
176  end if
177  end if
178  if (up > 1 .and. up < i1) then
179  if (this%times(up + 1) > t1 .and. &
180  this%times(up) <= t1) then
181  u = up
182  i1 = u
183  end if
184  end if
185  if (l == lp .and. u == up) then
186  this%selection = (/l, u/)
187  if (present(changed)) changed = .false.
188  return
189  end if
190  end if
191 
192  ! recompute bounding indices if needed
193  do i = i0, i1
194  t = this%times(i)
195  if (l < 0 .and. t > t0 .and. t <= t1) l = i
196  if (l > 0 .and. t <= t1) u = i
197  end do
198  this%selection = (/l, u/)
199  if (present(changed)) changed = l /= lp .or. u /= up
200 

◆ sort()

subroutine timeselectmodule::sort ( class(timeselecttype this)

Note that this routine does NOT remove duplicate times. Call increasing() to check for duplicates in the array.

Definition at line 269 of file TimeSelect.f90.

270  class(TimeSelectType) :: this
271  integer(I4B), allocatable :: indx(:)
272 
273  allocate (indx(size(this%times)))
274  call qsort(indx, this%times)
275  deallocate (indx)