MODFLOW 6  version 6.8.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, mempath, inunit, iout)
 Create a new TvkType object. More...
 
subroutine tvk_ar_set_pointers (this)
 Announce package and set pointers to variables. More...
 
subroutine tvk_apply_row_changes (this, n, node)
 Apply input K/K22/K33 column changes for period-data row n to node. 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_apply_row_changes()

subroutine tvkmodule::tvk_apply_row_changes ( class(tvktype this,
integer(i4b), intent(in)  n,
integer(i4b), intent(in)  node 
)
private

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

103  ! -- dummy
104  class(TvkType) :: this
105  integer(I4B), intent(in) :: n
106  integer(I4B), intent(in) :: node
107  ! -- local
108  character(len=LINELENGTH) :: cellstr
109  ! -- formats
110  character(len=*), parameter :: fmtvalchg = &
111  "(a, ' package: Setting ', a, ' value for cell ', a, ' at start of &
112  &stress period ', i0, ' = ', g12.5)"
113  !
114  ! -- K is processed before K22/K33 so that validate_change can use
115  ! -- the already-updated k11 value when ik22overk/ik33overk are set.
116  if (this%k11_src(n) /= dnodata) then
117  this%k11(node) = this%k11_src(n)
118  call this%validate_change(node, 'K')
119  if (this%iprpak /= 0) then
120  call this%dis%noder_to_string(node, cellstr)
121  write (this%iout, fmtvalchg) &
122  trim(adjustl(this%packName)), 'K', trim(cellstr), kper, this%k11(node)
123  end if
124  end if
125  !
126  if (this%k22_src(n) /= dnodata) then
127  this%k22(node) = this%k22_src(n)
128  call this%validate_change(node, 'K22')
129  if (this%iprpak /= 0) then
130  call this%dis%noder_to_string(node, cellstr)
131  write (this%iout, fmtvalchg) &
132  trim(adjustl(this%packName)), 'K22', trim(cellstr), kper, &
133  this%k22(node)
134  end if
135  end if
136  !
137  if (this%k33_src(n) /= dnodata) then
138  this%k33(node) = this%k33_src(n)
139  call this%validate_change(node, 'K33')
140  if (this%iprpak /= 0) then
141  call this%dis%noder_to_string(node, cellstr)
142  write (this%iout, fmtvalchg) &
143  trim(adjustl(this%packName)), 'K33', trim(cellstr), kper, &
144  this%k33(node)
145  end if
146  end if

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

73  ! -- dummy
74  class(TvkType) :: this
75  ! -- local
76  character(len=LENMEMPATH) :: npfMemoryPath
77  ! -- formats
78  character(len=*), parameter :: fmttvk = &
79  "(1x,/1x,'TVK -- TIME-VARYING K PACKAGE, VERSION 1, 08/18/2021', &
80  &' INPUT READ FROM MEMPATH ', A, //)"
81  !
82  write (this%iout, fmttvk) this%input_mempath
83  !
84  npfmemorypath = create_mem_path(this%name_model, 'NPF')
85  call mem_setptr(this%ik22overk, 'IK22OVERK', npfmemorypath)
86  call mem_setptr(this%ik33overk, 'IK33OVERK', npfmemorypath)
87  call mem_setptr(this%k11, 'K11', npfmemorypath)
88  call mem_setptr(this%k22, 'K22', npfmemorypath)
89  call mem_setptr(this%k33, 'K33', npfmemorypath)
90  call mem_setptr(this%kchangeper, 'KCHANGEPER', npfmemorypath)
91  call mem_setptr(this%kchangestp, 'KCHANGESTP', npfmemorypath)
92  call mem_setptr(this%nodekchange, 'NODEKCHANGE', npfmemorypath)
93  !
94  ! -- set input mempath pointers
95  call mem_setptr(this%k11_src, 'K', this%input_mempath)
96  call mem_setptr(this%k22_src, 'K22', this%input_mempath)
97  call mem_setptr(this%k33_src, 'K33', this%input_mempath)
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,
character(len=*), intent(in)  mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

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

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

56  ! -- dummy
57  type(TvkType), pointer, intent(out) :: tvk
58  character(len=*), intent(in) :: name_model
59  character(len=*), intent(in) :: mempath
60  integer(I4B), intent(in) :: inunit
61  integer(I4B), intent(in) :: iout
62  !
63  allocate (tvk)
64  call tvk%init(name_model, 'TVK', 'TVK', mempath, 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 238 of file gwf-tvk.f90.

239  ! -- dummy
240  class(TvkType) :: this
241  !
242  nullify (this%ik22overk)
243  nullify (this%ik33overk)
244  nullify (this%k11)
245  nullify (this%k22)
246  nullify (this%k33)
247  nullify (this%kchangeper)
248  nullify (this%kchangestp)
249  nullify (this%nodekchange)
250  nullify (this%k11_src)
251  nullify (this%k22_src)
252  nullify (this%k33_src)
253  call tvbase_da(this)
Here is the call graph for this function:

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

171  ! -- dummy variables
172  class(TvkType) :: this
173  ! -- local variables
174  integer(I4B) :: i
175  !
176  do i = 1, this%dis%nodes
177  this%nodekchange(i) = 0
178  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 154 of file gwf-tvk.f90.

155  ! -- dummy
156  class(TvkType) :: this
157  integer(I4B), intent(in) :: kper
158  integer(I4B), intent(in) :: kstp
159  !
160  this%kchangeper = kper
161  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 188 of file gwf-tvk.f90.

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