MODFLOW 6  version 6.8.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  apply_row_changes
 Apply input column changes for period-data row n to node. 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, mempath, inunit, iout)
 Initialize the TvBaseType object. More...
 
subroutine tvbase_allocate_scalars (this)
 Allocate scalar variables. More...
 
subroutine tvbase_source_options (this)
 Source common options from the input memory path. More...
 
subroutine tvbase_source_package_options (this)
 Source package-specific options from the input memory path. More...
 
subroutine ar (this, dis)
 Allocate and read static data for the package. More...
 
subroutine rp (this)
 Read and prepare stress period data for the package. More...
 
subroutine ad (this)
 Apply advanced values at each time step. More...
 
subroutine, public tvbase_da (this)
 Deallocate package memory. More...
 
integer(i4b) function tv_get_node (this, n)
 Return the reduced node number for CELLID row n; caller must bounds-check. 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

Definition at line 268 of file TvBase.f90.

269  ! -- dummy
270  class(TvBaseType) :: this
271  ! -- local variables
272  integer(I4B), pointer :: nbound
273  integer(I4B) :: n, node
274  !
275  ! -- no-op when timeseries aren't active.
276  if (.not. this%ts_active) return
277  !
278  call mem_setptr(nbound, 'NBOUND', this%input_mempath)
279  if (nbound > 0) then
280  ! -- Record that changes were made at the current time step
281  call this%set_changed_at(kper, kstp)
282  ! -- Reset node change flags
283  call this%reset_change_flags()
284  ! -- Apply row changes
285  do n = 1, nbound
286  node = this%tv_get_node(n)
287  if (node < 1 .or. node > this%dis%nodes) cycle
288  call this%apply_row_changes(n, node)
289  end do
290  end if
291  !
292  if (count_errors() > 0) then
293  call store_error_filename(this%input_fname)
294  end if
Here is the call graph for this function:

◆ ar()

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

Definition at line 203 of file TvBase.f90.

204  ! -- dummy
205  class(TvBaseType) :: this
206  class(DisBaseType), pointer, intent(in) :: dis
207  !
208  this%dis => dis
209  call this%ar_set_pointers()
210  !
211  call this%source_options()
212  !
213  ! -- set input mempath pointers
214  call mem_setptr(this%cellid, 'CELLID', this%input_mempath)
215  !
216  if (count_errors() > 0) then
217  call store_error_filename(this%input_fname)
218  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,
character(len=*), intent(in)  mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Allocate and initialize data members of the object.

Definition at line 127 of file TvBase.f90.

128  ! -- dummy
129  class(TvBaseType) :: this
130  character(len=*), intent(in) :: name_model
131  character(len=*), intent(in) :: pakname
132  character(len=*), intent(in) :: ftype
133  character(len=*), intent(in) :: mempath
134  integer(I4B), intent(in) :: inunit
135  integer(I4B), intent(in) :: iout
136  !
137  call this%set_names(1, name_model, pakname, ftype, mempath)
138  call this%tvbase_allocate_scalars()
139  this%inunit = inunit
140  this%iout = iout

◆ rp()

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

Definition at line 223 of file TvBase.f90.

224  ! -- dummy
225  class(TvBaseType) :: this
226  ! -- local variables
227  integer(I4B), pointer :: iper, nbound
228  integer(I4B) :: n, node
229  character(len=LINELENGTH) :: cellstr
230  !
231  ! -- check last loaded input period
232  call mem_setptr(iper, 'IPER', this%input_mempath)
233  if (iper /= kper) return
234  !
235  ! -- When timeseries are active, ad applies values at every time step
236  if (this%ts_active) return
237  !
238  call mem_setptr(nbound, 'NBOUND', this%input_mempath)
239  !
240  ! -- Reset per-node property change flags
241  call this%reset_change_flags()
242  !
243  do n = 1, nbound
244  node = this%tv_get_node(n)
245  if (node < 1 .or. node > this%dis%nodes) then
246  call this%dis%noder_to_string(node, cellstr)
247  write (errmsg, '(a,2(1x,a))') &
248  'CELLID', trim(cellstr), 'is not in the active model domain.'
249  call store_error(errmsg)
250  cycle
251  end if
252  !
253  call this%apply_row_changes(n, node)
254  end do
255  !
256  ! -- Record that changes were made at the current stress period / time step
257  if (nbound > 0) then
258  call this%set_changed_at(kper, kstp)
259  end if
260  !
261  if (count_errors() > 0) then
262  call store_error_filename(this%input_fname)
263  end if
Here is the call graph for this function:

◆ tv_get_node()

integer(i4b) function tvbasemodule::tv_get_node ( class(tvbasetype this,
integer(i4b), intent(in)  n 
)
private
Parameters
[in]nrow index in the period-data arrays

Definition at line 311 of file TvBase.f90.

312  ! -- dummy
313  class(TvBaseType) :: this
314  integer(I4B), intent(in) :: n !< row index in the period-data arrays
315  ! -- return
316  integer(I4B) :: node
317  ! -- local
318  integer(I4B) :: nodeu
319  !
320  if (this%dis%ndim == 1) then
321  nodeu = this%cellid(1, n)
322  elseif (this%dis%ndim == 2) then
323  nodeu = get_node(this%cellid(1, n), 1, &
324  this%cellid(2, n), &
325  this%dis%mshape(1), 1, &
326  this%dis%mshape(2))
327  else
328  nodeu = get_node(this%cellid(1, n), &
329  this%cellid(2, n), &
330  this%cellid(3, n), &
331  this%dis%mshape(1), &
332  this%dis%mshape(2), &
333  this%dis%mshape(3))
334  end if
335  node = this%dis%get_nodenumber(nodeu, 1)
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 147 of file TvBase.f90.

148  ! -- dummy
149  class(TvBaseType) :: this
150  !
151  ! -- Call standard NumericalPackageType allocate scalars
152  call this%NumericalPackageType%allocate_scalars()

◆ tvbase_da()

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

Deallocate package scalars and arrays.

Definition at line 301 of file TvBase.f90.

302  ! -- dummy
303  class(TvBaseType) :: this
304  !
305  nullify (this%cellid)
306  call this%NumericalPackageType%da()
Here is the caller graph for this function:

◆ tvbase_source_options()

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

Source common options and call derived package routine to source and log any package-specific options within the same block.

Definition at line 160 of file TvBase.f90.

161  ! -- modules
163  ! -- dummy
164  class(TvBaseType) :: this
165  ! -- locals
166  integer(I4B) :: isize
167  logical(LGP) :: found_print_input
168  !
169  write (this%iout, '(1x,a)') &
170  'PROCESSING '//trim(adjustl(this%packName))//' OPTIONS'
171  !
172  call mem_set_value(this%iprpak, 'PRINT_INPUT', this%input_mempath, &
173  found_print_input)
174  !
175  if (found_print_input) then
176  write (this%iout, '(4x,a)') 'TIME-VARYING INPUT WILL BE PRINTED.'
177  end if
178  !
179  call get_isize('TS6_FILENAME', this%input_mempath, isize)
180  if (isize > 0) this%ts_active = .true.
181  !
182  ! -- source package-specific options
183  call this%source_package_options()
184  !
185  write (this%iout, '(1x,a)') &
186  'END OF '//trim(adjustl(this%packName))//' OPTIONS'
Here is the call graph for this function:

◆ tvbase_source_package_options()

subroutine tvbasemodule::tvbase_source_package_options ( class(tvbasetype this)

Override in derived package to source and log package-specific options. Default implementation is a no-op.

Definition at line 194 of file TvBase.f90.

195  ! -- dummy
196  class(TvBaseType) :: this
197  !
198  ! -- no package-specific options in the base class