MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
TimeSelect.f90
Go to the documentation of this file.
1 !> @brief Specify times for some event to occur.
3 
4  use kindmodule, only: dp, i4b, lgp
5  use constantsmodule, only: dzero, done
7  use errorutilmodule, only: pstop
8  use sortmodule, only: qsort
9 
10  implicit none
11  public :: timeselecttype
12 
13  !> @brief Represents a series of instants at which some event should occur.
14  !!
15  !! Maintains an array of configured times which can be sliced to match e.g.
16  !! the current period & time step. Slicing can be performed manually, with
17  !! the select() routine, or automatically, with the advance() routine, for
18  !! a convenient view onto the applicable subset of the complete time array.
19  !!
20  !! Time selection uses the interval convention (t0, t1] (exclusive lower
21  !! bound, inclusive upper bound). This ensures times at exact time step
22  !! boundaries are captured by the earlier time step without duplication.
23  !!
24  !! Array storage can be expanded manually. Note: array expansion must take
25  !! place before selection; when expand() is called the selection is wiped.
26  !! Alternatively, the extend() routine will automatically expand the array
27  !! and sort it.
28  !!
29  !! Most use cases likely assume a strictly increasing time selection; this
30  !! can be checked with increasing(). Note that the sort() routine does not
31  !! check for duplicates, and should usually be followed by an increasing()
32  !! check before the time selection is used.
33  !<
35  real(dp), allocatable :: times(:)
36  integer(I4B) :: selection(2)
37  contains
38  procedure :: deallocate
39  procedure :: expand
40  procedure :: init
41  procedure :: increasing
42  procedure :: log
43  procedure :: select
44  procedure :: advance
45  procedure :: any
46  procedure :: count
47  procedure :: sort
48  procedure :: extend
49  end type timeselecttype
50 
51 contains
52 
53  !> @brief Deallocate the time selection object.
54  subroutine deallocate (this)
55  class(timeselecttype) :: this
56  deallocate (this%times)
57  end subroutine deallocate
58 
59  !> @brief Expand capacity by the given amount. Resets the current slice.
60  subroutine expand(this, increment)
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)/)
66  end subroutine expand
67 
68  !> @brief Initialize or clear the time selection object.
69  subroutine init(this)
70  class(timeselecttype) :: this
71 
72  if (allocated(this%times)) deallocate (this%times)
73  allocate (this%times(0))
74  this%selection = (/0, 0/)
75  end subroutine
76 
77  !> @brief Determine if times strictly increase.
78  !!
79  !! Returns true if the times array strictly increases,
80  !! as well as if the times array is empty, or not yet
81  !! allocated. Note that this function operates on the
82  !! entire times array, not the current selection. Note
83  !! also that this function conducts exact comparisons;
84  !! deduplication with tolerance must be done manually.
85  !<
86  function increasing(this) result(inc)
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
104  end function increasing
105 
106  !> @brief Show the current time selection, if any.
107  subroutine log(this, iout, verb)
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
121  end subroutine log
122 
123  !> @brief Select times in the interval (t0, t1] (exclusive lower, inclusive upper).
124  !!
125  !! Finds and stores the index of the first time strictly after the start time,
126  !! and of the last time at or before the end time. This interval convention
127  !! ensures that times at exact time step boundaries are captured by the later
128  !! (earlier in simulation time) step without duplication.
129  !!
130  !! Allows filtering the times for e.g. a particular stress period and time step.
131  !! Array indices are assumed to start at 1. If no times are found to fall within
132  !! the selection (i.e. the interval falls entirely between two consecutive times
133  !! or beyond the time range), indices are set to [-1, -1].
134  !!
135  !! The given start and end times are first checked against currently stored
136  !! indices to avoid recalculating them if possible, allowing multiple consuming
137  !! components (e.g., subdomain particle tracking solutions) to share the object
138  !! efficiently, provided all proceed through stress periods and time steps in
139  !! lockstep, i.e. they all solve any given period/step before any will proceed
140  !! to the next.
141  !<
142  subroutine select(this, t0, t1, changed)
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 
201  end subroutine
202 
203  !> @brief Update the selection to the current time step.
204  subroutine advance(this)
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)
227  end subroutine advance
228 
229  !> @brief Check if any times are currently selected.
230  !!
231  !! Indicates whether any times are selected for the
232  !! current time step.
233  !!
234  !! Note that this routine does NOT indicate whether
235  !! the times array has nonzero size; use the size
236  !! intrinsic for that.
237  !<
238  function any(this) result(a)
239  class(timeselecttype) :: this
240  logical(LGP) :: a
241 
242  a = all(this%selection > 0)
243  end function any
244 
245  !> @brief Return the number of times currently selected.
246  !!
247  !! Returns the number of times selected for the current
248  !! time step.
249  !!
250  !! Note that this routine does NOT return the total size
251  !! of the times array; use the size intrinsic for that.
252  !<
253  function count(this) result(n)
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
262  end function count
263 
264  !> @brief Sort the time selection in increasing order.
265  !!
266  !! Note that this routine does NOT remove duplicate times.
267  !! Call increasing() to check for duplicates in the array.
268  !<
269  subroutine sort(this)
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)
276  end subroutine sort
277 
278  !> @brief Extend the time selection with the given array.
279  !!
280  !! This routine sorts the selection after appending the
281  !! elements of the given array, but users should likely
282  !! still call increasing() to check for duplicate times.
283  !<
284  subroutine extend(this, a)
285  class(timeselecttype) :: this
286  real(DP) :: a(:)
287 
288  this%times = [this%times, a]
289  call this%sort()
290  end subroutine extend
291 
292 end module timeselectmodule
subroutine init()
Definition: GridSorting.f90:24
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
This module defines variable data types.
Definition: kind.f90:8
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
Specify times for some event to occur.
Definition: TimeSelect.f90:2
logical(lgp) function increasing(this)
Determine if times strictly increase.
Definition: TimeSelect.f90:87
subroutine advance(this)
Update the selection to the current time step.
Definition: TimeSelect.f90:205
integer(i4b) function count(this)
Return the number of times currently selected.
Definition: TimeSelect.f90:254
subroutine log(this, iout, verb)
Show the current time selection, if any.
Definition: TimeSelect.f90:108
subroutine deallocate(this)
Deallocate the time selection object.
Definition: TimeSelect.f90:55
logical(lgp) function any(this)
Check if any times are currently selected.
Definition: TimeSelect.f90:239
subroutine sort(this)
Sort the time selection in increasing order.
Definition: TimeSelect.f90:270
subroutine expand(this, increment)
Expand capacity by the given amount. Resets the current slice.
Definition: TimeSelect.f90:61
subroutine select(this, t0, t1, changed)
Select times in the interval (t0, t1] (exclusive lower, inclusive upper).
Definition: TimeSelect.f90:143
subroutine extend(this, a)
Extend the time selection with the given array.
Definition: TimeSelect.f90:285
Represents a series of instants at which some event should occur.
Definition: TimeSelect.f90:34