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