27 integer(I4B),
pointer :: ik22overk => null()
28 integer(I4B),
pointer :: ik33overk => null()
29 real(dp),
dimension(:),
pointer,
contiguous :: k11 => null()
30 real(dp),
dimension(:),
pointer,
contiguous :: k22 => null()
31 real(dp),
dimension(:),
pointer,
contiguous :: k33 => null()
32 integer(I4B),
pointer :: kchangeper => null()
33 integer(I4B),
pointer :: kchangestp => null()
34 integer(I4B),
dimension(:),
pointer,
contiguous :: nodekchange => null()
35 real(dp),
dimension(:),
pointer,
contiguous :: k11_src => null()
36 real(dp),
dimension(:),
pointer,
contiguous :: k22_src => null()
37 real(dp),
dimension(:),
pointer,
contiguous :: k33_src => null()
55 subroutine tvk_cr(tvk, name_model, mempath, inunit, iout)
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
64 call tvk%init(name_model,
'TVK',
'TVK', mempath, inunit, iout)
76 character(len=LENMEMPATH) :: npfMemoryPath
78 character(len=*),
parameter :: fmttvk = &
79 "(1x,/1x,'TVK -- TIME-VARYING K PACKAGE, VERSION 1, 08/18/2021', &
80 &' INPUT READ FROM MEMPATH ', A, //)"
82 write (this%iout, fmttvk) this%input_mempath
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)
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)
105 integer(I4B),
intent(in) :: n
106 integer(I4B),
intent(in) :: node
108 character(len=LINELENGTH) :: cellstr
110 character(len=*),
parameter :: fmtvalchg = &
111 "(a, ' package: Setting ', a, ' value for cell ', a, ' at start of &
112 &stress period ', i0, ' = ', g12.5)"
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)
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, &
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, &
157 integer(I4B),
intent(in) :: kper
158 integer(I4B),
intent(in) :: kstp
160 this%kchangeper = kper
161 this%kchangestp = kstp
176 do i = 1, this%dis%nodes
177 this%nodekchange(i) = 0
191 integer(I4B),
intent(in) :: n
192 character(len=*),
intent(in) :: varName
194 character(len=LINELENGTH) :: cellstr
196 character(len=*),
parameter :: fmtkerr = &
197 "(1x, a, ' changed hydraulic property ', a, ' is <= 0 for cell ', a, &
201 this%nodekchange(n) = 1
204 if (varname ==
'K')
then
205 if (this%k11(n) <=
dzero)
then
206 call this%dis%noder_to_string(n, cellstr)
208 trim(adjustl(this%packName)),
'K', trim(cellstr), this%k11(n)
211 elseif (varname ==
'K22')
then
212 if (this%ik22overk == 1)
then
213 this%k22(n) = this%k22(n) * this%k11(n)
215 if (this%k22(n) <=
dzero)
then
216 call this%dis%noder_to_string(n, cellstr)
218 trim(adjustl(this%packName)),
'K22', trim(cellstr), this%k22(n)
221 elseif (varname ==
'K33')
then
222 if (this%ik33overk == 1)
then
223 this%k33(n) = this%k33(n) * this%k11(n)
225 if (this%k33(n) <=
dzero)
then
226 call this%dis%noder_to_string(n, cellstr)
228 trim(adjustl(this%packName)),
'K33', trim(cellstr), this%k33(n)
242 nullify (this%ik22overk)
243 nullify (this%ik33overk)
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)
Apply input column changes for period-data row n to node.
Announce package and set pointers to variables.
Clear all per-node change flags.
Mark property changes as having occurred at (kper, kstp)
Check that a given property value is valid.
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
real(dp), parameter dnodata
real no data constant
real(dp), parameter dzero
real constant zero
integer(i4b), parameter lenmempath
maximum length of the memory path
This module defines variable data types.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kper
current stress period number
This module contains common time-varying property functionality.
subroutine, public tvbase_da(this)
Deallocate package memory.
This module contains time-varying conductivity package methods.
subroutine tvk_set_changed_at(this, kper, kstp)
Mark property changes as having occurred at (kper, kstp)
subroutine tvk_da(this)
Deallocate package memory.
subroutine tvk_ar_set_pointers(this)
Announce package and set pointers to variables.
subroutine tvk_apply_row_changes(this, n, node)
Apply input K/K22/K33 column changes for period-data row n to node.
subroutine tvk_validate_change(this, n, varName)
Check that a given property value is valid.
subroutine tvk_reset_change_flags(this)
Clear all per-node change flags.
subroutine, public tvk_cr(tvk, name_model, mempath, inunit, iout)
Create a new TvkType object.