MODFLOW 6  version 6.7.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
17 
18  implicit none
19 
20  private
21 
22  public :: tvktype
23  public :: tvk_cr
24 
25  type, extends(tvbasetype) :: tvktype
26  integer(I4B), pointer :: ik22overk => null() !< NPF flag that k22 is specified as anisotropy ratio
27  integer(I4B), pointer :: ik33overk => null() !< NPF flag that k33 is specified as anisotropy ratio
28  real(dp), dimension(:), pointer, contiguous :: k11 => null() !< NPF hydraulic conductivity; if anisotropic, then this is Kx prior to rotation
29  real(dp), dimension(:), pointer, contiguous :: k22 => null() !< NPF hydraulic conductivity; if specified then this is Ky prior to rotation
30  real(dp), dimension(:), pointer, contiguous :: k33 => null() !< NPF hydraulic conductivity; if specified then this is Kz prior to rotation
31  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)
32  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)
33  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)
34 
35  contains
36 
37  procedure :: da => tvk_da
39  procedure :: read_option => tvk_read_option
44  end type tvktype
45 
46 contains
47 
48  !> @brief Create a new TvkType object
49  !!
50  !! Create a new time-varying conductivity (TvkType) object.
51  !<
52  subroutine tvk_cr(tvk, name_model, inunit, iout)
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  end subroutine tvk_cr
62 
63  !> @brief Announce package and set pointers to variables
64  !!
65  !! Announce package version and set array and variable pointers from the NPF
66  !! package for access by TVK.
67  !<
68  subroutine tvk_ar_set_pointers(this)
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)
92  end subroutine tvk_ar_set_pointers
93 
94  !> @brief Read a TVK-specific option from the OPTIONS block
95  !!
96  !! Process a single TVK-specific option. Used when reading the OPTIONS block
97  !! of the TVK package input file.
98  !<
99  function tvk_read_option(this, keyword) result(success)
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.
108  end function tvk_read_option
109 
110  !> @brief Get an array value pointer given a variable name and node index
111  !!
112  !! Return a pointer to the given node's value in the appropriate NPF array
113  !! based on the given variable name string.
114  !<
115  function tvk_get_pointer_to_value(this, n, varName) result(bndElem)
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
133  end function tvk_get_pointer_to_value
134 
135  !> @brief Mark property changes as having occurred at (kper, kstp)
136  !!
137  !! Deferred procedure implementation called by the TvBaseType code when a
138  !! property value change occurs at (kper, kstp).
139  !<
140  subroutine tvk_set_changed_at(this, kper, kstp)
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
148  end subroutine tvk_set_changed_at
149 
150  !> @brief Clear all per-node change flags
151  !!
152  !! Deferred procedure implementation called by the TvBaseType code when a
153  !! new time step commences, indicating that any previously set per-node
154  !! property value change flags should be reset.
155  !<
156  subroutine tvk_reset_change_flags(this)
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
166  end subroutine tvk_reset_change_flags
167 
168  !> @brief Check that a given property value is valid
169  !!
170  !! Deferred procedure implementation called by the TvBaseType code after a
171  !! property value change occurs. Check if the property value of the given
172  !! variable at the given node is invalid, and log an error if so. Update
173  !! K22 and K33 values appropriately when specified as anisotropy.
174  !<
175  subroutine tvk_validate_change(this, n, varName)
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
219  end subroutine tvk_validate_change
220 
221  !> @brief Deallocate package memory
222  !!
223  !! Deallocate TVK package scalars and arrays.
224  !<
225  subroutine tvk_da(this)
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)
241  end subroutine tvk_da
242 
243 end module tvkmodule
Announce package and set pointers to variables.
Definition: TvBase.f90:54
Get an array value pointer given a variable name and node index.
Definition: TvBase.f90:83
Announce package and set pointers to variables.
Definition: TvBase.f90:67
Clear all per-node change flags.
Definition: TvBase.f90:118
Mark property changes as having occurred at (kper, kstp)
Definition: TvBase.f90:101
Check that a given property value is valid.
Definition: TvBase.f90:133
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 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
This module contains common time-varying property functionality.
Definition: TvBase.f90:8
subroutine, public tvbase_da(this)
Deallocate package memory.
Definition: TvBase.f90:460
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:141
subroutine tvk_da(this)
Deallocate package memory.
Definition: gwf-tvk.f90:226
logical function tvk_read_option(this, keyword)
Read a TVK-specific option from the OPTIONS block.
Definition: gwf-tvk.f90:100
subroutine, public tvk_cr(tvk, name_model, inunit, iout)
Create a new TvkType object.
Definition: gwf-tvk.f90:53
subroutine tvk_ar_set_pointers(this)
Announce package and set pointers to variables.
Definition: gwf-tvk.f90:69
real(dp) function, pointer tvk_get_pointer_to_value(this, n, varName)
Get an array value pointer given a variable name and node index.
Definition: gwf-tvk.f90:116
subroutine tvk_validate_change(this, n, varName)
Check that a given property value is valid.
Definition: gwf-tvk.f90:176
subroutine tvk_reset_change_flags(this)
Clear all per-node change flags.
Definition: gwf-tvk.f90:157