MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 422 of file TvBase.f90.

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

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
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)
166  !
167  return

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

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

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

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

◆ tvbase_da()

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

Deallocate package scalars and arrays.

Definition at line 471 of file TvBase.f90.

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