MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 71 of file gwf-tvk.f90.

72  ! -- dummy
73  class(TvkType) :: this
74  ! -- local
75  character(len=LENMEMPATH) :: npfMemoryPath
76  ! -- formats
77  character(len=*), parameter :: fmttvk = &
78  "(1x,/1x,'TVK -- TIME-VARYING K PACKAGE, VERSION 1, 08/18/2021', &
79  &' INPUT READ FROM UNIT ', i0, //)"
80  !
81  ! -- Print a message identifying the TVK package
82  write (this%iout, fmttvk) this%inunit
83  !
84  ! -- Set pointers to other package variables
85  ! -- NPF
86  npfmemorypath = create_mem_path(this%name_model, 'NPF')
87  call mem_setptr(this%ik22overk, 'IK22OVERK', npfmemorypath)
88  call mem_setptr(this%ik33overk, 'IK33OVERK', npfmemorypath)
89  call mem_setptr(this%k11, 'K11', npfmemorypath)
90  call mem_setptr(this%k22, 'K22', npfmemorypath)
91  call mem_setptr(this%k33, 'K33', npfmemorypath)
92  call mem_setptr(this%kchangeper, 'KCHANGEPER', npfmemorypath)
93  call mem_setptr(this%kchangestp, 'KCHANGESTP', npfmemorypath)
94  call mem_setptr(this%nodekchange, 'NODEKCHANGE', npfmemorypath)
95  !
96  ! -- Return
97  return
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)
61  !
62  ! -- Return
63  return
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 246 of file gwf-tvk.f90.

247  ! -- dummy
248  class(TvkType) :: this
249  !
250  ! -- Nullify pointers to other package variables
251  nullify (this%ik22overk)
252  nullify (this%ik33overk)
253  nullify (this%k11)
254  nullify (this%k22)
255  nullify (this%k33)
256  nullify (this%kchangeper)
257  nullify (this%kchangestp)
258  nullify (this%nodekchange)
259  !
260  ! -- Deallocate parent
261  call tvbase_da(this)
262  !
263  ! -- Return
264  return
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 124 of file gwf-tvk.f90.

125  ! -- dummy
126  class(TvkType) :: this
127  integer(I4B), intent(in) :: n
128  character(len=*), intent(in) :: varName
129  ! -- return
130  real(DP), pointer :: bndElem
131  !
132  select case (varname)
133  case ('K')
134  bndelem => this%k11(n)
135  case ('K22')
136  bndelem => this%k22(n)
137  case ('K33')
138  bndelem => this%k33(n)
139  case default
140  bndelem => null()
141  end select
142  !
143  ! -- Return
144  return

◆ 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 105 of file gwf-tvk.f90.

106  ! -- dummy
107  class(TvkType) :: this
108  character(len=*), intent(in) :: keyword
109  ! -- return
110  logical :: success
111  !
112  ! -- There are no TVK-specific options, so just return false
113  success = .false.
114  !
115  ! -- Return
116  return

◆ 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 171 of file gwf-tvk.f90.

172  ! -- dummy variables
173  class(TvkType) :: this
174  ! -- local variables
175  integer(I4B) :: i
176  !
177  ! -- Clear NPF's nodekchange array
178  do i = 1, this%dis%nodes
179  this%nodekchange(i) = 0
180  end do
181  !
182  ! -- Return
183  return

◆ 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 152 of file gwf-tvk.f90.

153  ! -- dummy
154  class(TvkType) :: this
155  integer(I4B), intent(in) :: kper
156  integer(I4B), intent(in) :: kstp
157  !
158  this%kchangeper = kper
159  this%kchangestp = kstp
160  !
161  ! -- Return
162  return

◆ 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 193 of file gwf-tvk.f90.

194  ! -- dummy
195  class(TvkType) :: this
196  integer(I4B), intent(in) :: n
197  character(len=*), intent(in) :: varName
198  ! -- local
199  character(len=LINELENGTH) :: cellstr
200  ! -- formats
201  character(len=*), parameter :: fmtkerr = &
202  "(1x, a, ' changed hydraulic property ', a, ' is <= 0 for cell ', a, &
203  &' ', 1pg15.6)"
204  !
205  ! -- Mark the node as being changed this time step
206  this%nodekchange(n) = 1
207  !
208  ! -- Check the changed value is ok
209  if (varname == 'K') then
210  if (this%k11(n) <= dzero) then
211  call this%dis%noder_to_string(n, cellstr)
212  write (errmsg, fmtkerr) &
213  trim(adjustl(this%packName)), 'K', trim(cellstr), this%k11(n)
214  call store_error(errmsg)
215  end if
216  elseif (varname == 'K22') then
217  if (this%ik22overk == 1) then
218  this%k22(n) = this%k22(n) * this%k11(n)
219  end if
220  if (this%k22(n) <= dzero) then
221  call this%dis%noder_to_string(n, cellstr)
222  write (errmsg, fmtkerr) &
223  trim(adjustl(this%packName)), 'K22', trim(cellstr), this%k22(n)
224  call store_error(errmsg)
225  end if
226  elseif (varname == 'K33') then
227  if (this%ik33overk == 1) then
228  this%k33(n) = this%k33(n) * this%k11(n)
229  end if
230  if (this%k33(n) <= dzero) then
231  call this%dis%noder_to_string(n, cellstr)
232  write (errmsg, fmtkerr) &
233  trim(adjustl(this%packName)), 'K33', trim(cellstr), this%k33(n)
234  call store_error(errmsg)
235  end if
236  end if
237  !
238  ! -- Return
239  return
Here is the call graph for this function: