MODFLOW 6  version 6.6.0.dev0
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 between t0 and t1 (inclusive). 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 197 of file TimeSelect.f90.

198  ! modules
199  use tdismodule, only: kper, kstp, nper, nstp, totimc, delt
200  ! dummy
201  class(TimeSelectType) :: this
202  ! local
203  real(DP) :: l, u
204 
205  l = minval(this%times)
206  u = maxval(this%times)
207  if (.not. (kper == 1 .and. kstp == 1)) l = totimc
208  if (.not. (kper == nper .and. kstp == nstp(kper))) u = totimc + delt
209  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 221 of file TimeSelect.f90.

222  class(TimeSelectType) :: this
223  logical(LGP) :: a
224 
225  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 236 of file TimeSelect.f90.

237  class(TimeSelectType) :: this
238  integer(I4B) :: n
239 
240  if (this%any()) then
241  n = this%selection(2) - this%selection(1)
242  else
243  n = 0
244  end if

◆ deallocate()

subroutine timeselectmodule::deallocate ( class(timeselecttype this)

Definition at line 50 of file TimeSelect.f90.

51  class(TimeSelectType) :: this
52  deallocate (this%times)

◆ expand()

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

Definition at line 56 of file TimeSelect.f90.

57  class(TimeSelectType) :: this
58  integer(I4B), optional, intent(in) :: increment
59 
60  call expandarray(this%times, increment=increment)
61  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 267 of file TimeSelect.f90.

268  class(TimeSelectType) :: this
269  real(DP) :: a(:)
270 
271  this%times = [this%times, a]
272  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 82 of file TimeSelect.f90.

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

◆ init()

subroutine timeselectmodule::init ( class(timeselecttype this)

Definition at line 65 of file TimeSelect.f90.

66  class(TimeSelectType) :: this
67 
68  if (allocated(this%times)) deallocate (this%times)
69  allocate (this%times(0))
70  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 103 of file TimeSelect.f90.

104  ! dummy
105  class(TimeSelectType) :: this !< this instance
106  integer(I4B), intent(in) :: iout !< output unit
107  character(len=*), intent(in) :: verb !< selection name
108  ! formats
109  character(len=*), parameter :: fmt = &
110  &"(6x,'THE FOLLOWING TIMES WILL BE ',A,': ',50(G0,' '))"
111 
112  if (this%any()) then
113  write (iout, fmt) verb, this%times(this%selection(1):this%selection(2))
114  else
115  write (iout, "(a,1x,a)") 'NO TIMES WILL BE', verb
116  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 at the same instant as or following the start time, and of the last time at the same instant as or preceding the end time. 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. it 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 136 of file TimeSelect.f90.

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

◆ 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 252 of file TimeSelect.f90.

253  class(TimeSelectType) :: this
254  integer(I4B), allocatable :: indx(:)
255 
256  allocate (indx(size(this%times)))
257  call qsort(indx, this%times)
258  deallocate (indx)