MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
tvkmodule Module Reference

This module contains time-varying conductivity package methods. More...

Data Types

type  tvktype
 

Functions/Subroutines

subroutine, public tvk_cr (tvk, name_model, inunit, iout)
 Create a new TvkType object. More...
 
subroutine tvk_ar_set_pointers (this)
 Announce package and set pointers to variables. More...
 
logical function tvk_read_option (this, keyword)
 Read a TVK-specific option from the OPTIONS block. More...
 
real(dp) function, pointer tvk_get_pointer_to_value (this, n, varName)
 Get an array value pointer given a variable name and node index. More...
 
subroutine tvk_set_changed_at (this, kper, kstp)
 Mark property changes as having occurred at (kper, kstp) More...
 
subroutine tvk_reset_change_flags (this)
 Clear all per-node change flags. More...
 
subroutine tvk_validate_change (this, n, varName)
 Check that a given property value is valid. More...
 
subroutine tvk_da (this)
 Deallocate package memory. More...
 

Detailed Description

This module contains the methods used to allow hydraulic conductivity parameters in the NPF package (K11, K22, K33) to be varied throughout a simulation.

Function/Subroutine Documentation

◆ tvk_ar_set_pointers()

subroutine tvkmodule::tvk_ar_set_pointers ( class(tvktype this)
private

Announce package version and set array and variable pointers from the NPF package for access by TVK.

Definition at line 68 of file gwf-tvk.f90.

69  ! -- dummy
70  class(TvkType) :: this
71  ! -- local
72  character(len=LENMEMPATH) :: npfMemoryPath
73  ! -- formats
74  character(len=*), parameter :: fmttvk = &
75  "(1x,/1x,'TVK -- TIME-VARYING K PACKAGE, VERSION 1, 08/18/2021', &
76  &' INPUT READ FROM UNIT ', i0, //)"
77  !
78  ! -- Print a message identifying the TVK package
79  write (this%iout, fmttvk) this%inunit
80  !
81  ! -- Set pointers to other package variables
82  ! -- NPF
83  npfmemorypath = create_mem_path(this%name_model, 'NPF')
84  call mem_setptr(this%ik22overk, 'IK22OVERK', npfmemorypath)
85  call mem_setptr(this%ik33overk, 'IK33OVERK', npfmemorypath)
86  call mem_setptr(this%k11, 'K11', npfmemorypath)
87  call mem_setptr(this%k22, 'K22', npfmemorypath)
88  call mem_setptr(this%k33, 'K33', npfmemorypath)
89  call mem_setptr(this%kchangeper, 'KCHANGEPER', npfmemorypath)
90  call mem_setptr(this%kchangestp, 'KCHANGESTP', npfmemorypath)
91  call mem_setptr(this%nodekchange, 'NODEKCHANGE', npfmemorypath)
Here is the call graph for this function:

◆ tvk_cr()

subroutine, public tvkmodule::tvk_cr ( type(tvktype), intent(out), pointer  tvk,
character(len=*), intent(in)  name_model,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Create a new time-varying conductivity (TvkType) object.

Definition at line 52 of file gwf-tvk.f90.

53  ! -- dummy
54  type(TvkType), pointer, intent(out) :: tvk
55  character(len=*), intent(in) :: name_model
56  integer(I4B), intent(in) :: inunit
57  integer(I4B), intent(in) :: iout
58  !
59  allocate (tvk)
60  call tvk%init(name_model, 'TVK', 'TVK', inunit, iout)
Here is the caller graph for this function:

◆ tvk_da()

subroutine tvkmodule::tvk_da ( class(tvktype this)
private

Deallocate TVK package scalars and arrays.

Definition at line 225 of file gwf-tvk.f90.

226  ! -- dummy
227  class(TvkType) :: this
228  !
229  ! -- Nullify pointers to other package variables
230  nullify (this%ik22overk)
231  nullify (this%ik33overk)
232  nullify (this%k11)
233  nullify (this%k22)
234  nullify (this%k33)
235  nullify (this%kchangeper)
236  nullify (this%kchangestp)
237  nullify (this%nodekchange)
238  !
239  ! -- Deallocate parent
240  call tvbase_da(this)
Here is the call graph for this function:

◆ tvk_get_pointer_to_value()

real(dp) function, pointer tvkmodule::tvk_get_pointer_to_value ( class(tvktype this,
integer(i4b), intent(in)  n,
character(len=*), intent(in)  varName 
)
private

Return a pointer to the given node's value in the appropriate NPF array based on the given variable name string.

Definition at line 115 of file gwf-tvk.f90.

116  ! -- dummy
117  class(TvkType) :: this
118  integer(I4B), intent(in) :: n
119  character(len=*), intent(in) :: varName
120  ! -- return
121  real(DP), pointer :: bndElem
122  !
123  select case (varname)
124  case ('K')
125  bndelem => this%k11(n)
126  case ('K22')
127  bndelem => this%k22(n)
128  case ('K33')
129  bndelem => this%k33(n)
130  case default
131  bndelem => null()
132  end select

◆ tvk_read_option()

logical function tvkmodule::tvk_read_option ( class(tvktype this,
character(len=*), intent(in)  keyword 
)
private

Process a single TVK-specific option. Used when reading the OPTIONS block of the TVK package input file.

Definition at line 99 of file gwf-tvk.f90.

100  ! -- dummy
101  class(TvkType) :: this
102  character(len=*), intent(in) :: keyword
103  ! -- return
104  logical :: success
105  !
106  ! -- There are no TVK-specific options, so just return false
107  success = .false.

◆ tvk_reset_change_flags()

subroutine tvkmodule::tvk_reset_change_flags ( class(tvktype this)
private

Deferred procedure implementation called by the TvBaseType code when a new time step commences, indicating that any previously set per-node property value change flags should be reset.

Definition at line 156 of file gwf-tvk.f90.

157  ! -- dummy variables
158  class(TvkType) :: this
159  ! -- local variables
160  integer(I4B) :: i
161  !
162  ! -- Clear NPF's nodekchange array
163  do i = 1, this%dis%nodes
164  this%nodekchange(i) = 0
165  end do

◆ tvk_set_changed_at()

subroutine tvkmodule::tvk_set_changed_at ( class(tvktype this,
integer(i4b), intent(in)  kper,
integer(i4b), intent(in)  kstp 
)
private

Deferred procedure implementation called by the TvBaseType code when a property value change occurs at (kper, kstp).

Definition at line 140 of file gwf-tvk.f90.

141  ! -- dummy
142  class(TvkType) :: this
143  integer(I4B), intent(in) :: kper
144  integer(I4B), intent(in) :: kstp
145  !
146  this%kchangeper = kper
147  this%kchangestp = kstp

◆ tvk_validate_change()

subroutine tvkmodule::tvk_validate_change ( class(tvktype this,
integer(i4b), intent(in)  n,
character(len=*), intent(in)  varName 
)
private

Deferred procedure implementation called by the TvBaseType code after a property value change occurs. Check if the property value of the given variable at the given node is invalid, and log an error if so. Update K22 and K33 values appropriately when specified as anisotropy.

Definition at line 175 of file gwf-tvk.f90.

176  ! -- dummy
177  class(TvkType) :: this
178  integer(I4B), intent(in) :: n
179  character(len=*), intent(in) :: varName
180  ! -- local
181  character(len=LINELENGTH) :: cellstr
182  ! -- formats
183  character(len=*), parameter :: fmtkerr = &
184  "(1x, a, ' changed hydraulic property ', a, ' is <= 0 for cell ', a, &
185  &' ', 1pg15.6)"
186  !
187  ! -- Mark the node as being changed this time step
188  this%nodekchange(n) = 1
189  !
190  ! -- Check the changed value is ok
191  if (varname == 'K') then
192  if (this%k11(n) <= dzero) then
193  call this%dis%noder_to_string(n, cellstr)
194  write (errmsg, fmtkerr) &
195  trim(adjustl(this%packName)), 'K', trim(cellstr), this%k11(n)
196  call store_error(errmsg)
197  end if
198  elseif (varname == 'K22') then
199  if (this%ik22overk == 1) then
200  this%k22(n) = this%k22(n) * this%k11(n)
201  end if
202  if (this%k22(n) <= dzero) then
203  call this%dis%noder_to_string(n, cellstr)
204  write (errmsg, fmtkerr) &
205  trim(adjustl(this%packName)), 'K22', trim(cellstr), this%k22(n)
206  call store_error(errmsg)
207  end if
208  elseif (varname == 'K33') then
209  if (this%ik33overk == 1) then
210  this%k33(n) = this%k33(n) * this%k11(n)
211  end if
212  if (this%k33(n) <= dzero) then
213  call this%dis%noder_to_string(n, cellstr)
214  write (errmsg, fmtkerr) &
215  trim(adjustl(this%packName)), 'K33', trim(cellstr), this%k33(n)
216  call store_error(errmsg)
217  end if
218  end if
Here is the call graph for this function: