MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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  !
167  return
168  end subroutine init
169 
170  !> @brief Allocate scalar variables
171  !!
172  !! Allocate scalar data members of the object.
173  !!
174  !<
175  subroutine tvbase_allocate_scalars(this)
176  ! -- dummy variables
177  class(tvbasetype) :: this
178  !
179  ! -- Call standard NumericalPackageType allocate scalars
180  call this%NumericalPackageType%allocate_scalars()
181  !
182  ! -- Allocate time series manager
183  allocate (this%tsmanager)
184  !
185  return
186  end subroutine tvbase_allocate_scalars
187 
188  !> @brief Allocate and read method for package
189  !!
190  !! Allocate and read static data for the package.
191  !!
192  !<
193  subroutine ar(this, dis)
194  ! -- dummy variables
195  class(tvbasetype) :: this
196  class(disbasetype), pointer, intent(in) :: dis
197  !
198  ! -- Set pointers to other package variables and announce package
199  this%dis => dis
200  call this%ar_set_pointers()
201  !
202  ! -- Create time series manager
203  call tsmanager_cr(this%tsmanager, &
204  this%iout, &
205  removetslinksoncompletion=.true., &
206  extendtstoendofsimulation=.true.)
207  !
208  ! -- Read options
209  call this%read_options()
210  !
211  ! -- Now that time series will have been read, need to call the df routine
212  ! -- to define the manager
213  call this%tsmanager%tsmanager_df()
214  !
215  ! -- Terminate if any errors were encountered
216  if (count_errors() > 0) then
217  call this%parser%StoreErrorUnit()
218  call ustop()
219  end if
220  !
221  return
222  end subroutine ar
223 
224  !> @brief Read OPTIONS block of package input file
225  !!
226  !! Reads the OPTIONS block of the package's input file, deferring to the
227  !! derived type to process any package-specific keywords.
228  !!
229  !<
230  subroutine read_options(this)
231  ! -- dummy variables
232  class(tvbasetype) :: this
233  ! -- local variables
234  character(len=LINELENGTH) :: keyword
235  character(len=MAXCHARLEN) :: fname
236  logical :: isfound
237  logical :: endOfBlock
238  integer(I4B) :: ierr
239  ! -- formats
240  character(len=*), parameter :: fmtts = &
241  &"(4x, 'TIME-SERIES DATA WILL BE READ FROM FILE: ', a)"
242  !
243  ! -- Get options block
244  call this%parser%GetBlock('OPTIONS', isfound, ierr, &
245  blockrequired=.false., supportopenclose=.true.)
246  !
247  ! -- Parse options block if detected
248  if (isfound) then
249  write (this%iout, '(1x,a)') &
250  'PROCESSING '//trim(adjustl(this%packName))//' OPTIONS'
251  do
252  call this%parser%GetNextLine(endofblock)
253  if (endofblock) then
254  exit
255  end if
256  call this%parser%GetStringCaps(keyword)
257  select case (keyword)
258  case ('PRINT_INPUT')
259  this%iprpak = 1
260  write (this%iout, '(4x,a)') 'TIME-VARYING INPUT WILL BE PRINTED.'
261  case ('TS6')
262  !
263  ! -- Add a time series file
264  call this%parser%GetStringCaps(keyword)
265  if (trim(adjustl(keyword)) /= 'FILEIN') then
266  errmsg = &
267  'TS6 keyword must be followed by "FILEIN" then by filename.'
268  call store_error(errmsg)
269  call this%parser%StoreErrorUnit()
270  call ustop()
271  end if
272  call this%parser%GetString(fname)
273  write (this%iout, fmtts) trim(fname)
274  call this%tsmanager%add_tsfile(fname, this%inunit)
275  case default
276  !
277  ! -- Defer to subtype to read the option;
278  ! -- if the subtype can't handle it, report an error
279  if (.not. this%read_option(keyword)) then
280  write (errmsg, '(a,3(1x,a),a)') &
281  'Unknown', trim(adjustl(this%packName)), "option '", &
282  trim(keyword), "'."
283  call store_error(errmsg)
284  end if
285  end select
286  end do
287  write (this%iout, '(1x,a)') &
288  'END OF '//trim(adjustl(this%packName))//' OPTIONS'
289  end if
290  !
291  return
292  end subroutine read_options
293 
294  !> @brief Read and prepare method for package
295  !!
296  !! Read and prepare stress period data for the package.
297  !!
298  !<
299  subroutine rp(this)
300  ! -- dummy variables
301  class(tvbasetype) :: this
302  ! -- local variables
303  character(len=LINELENGTH) :: line, cellid, varName, text
304  logical :: isfound, endOfBlock, haveChanges
305  integer(I4B) :: ierr, node
306  real(DP), pointer :: bndElem => null()
307  ! -- formats
308  character(len=*), parameter :: fmtblkerr = &
309  &"('Looking for BEGIN PERIOD iper. Found ', a, ' instead.')"
310  character(len=*), parameter :: fmtvalchg = &
311  "(a, ' package: Setting ', a, ' value for cell ', a, ' at start of &
312  &stress period ', i0, ' = ', g12.5)"
313  !
314  if (this%inunit == 0) return
315  !
316  ! -- Get stress period data
317  if (this%ionper < kper) then
318  !
319  ! -- Get PERIOD block
320  call this%parser%GetBlock('PERIOD', isfound, ierr, &
321  supportopenclose=.true., &
322  blockrequired=.false.)
323  if (isfound) then
324  !
325  ! -- Read ionper and check for increasing period numbers
326  call this%read_check_ionper()
327  else
328  !
329  ! -- PERIOD block not found
330  if (ierr < 0) then
331  ! -- End of file found; data applies for remainder of simulation.
332  this%ionper = nper + 1
333  else
334  ! -- Found invalid block
335  call this%parser%GetCurrentLine(line)
336  write (errmsg, fmtblkerr) adjustl(trim(line))
337  call store_error(errmsg)
338  end if
339  end if
340  end if
341  !
342  ! -- Read data if ionper == kper
343  if (this%ionper == kper) then
344  !
345  ! -- Reset per-node property change flags
346  call this%reset_change_flags()
347  !
348  havechanges = .false.
349  do
350  call this%parser%GetNextLine(endofblock)
351  if (endofblock) then
352  exit
353  end if
354  !
355  ! -- Read cell ID
356  call this%parser%GetCellid(this%dis%ndim, cellid)
357  node = this%dis%noder_from_cellid(cellid, this%parser%iuactive, &
358  this%iout)
359  !
360  ! -- Validate cell ID
361  if (node < 1 .or. node > this%dis%nodes) then
362  write (errmsg, '(a,2(1x,a))') &
363  'CELLID', cellid, 'is not in the active model domain.'
364  call store_error(errmsg)
365  cycle
366  end if
367  !
368  ! -- Read variable name
369  call this%parser%GetStringCaps(varname)
370  !
371  ! -- Get a pointer to the property value given by varName for the node
372  ! -- with the specified cell ID
373  bndelem => this%get_pointer_to_value(node, varname)
374  if (.not. associated(bndelem)) then
375  write (errmsg, '(a,3(1x,a),a)') &
376  'Unknown', trim(adjustl(this%packName)), "variable '", &
377  trim(varname), "'."
378  call store_error(errmsg)
379  cycle
380  end if
381  !
382  ! -- Read and apply value or time series link
383  call this%parser%GetString(text)
384  call read_value_or_time_series_adv(text, node, 0, bndelem, &
385  this%packName, 'BND', &
386  this%tsmanager, this%iprpak, &
387  varname)
388  !
389  ! -- Report value change
390  if (this%iprpak /= 0) then
391  write (this%iout, fmtvalchg) &
392  trim(adjustl(this%packName)), trim(varname), trim(cellid), &
393  kper, bndelem
394  end if
395  !
396  ! -- Validate the new property value
397  call this%validate_change(node, varname)
398  havechanges = .true.
399  end do
400  !
401  ! -- Record that any changes were made at the first time step of the
402  ! -- stress period
403  if (havechanges) then
404  call this%set_changed_at(kper, 1)
405  end if
406  end if
407  !
408  ! -- Terminate if errors were encountered in the PERIOD block
409  if (count_errors() > 0) then
410  call this%parser%StoreErrorUnit()
411  call ustop()
412  end if
413  !
414  return
415  end subroutine rp
416 
417  !> @brief Advance the package
418  !!
419  !! Advance data for a new time step.
420  !!
421  !<
422  subroutine ad(this)
423  ! -- dummy variables
424  class(tvbasetype) :: this
425  ! -- local variables
426  integer(I4B) :: i, n, numlinks
427  type(timeserieslinktype), pointer :: tsLink
428  !
429  ! -- Advance the time series manager;
430  ! -- this will apply any time series changes to property values
431  call this%tsmanager%ad()
432  !
433  ! -- If there are no time series property changes,
434  ! -- there is nothing else to be done
435  numlinks = this%tsmanager%CountLinks('BND')
436  if (numlinks <= 0) then
437  return
438  end if
439  !
440  ! -- Record that changes were made at the current time step
441  ! -- (as we have at least one active time series link)
442  call this%set_changed_at(kper, kstp)
443  !
444  ! -- Reset node K change flags at all time steps except the first of each
445  ! -- period (the first is done in rp(), to allow non-time series changes)
446  if (kstp /= 1) then
447  call this%reset_change_flags()
448  end if
449  !
450  ! -- Validate all new property values
451  do i = 1, numlinks
452  tslink => this%tsmanager%GetLink('BND', i)
453  n = tslink%iRow
454  call this%validate_change(n, tslink%Text)
455  end do
456  !
457  ! -- Terminate if there were errors
458  if (count_errors() > 0) then
459  call this%parser%StoreErrorUnit()
460  call ustop()
461  end if
462  !
463  return
464  end subroutine ad
465 
466  !> @brief Deallocate package memory
467  !!
468  !! Deallocate package scalars and arrays.
469  !!
470  !<
471  subroutine tvbase_da(this)
472  ! -- dummy variables
473  class(tvbasetype) :: this
474  !
475  ! -- Deallocate time series manager
476  deallocate (this%tsmanager)
477  !
478  ! -- Deallocate parent
479  call this%NumericalPackageType%da()
480  !
481  return
482  end subroutine tvbase_da
483 
484 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:44
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:46
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:315
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:231
subroutine, public tvbase_da(this)
Deallocate package memory.
Definition: TvBase.f90:472
subroutine tvbase_allocate_scalars(this)
Allocate scalar variables.
Definition: TvBase.f90:176
subroutine rp(this)
Read and prepare method for package.
Definition: TvBase.f90:300
subroutine ad(this)
Advance the package.
Definition: TvBase.f90:423
subroutine ar(this, dis)
Allocate and read method for package.
Definition: TvBase.f90:194