MODFLOW 6  version 6.7.0.dev1
USGS Modular Hydrologic Model
tvbasemodule Module Reference

This module contains common time-varying property functionality. More...

Data Types

type  tvbasetype
 
interface  ar_set_pointers
 Announce package and set pointers to variables. More...
 
interface  read_option
 Announce package and set pointers to variables. More...
 
interface  get_pointer_to_value
 Get an array value pointer given a variable name and node index. More...
 
interface  set_changed_at
 Mark property changes as having occurred at (kper, kstp) More...
 
interface  reset_change_flags
 Clear all per-node change flags. More...
 
interface  validate_change
 Check that a given property value is valid. More...
 

Functions/Subroutines

subroutine init (this, name_model, pakname, ftype, inunit, iout)
 Initialize the TvBaseType object. More...
 
subroutine tvbase_allocate_scalars (this)
 Allocate scalar variables. More...
 
subroutine ar (this, dis)
 Allocate and read method for package. More...
 
subroutine read_options (this)
 Read OPTIONS block of package input file. More...
 
subroutine rp (this)
 Read and prepare method for package. More...
 
subroutine ad (this)
 Advance the package. More...
 
subroutine, public tvbase_da (this)
 Deallocate package memory. More...
 

Detailed Description

This module contains methods implementing functionality common to both time-varying hydraulic conductivity (TVK) and time-varying storage (TVS) packages.

Function/Subroutine Documentation

◆ ad()

subroutine tvbasemodule::ad ( class(tvbasetype this)
private

Advance data for a new time step.

Definition at line 400 of file TvBase.f90.

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
Here is the call graph for this function:

◆ ar()

subroutine tvbasemodule::ar ( class(tvbasetype this,
class(disbasetype), intent(in), pointer  dis 
)
private

Allocate and read static data for the package.

Definition at line 180 of file TvBase.f90.

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
Here is the call graph for this function:

◆ init()

subroutine tvbasemodule::init ( class(tvbasetype this,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  ftype,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Allocate and initialize data members of the object.

Definition at line 145 of file TvBase.f90.

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)

◆ read_options()

subroutine tvbasemodule::read_options ( class(tvbasetype this)
private

Reads the OPTIONS block of the package's input file, deferring to the derived type to process any package-specific keywords.

Definition at line 214 of file TvBase.f90.

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
Here is the call graph for this function:

◆ rp()

subroutine tvbasemodule::rp ( class(tvbasetype this)
private

Read and prepare stress period data for the package.

Definition at line 280 of file TvBase.f90.

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
Here is the call graph for this function:

◆ tvbase_allocate_scalars()

subroutine tvbasemodule::tvbase_allocate_scalars ( class(tvbasetype this)
private

Allocate scalar data members of the object.

Definition at line 165 of file TvBase.f90.

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)

◆ tvbase_da()

subroutine, public tvbasemodule::tvbase_da ( class(tvbasetype this)

Deallocate package scalars and arrays.

Definition at line 446 of file TvBase.f90.

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()
Here is the caller graph for this function: