MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwf-tvk.f90
Go to the documentation of this file.
1 !> @brief This module contains time-varying conductivity package methods
2 !!
3 !! This module contains the methods used to allow hydraulic conductivity
4 !! parameters in the NPF package (K11, K22, K33) to be varied throughout
5 !! a simulation.
6 !!
7 !<
8 module tvkmodule
9  use basedismodule, only: disbasetype
11  use kindmodule, only: i4b, dp
14  use simmodule, only: store_error
15  use simvariablesmodule, only: errmsg
16  use tdismodule, only: kper
18 
19  implicit none
20 
21  private
22 
23  public :: tvktype
24  public :: tvk_cr
25 
26  type, extends(tvbasetype) :: tvktype
27  integer(I4B), pointer :: ik22overk => null() !< NPF flag that k22 is specified as anisotropy ratio
28  integer(I4B), pointer :: ik33overk => null() !< NPF flag that k33 is specified as anisotropy ratio
29  real(dp), dimension(:), pointer, contiguous :: k11 => null() !< NPF hydraulic conductivity; if anisotropic, then this is Kx prior to rotation
30  real(dp), dimension(:), pointer, contiguous :: k22 => null() !< NPF hydraulic conductivity; if specified then this is Ky prior to rotation
31  real(dp), dimension(:), pointer, contiguous :: k33 => null() !< NPF hydraulic conductivity; if specified then this is Kz prior to rotation
32  integer(I4B), pointer :: kchangeper => null() !< NPF last stress period in which any node K (or K22, or K33) values were changed (0 if unchanged from start of simulation)
33  integer(I4B), pointer :: kchangestp => null() !< NPF last time step in which any node K (or K22, or K33) values were changed (0 if unchanged from start of simulation)
34  integer(I4B), dimension(:), pointer, contiguous :: nodekchange => null() !< NPF grid array of flags indicating for each node whether its K (or K22, or K33) value changed (1) at (kchangeper, kchangestp) or not (0)
35  real(dp), dimension(:), pointer, contiguous :: k11_src => null() !< input K values
36  real(dp), dimension(:), pointer, contiguous :: k22_src => null() !< input K22 values
37  real(dp), dimension(:), pointer, contiguous :: k33_src => null() !< input K33 values
38 
39  contains
40 
41  procedure :: da => tvk_da
47  end type tvktype
48 
49 contains
50 
51  !> @brief Create a new TvkType object
52  !!
53  !! Create a new time-varying conductivity (TvkType) object.
54  !<
55  subroutine tvk_cr(tvk, name_model, mempath, inunit, iout)
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)
65  end subroutine tvk_cr
66 
67  !> @brief Announce package and set pointers to variables
68  !!
69  !! Announce package version and set array and variable pointers from the NPF
70  !! package for access by TVK.
71  !<
72  subroutine tvk_ar_set_pointers(this)
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)
98  end subroutine tvk_ar_set_pointers
99 
100  !> @brief Apply input K/K22/K33 column changes for period-data row n to node.
101  !<
102  subroutine tvk_apply_row_changes(this, n, node)
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
147  end subroutine tvk_apply_row_changes
148 
149  !> @brief Mark property changes as having occurred at (kper, kstp)
150  !!
151  !! Deferred procedure implementation called by the TvBaseType code when a
152  !! property value change occurs at (kper, kstp).
153  !<
154  subroutine tvk_set_changed_at(this, kper, kstp)
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
162  end subroutine tvk_set_changed_at
163 
164  !> @brief Clear all per-node change flags
165  !!
166  !! Deferred procedure implementation called by the TvBaseType code when a
167  !! new time step commences, indicating that any previously set per-node
168  !! property value change flags should be reset.
169  !<
170  subroutine tvk_reset_change_flags(this)
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
179  end subroutine tvk_reset_change_flags
180 
181  !> @brief Check that a given property value is valid
182  !!
183  !! Deferred procedure implementation called by the TvBaseType code after a
184  !! property value change occurs. Check if the property value of the given
185  !! variable at the given node is invalid, and log an error if so. Update
186  !! K22 and K33 values appropriately when specified as anisotropy.
187  !<
188  subroutine tvk_validate_change(this, n, varName)
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
232  end subroutine tvk_validate_change
233 
234  !> @brief Deallocate package memory
235  !!
236  !! Deallocate TVK package scalars and arrays.
237  !<
238  subroutine tvk_da(this)
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)
254  end subroutine tvk_da
255 
256 end module tvkmodule
Apply input column changes for period-data row n to node.
Definition: TvBase.f90:63
Announce package and set pointers to variables.
Definition: TvBase.f90:54
Clear all per-node change flags.
Definition: TvBase.f90:94
Mark property changes as having occurred at (kper, kstp)
Definition: TvBase.f90:78
Check that a given property value is valid.
Definition: TvBase.f90:109
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
This module contains common time-varying property functionality.
Definition: TvBase.f90:8
subroutine, public tvbase_da(this)
Deallocate package memory.
Definition: TvBase.f90:302
This module contains time-varying conductivity package methods.
Definition: gwf-tvk.f90:8
subroutine tvk_set_changed_at(this, kper, kstp)
Mark property changes as having occurred at (kper, kstp)
Definition: gwf-tvk.f90:155
subroutine tvk_da(this)
Deallocate package memory.
Definition: gwf-tvk.f90:239
subroutine tvk_ar_set_pointers(this)
Announce package and set pointers to variables.
Definition: gwf-tvk.f90:73
subroutine tvk_apply_row_changes(this, n, node)
Apply input K/K22/K33 column changes for period-data row n to node.
Definition: gwf-tvk.f90:103
subroutine tvk_validate_change(this, n, varName)
Check that a given property value is valid.
Definition: gwf-tvk.f90:189
subroutine tvk_reset_change_flags(this)
Clear all per-node change flags.
Definition: gwf-tvk.f90:171
subroutine, public tvk_cr(tvk, name_model, mempath, inunit, iout)
Create a new TvkType object.
Definition: gwf-tvk.f90:56