MODFLOW 6  version 6.7.0.dev0
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)
 Initialise 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 412 of file TvBase.f90.

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
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 189 of file TvBase.f90.

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
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 152 of file TvBase.f90.

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)

◆ 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 224 of file TvBase.f90.

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
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 291 of file TvBase.f90.

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
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 173 of file TvBase.f90.

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)

◆ tvbase_da()

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

Deallocate package scalars and arrays.

Definition at line 459 of file TvBase.f90.

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