MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
TvBase.f90
Go to the documentation of this file.
1 !> @brief This module contains common time-varying property functionality
2 !!
3 !! This module contains methods implementing functionality common to both
4 !! time-varying hydraulic conductivity (TVK) and time-varying storage (TVS)
5 !! packages.
6 !!
7 !<
9  use basedismodule, only: disbasetype
11  use kindmodule, only: i4b, dp
14  use simvariablesmodule, only: errmsg
15  use tdismodule, only: kper, nper, kstp
19 
20  implicit none
21 
22  private
23 
24  public :: tvbasetype
25  public :: tvbase_da
26 
27  type, abstract, extends(numericalpackagetype) :: tvbasetype
28  type(timeseriesmanagertype), pointer, private :: tsmanager => null()
29  contains
30  procedure :: init
31  procedure :: ar
32  procedure :: rp
33  procedure :: ad
34  procedure :: da => tvbase_da
35  procedure, private :: tvbase_allocate_scalars
36  procedure, private :: read_options
37  procedure(ar_set_pointers), deferred :: ar_set_pointers
38  procedure(read_option), deferred :: read_option
40  procedure(set_changed_at), deferred :: set_changed_at
41  procedure(reset_change_flags), deferred :: reset_change_flags
42  procedure(validate_change), deferred :: validate_change
43  end type tvbasetype
44 
45  abstract interface
46 
47  !> @brief Announce package and set pointers to variables
48  !!
49  !! Deferred procedure called by the TvBaseType code to announce the
50  !! specific package version and set any required array and variable
51  !! pointers from other packages.
52  !!
53  !<
54  subroutine ar_set_pointers(this)
55  ! -- modules
56  import tvbasetype
57  ! -- dummy variables
58  class(tvbasetype) :: this
59  end subroutine
60 
61  !> @brief Announce package and set pointers to variables
62  !!
63  !! Deferred procedure called by the TvBaseType code to process a single
64  !! keyword from the OPTIONS block of the package input file.
65  !!
66  !<
67  function read_option(this, keyword) result(success)
68  ! -- modules
69  import tvbasetype
70  ! -- dummy variables
71  class(tvbasetype) :: this
72  character(len=*), intent(in) :: keyword
73  ! -- return
74  logical :: success
75  end function
76 
77  !> @brief Get an array value pointer given a variable name and node index
78  !!
79  !! Deferred procedure called by the TvBaseType code to retrieve a pointer
80  !! to a given node's value for a given named variable.
81  !!
82  !<
83  function get_pointer_to_value(this, n, varName) result(bndElem)
84  ! -- modules
85  use kindmodule, only: i4b, dp
86  import tvbasetype
87  ! -- dummy variables
88  class(tvbasetype) :: this
89  integer(I4B), intent(in) :: n
90  character(len=*), intent(in) :: varname
91  ! -- return
92  real(dp), pointer :: bndelem
93  end function
94 
95  !> @brief Mark property changes as having occurred at (kper, kstp)
96  !!
97  !! Deferred procedure called by the TvBaseType code when a property value
98  !! change occurs at (kper, kstp).
99  !!
100  !<
101  subroutine set_changed_at(this, kper, kstp)
102  ! -- modules
103  use kindmodule, only: i4b
104  import tvbasetype
105  ! -- dummy variables
106  class(tvbasetype) :: this
107  integer(I4B), intent(in) :: kper
108  integer(I4B), intent(in) :: kstp
109  end subroutine
110 
111  !> @brief Clear all per-node change flags
112  !!
113  !! Deferred procedure called by the TvBaseType code when a new time step
114  !! commences, indicating that any previously set per-node property value
115  !! change flags should be reset.
116  !!
117  !<
118  subroutine reset_change_flags(this)
119  ! -- modules
120  import tvbasetype
121  ! -- dummy variables
122  class(tvbasetype) :: this
123  end subroutine
124 
125  !> @brief Check that a given property value is valid
126  !!
127  !! Deferred procedure called by the TvBaseType code after a property value
128  !! change occurs to perform any required validity checks on the value of
129  !! the given variable at the given node. Perform any required updates to
130  !! the property value if it is valid, or log an error if not.
131  !!
132  !<
133  subroutine validate_change(this, n, varName)
134  ! -- modules
135  use kindmodule, only: i4b
136  import tvbasetype
137  ! -- dummy variables
138  class(tvbasetype) :: this
139  integer(I4B), intent(in) :: n
140  character(len=*), intent(in) :: varName
141  end subroutine
142 
143  end interface
144 
145 contains
146 
147  !> @brief Initialise the TvBaseType object
148  !!
149  !! Allocate and initialize data members of the object.
150  !!
151  !<
152  subroutine init(this, name_model, pakname, ftype, inunit, iout)
153  ! -- dummy variables
154  class(tvbasetype) :: this
155  character(len=*), intent(in) :: name_model
156  character(len=*), intent(in) :: pakname
157  character(len=*), intent(in) :: ftype
158  integer(I4B), intent(in) :: inunit
159  integer(I4B), intent(in) :: iout
160  !
161  call this%set_names(1, name_model, pakname, ftype)
162  call this%tvbase_allocate_scalars()
163  this%inunit = inunit
164  this%iout = iout
165  call this%parser%Initialize(this%inunit, this%iout)
166  end subroutine init
167 
168  !> @brief Allocate scalar variables
169  !!
170  !! Allocate scalar data members of the object.
171  !!
172  !<
173  subroutine tvbase_allocate_scalars(this)
174  ! -- dummy variables
175  class(tvbasetype) :: this
176  !
177  ! -- Call standard NumericalPackageType allocate scalars
178  call this%NumericalPackageType%allocate_scalars()
179  !
180  ! -- Allocate time series manager
181  allocate (this%tsmanager)
182  end subroutine tvbase_allocate_scalars
183 
184  !> @brief Allocate and read method for package
185  !!
186  !! Allocate and read static data for the package.
187  !!
188  !<
189  subroutine ar(this, dis)
190  ! -- dummy variables
191  class(tvbasetype) :: this
192  class(disbasetype), pointer, intent(in) :: dis
193  !
194  ! -- Set pointers to other package variables and announce package
195  this%dis => dis
196  call this%ar_set_pointers()
197  !
198  ! -- Create time series manager
199  call tsmanager_cr(this%tsmanager, &
200  this%iout, &
201  removetslinksoncompletion=.true., &
202  extendtstoendofsimulation=.true.)
203  !
204  ! -- Read options
205  call this%read_options()
206  !
207  ! -- Now that time series will have been read, need to call the df routine
208  ! -- to define the manager
209  call this%tsmanager%tsmanager_df()
210  !
211  ! -- Terminate if any errors were encountered
212  if (count_errors() > 0) then
213  call this%parser%StoreErrorUnit()
214  call ustop()
215  end if
216  end subroutine ar
217 
218  !> @brief Read OPTIONS block of package input file
219  !!
220  !! Reads the OPTIONS block of the package's input file, deferring to the
221  !! derived type to process any package-specific keywords.
222  !!
223  !<
224  subroutine read_options(this)
225  ! -- dummy variables
226  class(tvbasetype) :: this
227  ! -- local variables
228  character(len=LINELENGTH) :: keyword
229  character(len=MAXCHARLEN) :: fname
230  logical :: isfound
231  logical :: endOfBlock
232  integer(I4B) :: ierr
233  ! -- formats
234  character(len=*), parameter :: fmtts = &
235  &"(4x, 'TIME-SERIES DATA WILL BE READ FROM FILE: ', a)"
236  !
237  ! -- Get options block
238  call this%parser%GetBlock('OPTIONS', isfound, ierr, &
239  blockrequired=.false., supportopenclose=.true.)
240  !
241  ! -- Parse options block if detected
242  if (isfound) then
243  write (this%iout, '(1x,a)') &
244  'PROCESSING '//trim(adjustl(this%packName))//' OPTIONS'
245  do
246  call this%parser%GetNextLine(endofblock)
247  if (endofblock) then
248  exit
249  end if
250  call this%parser%GetStringCaps(keyword)
251  select case (keyword)
252  case ('PRINT_INPUT')
253  this%iprpak = 1
254  write (this%iout, '(4x,a)') 'TIME-VARYING INPUT WILL BE PRINTED.'
255  case ('TS6')
256  !
257  ! -- Add a time series file
258  call this%parser%GetStringCaps(keyword)
259  if (trim(adjustl(keyword)) /= 'FILEIN') then
260  errmsg = &
261  'TS6 keyword must be followed by "FILEIN" then by filename.'
262  call store_error(errmsg)
263  call this%parser%StoreErrorUnit()
264  call ustop()
265  end if
266  call this%parser%GetString(fname)
267  write (this%iout, fmtts) trim(fname)
268  call this%tsmanager%add_tsfile(fname, this%inunit)
269  case default
270  !
271  ! -- Defer to subtype to read the option;
272  ! -- if the subtype can't handle it, report an error
273  if (.not. this%read_option(keyword)) then
274  write (errmsg, '(a,3(1x,a),a)') &
275  'Unknown', trim(adjustl(this%packName)), "option '", &
276  trim(keyword), "'."
277  call store_error(errmsg)
278  end if
279  end select
280  end do
281  write (this%iout, '(1x,a)') &
282  'END OF '//trim(adjustl(this%packName))//' OPTIONS'
283  end if
284  end subroutine read_options
285 
286  !> @brief Read and prepare method for package
287  !!
288  !! Read and prepare stress period data for the package.
289  !!
290  !<
291  subroutine rp(this)
292  ! -- dummy variables
293  class(tvbasetype) :: this
294  ! -- local variables
295  character(len=LINELENGTH) :: line, cellid, varName, text
296  logical :: isfound, endOfBlock, haveChanges
297  integer(I4B) :: ierr, node
298  real(DP), pointer :: bndElem => null()
299  ! -- formats
300  character(len=*), parameter :: fmtblkerr = &
301  &"('Looking for BEGIN PERIOD iper. Found ', a, ' instead.')"
302  character(len=*), parameter :: fmtvalchg = &
303  "(a, ' package: Setting ', a, ' value for cell ', a, ' at start of &
304  &stress period ', i0, ' = ', g12.5)"
305  !
306  if (this%inunit == 0) return
307  !
308  ! -- Get stress period data
309  if (this%ionper < kper) then
310  !
311  ! -- Get PERIOD block
312  call this%parser%GetBlock('PERIOD', isfound, ierr, &
313  supportopenclose=.true., &
314  blockrequired=.false.)
315  if (isfound) then
316  !
317  ! -- Read ionper and check for increasing period numbers
318  call this%read_check_ionper()
319  else
320  !
321  ! -- PERIOD block not found
322  if (ierr < 0) then
323  ! -- End of file found; data applies for remainder of simulation.
324  this%ionper = nper + 1
325  else
326  ! -- Found invalid block
327  call this%parser%GetCurrentLine(line)
328  write (errmsg, fmtblkerr) adjustl(trim(line))
329  call store_error(errmsg)
330  end if
331  end if
332  end if
333  !
334  ! -- Read data if ionper == kper
335  if (this%ionper == kper) then
336  !
337  ! -- Reset per-node property change flags
338  call this%reset_change_flags()
339  !
340  havechanges = .false.
341  do
342  call this%parser%GetNextLine(endofblock)
343  if (endofblock) then
344  exit
345  end if
346  !
347  ! -- Read cell ID
348  call this%parser%GetCellid(this%dis%ndim, cellid)
349  node = this%dis%noder_from_cellid(cellid, this%parser%iuactive, &
350  this%iout)
351  !
352  ! -- Validate cell ID
353  if (node < 1 .or. node > this%dis%nodes) then
354  write (errmsg, '(a,2(1x,a))') &
355  'CELLID', cellid, 'is not in the active model domain.'
356  call store_error(errmsg)
357  cycle
358  end if
359  !
360  ! -- Read variable name
361  call this%parser%GetStringCaps(varname)
362  !
363  ! -- Get a pointer to the property value given by varName for the node
364  ! -- with the specified cell ID
365  bndelem => this%get_pointer_to_value(node, varname)
366  if (.not. associated(bndelem)) then
367  write (errmsg, '(a,3(1x,a),a)') &
368  'Unknown', trim(adjustl(this%packName)), "variable '", &
369  trim(varname), "'."
370  call store_error(errmsg)
371  cycle
372  end if
373  !
374  ! -- Read and apply value or time series link
375  call this%parser%GetString(text)
376  call read_value_or_time_series_adv(text, node, 0, bndelem, &
377  this%packName, 'BND', &
378  this%tsmanager, this%iprpak, &
379  varname)
380  !
381  ! -- Report value change
382  if (this%iprpak /= 0) then
383  write (this%iout, fmtvalchg) &
384  trim(adjustl(this%packName)), trim(varname), trim(cellid), &
385  kper, bndelem
386  end if
387  !
388  ! -- Validate the new property value
389  call this%validate_change(node, varname)
390  havechanges = .true.
391  end do
392  !
393  ! -- Record that any changes were made at the first time step of the
394  ! -- stress period
395  if (havechanges) then
396  call this%set_changed_at(kper, 1)
397  end if
398  end if
399  !
400  ! -- Terminate if errors were encountered in the PERIOD block
401  if (count_errors() > 0) then
402  call this%parser%StoreErrorUnit()
403  call ustop()
404  end if
405  end subroutine rp
406 
407  !> @brief Advance the package
408  !!
409  !! Advance data for a new time step.
410  !!
411  !<
412  subroutine ad(this)
413  ! -- dummy variables
414  class(tvbasetype) :: this
415  ! -- local variables
416  integer(I4B) :: i, n, numlinks
417  type(timeserieslinktype), pointer :: tsLink
418  !
419  ! -- Advance the time series manager;
420  ! -- this will apply any time series changes to property values
421  call this%tsmanager%ad()
422  !
423  ! -- If there are no time series property changes,
424  ! -- there is nothing else to be done
425  numlinks = this%tsmanager%CountLinks('BND')
426  if (numlinks <= 0) then
427  return
428  end if
429  !
430  ! -- Record that changes were made at the current time step
431  ! -- (as we have at least one active time series link)
432  call this%set_changed_at(kper, kstp)
433  !
434  ! -- Reset node K change flags at all time steps except the first of each
435  ! -- period (the first is done in rp(), to allow non-time series changes)
436  if (kstp /= 1) then
437  call this%reset_change_flags()
438  end if
439  !
440  ! -- Validate all new property values
441  do i = 1, numlinks
442  tslink => this%tsmanager%GetLink('BND', i)
443  n = tslink%iRow
444  call this%validate_change(n, tslink%Text)
445  end do
446  !
447  ! -- Terminate if there were errors
448  if (count_errors() > 0) then
449  call this%parser%StoreErrorUnit()
450  call ustop()
451  end if
452  end subroutine ad
453 
454  !> @brief Deallocate package memory
455  !!
456  !! Deallocate package scalars and arrays.
457  !!
458  !<
459  subroutine tvbase_da(this)
460  ! -- dummy variables
461  class(tvbasetype) :: this
462  !
463  ! -- Deallocate time series manager
464  deallocate (this%tsmanager)
465  !
466  ! -- Deallocate parent
467  call this%NumericalPackageType%da()
468  end subroutine tvbase_da
469 
470 end module tvbasemodule
subroutine init()
Definition: GridSorting.f90:24
Announce package and set pointers to variables.
Definition: TvBase.f90:54
Get an array value pointer given a variable name and node index.
Definition: TvBase.f90:83
Announce package and set pointers to variables.
Definition: TvBase.f90:67
Clear all per-node change flags.
Definition: TvBase.f90:118
Mark property changes as having occurred at (kper, kstp)
Definition: TvBase.f90:101
Check that a given property value is valid.
Definition: TvBase.f90:133
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
This module defines variable data types.
Definition: kind.f90:8
This module contains the base numerical package type.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public ustop(stopmess, ioutlocal)
Stop the simulation.
Definition: Sim.f90:312
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
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
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21
subroutine, public read_value_or_time_series_adv(textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, varName)
Call this subroutine from advanced packages to define timeseries link for a variable (varName).
subroutine, public tsmanager_cr(this, iout, removeTsLinksOnCompletion, extendTsToEndOfSimulation)
Create the tsmanager.
This module contains common time-varying property functionality.
Definition: TvBase.f90:8
subroutine read_options(this)
Read OPTIONS block of package input file.
Definition: TvBase.f90:225
subroutine, public tvbase_da(this)
Deallocate package memory.
Definition: TvBase.f90:460
subroutine tvbase_allocate_scalars(this)
Allocate scalar variables.
Definition: TvBase.f90:174
subroutine rp(this)
Read and prepare method for package.
Definition: TvBase.f90:292
subroutine ad(this)
Advance the package.
Definition: TvBase.f90:413
subroutine ar(this, dis)
Allocate and read method for package.
Definition: TvBase.f90:190