MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
TvBase.f90
Go to the documentation of this file.
1 !> @brief This module contains common time-varying property functionality
2 !!
3 !! This module contains methods implementing functionality common to both
4 !! time-varying hydraulic conductivity (TVK) and time-varying storage (TVS)
5 !! packages.
6 !!
7 !<
9  use basedismodule, only: disbasetype
11  use geomutilmodule, only: get_node
12  use kindmodule, only: i4b, dp, lgp
15  use simvariablesmodule, only: errmsg
16  use tdismodule, only: kper, nper, kstp
18 
19  implicit none
20 
21  private
22 
23  public :: tvbasetype
24  public :: tvbase_da
25 
26  type, abstract, extends(numericalpackagetype) :: tvbasetype
27  integer(I4B), dimension(:, :), pointer, contiguous :: cellid => null() !< input cellids of time varying value
28  logical(LGP) :: ts_active = .false. !< is timeseries active in package
29  contains
30  procedure :: init
31  procedure :: ar
32  procedure :: rp
33  procedure :: ad
34  procedure :: da => tvbase_da
35  procedure, private :: tvbase_allocate_scalars
36  procedure, private :: tv_get_node
37  procedure(ar_set_pointers), deferred :: ar_set_pointers
38  procedure :: source_options => tvbase_source_options
39  procedure :: source_package_options => tvbase_source_package_options
40  procedure(apply_row_changes), deferred :: apply_row_changes
41  procedure(set_changed_at), deferred :: set_changed_at
42  procedure(reset_change_flags), deferred :: reset_change_flags
43  procedure(validate_change), deferred :: validate_change
44  end type tvbasetype
45 
46  abstract interface
47 
48  !> @brief Announce package and set pointers to variables
49  !!
50  !! Deferred procedure called by the TvBaseType code to announce the
51  !! specific package version and set any required array and variable
52  !! pointers from other packages.
53  !<
54  subroutine ar_set_pointers(this)
55  ! -- modules
56  import tvbasetype
57  ! -- dummy
58  class(tvbasetype) :: this
59  end subroutine
60 
61  !> @brief Apply input column changes for period-data row n to node.
62  !<
63  subroutine apply_row_changes(this, n, node)
64  ! -- modules
65  use kindmodule, only: i4b
66  import tvbasetype
67  ! -- dummy
68  class(tvbasetype) :: this
69  integer(I4B), intent(in) :: n !< row index in period data arrays
70  integer(I4B), intent(in) :: node !< reduced node number
71  end subroutine
72 
73  !> @brief Mark property changes as having occurred at (kper, kstp)
74  !!
75  !! Deferred procedure called by the TvBaseType code when a property value
76  !! change occurs at (kper, kstp).
77  !<
78  subroutine set_changed_at(this, kper, kstp)
79  ! -- modules
80  use kindmodule, only: i4b
81  import tvbasetype
82  ! -- dummy
83  class(tvbasetype) :: this
84  integer(I4B), intent(in) :: kper
85  integer(I4B), intent(in) :: kstp
86  end subroutine
87 
88  !> @brief Clear all per-node change flags
89  !!
90  !! Deferred procedure called by the TvBaseType code when a new time step
91  !! commences, indicating that any previously set per-node property value
92  !! change flags should be reset.
93  !<
94  subroutine reset_change_flags(this)
95  ! -- modules
96  import tvbasetype
97  ! -- dummy
98  class(tvbasetype) :: this
99  end subroutine
100 
101  !> @brief Check that a given property value is valid
102  !!
103  !! Deferred procedure called by each derived type's apply_row_changes
104  !! implementation after a property value change occurs. Performs any
105  !! required validity checks on the value of the given variable at the
106  !! given node. Perform any required updates to the property value if it
107  !! is valid, or log an error if not.
108  !<
109  subroutine validate_change(this, n, varName)
110  ! -- modules
111  use kindmodule, only: i4b
112  import tvbasetype
113  ! -- dummy
114  class(tvbasetype) :: this
115  integer(I4B), intent(in) :: n
116  character(len=*), intent(in) :: varName
117  end subroutine
118 
119  end interface
120 
121 contains
122 
123  !> @brief Initialize the TvBaseType object
124  !!
125  !! Allocate and initialize data members of the object.
126  !<
127  subroutine init(this, name_model, pakname, ftype, mempath, inunit, iout)
128  ! -- dummy
129  class(tvbasetype) :: this
130  character(len=*), intent(in) :: name_model
131  character(len=*), intent(in) :: pakname
132  character(len=*), intent(in) :: ftype
133  character(len=*), intent(in) :: mempath
134  integer(I4B), intent(in) :: inunit
135  integer(I4B), intent(in) :: iout
136  !
137  call this%set_names(1, name_model, pakname, ftype, mempath)
138  call this%tvbase_allocate_scalars()
139  this%inunit = inunit
140  this%iout = iout
141  end subroutine init
142 
143  !> @brief Allocate scalar variables
144  !!
145  !! Allocate scalar data members of the object.
146  !<
147  subroutine tvbase_allocate_scalars(this)
148  ! -- dummy
149  class(tvbasetype) :: this
150  !
151  ! -- Call standard NumericalPackageType allocate scalars
152  call this%NumericalPackageType%allocate_scalars()
153  end subroutine tvbase_allocate_scalars
154 
155  !> @brief Source common options from the input memory path.
156  !!
157  !! Source common options and call derived package routine to source
158  !! and log any package-specific options within the same block.
159  !<
160  subroutine tvbase_source_options(this)
161  ! -- modules
163  ! -- dummy
164  class(tvbasetype) :: this
165  ! -- locals
166  integer(I4B) :: isize
167  logical(LGP) :: found_print_input
168  !
169  write (this%iout, '(1x,a)') &
170  'PROCESSING '//trim(adjustl(this%packName))//' OPTIONS'
171  !
172  call mem_set_value(this%iprpak, 'PRINT_INPUT', this%input_mempath, &
173  found_print_input)
174  !
175  if (found_print_input) then
176  write (this%iout, '(4x,a)') 'TIME-VARYING INPUT WILL BE PRINTED.'
177  end if
178  !
179  call get_isize('TS6_FILENAME', this%input_mempath, isize)
180  if (isize > 0) this%ts_active = .true.
181  !
182  ! -- source package-specific options
183  call this%source_package_options()
184  !
185  write (this%iout, '(1x,a)') &
186  'END OF '//trim(adjustl(this%packName))//' OPTIONS'
187  end subroutine tvbase_source_options
188 
189  !> @brief Source package-specific options from the input memory path.
190  !!
191  !! Override in derived package to source and log package-specific options.
192  !! Default implementation is a no-op.
193  !<
195  ! -- dummy
196  class(tvbasetype) :: this
197  !
198  ! -- no package-specific options in the base class
199  end subroutine tvbase_source_package_options
200 
201  !> @brief Allocate and read static data for the package.
202  !<
203  subroutine ar(this, dis)
204  ! -- dummy
205  class(tvbasetype) :: this
206  class(disbasetype), pointer, intent(in) :: dis
207  !
208  this%dis => dis
209  call this%ar_set_pointers()
210  !
211  call this%source_options()
212  !
213  ! -- set input mempath pointers
214  call mem_setptr(this%cellid, 'CELLID', this%input_mempath)
215  !
216  if (count_errors() > 0) then
217  call store_error_filename(this%input_fname)
218  end if
219  end subroutine ar
220 
221  !> @brief Read and prepare stress period data for the package.
222  !<
223  subroutine rp(this)
224  ! -- dummy
225  class(tvbasetype) :: this
226  ! -- local variables
227  integer(I4B), pointer :: iper, nbound
228  integer(I4B) :: n, node
229  character(len=LINELENGTH) :: cellstr
230  !
231  ! -- check last loaded input period
232  call mem_setptr(iper, 'IPER', this%input_mempath)
233  if (iper /= kper) return
234  !
235  ! -- When timeseries are active, ad applies values at every time step
236  if (this%ts_active) return
237  !
238  call mem_setptr(nbound, 'NBOUND', this%input_mempath)
239  !
240  ! -- Reset per-node property change flags
241  call this%reset_change_flags()
242  !
243  do n = 1, nbound
244  node = this%tv_get_node(n)
245  if (node < 1 .or. node > this%dis%nodes) then
246  call this%dis%noder_to_string(node, cellstr)
247  write (errmsg, '(a,2(1x,a))') &
248  'CELLID', trim(cellstr), 'is not in the active model domain.'
249  call store_error(errmsg)
250  cycle
251  end if
252  !
253  call this%apply_row_changes(n, node)
254  end do
255  !
256  ! -- Record that changes were made at the current stress period / time step
257  if (nbound > 0) then
258  call this%set_changed_at(kper, kstp)
259  end if
260  !
261  if (count_errors() > 0) then
262  call store_error_filename(this%input_fname)
263  end if
264  end subroutine rp
265 
266  !> @brief Apply advanced values at each time step.
267  !<
268  subroutine ad(this)
269  ! -- dummy
270  class(tvbasetype) :: this
271  ! -- local variables
272  integer(I4B), pointer :: nbound
273  integer(I4B) :: n, node
274  !
275  ! -- no-op when timeseries aren't active.
276  if (.not. this%ts_active) return
277  !
278  call mem_setptr(nbound, 'NBOUND', this%input_mempath)
279  if (nbound > 0) then
280  ! -- Record that changes were made at the current time step
281  call this%set_changed_at(kper, kstp)
282  ! -- Reset node change flags
283  call this%reset_change_flags()
284  ! -- Apply row changes
285  do n = 1, nbound
286  node = this%tv_get_node(n)
287  if (node < 1 .or. node > this%dis%nodes) cycle
288  call this%apply_row_changes(n, node)
289  end do
290  end if
291  !
292  if (count_errors() > 0) then
293  call store_error_filename(this%input_fname)
294  end if
295  end subroutine ad
296 
297  !> @brief Deallocate package memory
298  !!
299  !! Deallocate package scalars and arrays.
300  !<
301  subroutine tvbase_da(this)
302  ! -- dummy
303  class(tvbasetype) :: this
304  !
305  nullify (this%cellid)
306  call this%NumericalPackageType%da()
307  end subroutine tvbase_da
308 
309  !> @brief Return the reduced node number for CELLID row n; caller must bounds-check.
310  !<
311  function tv_get_node(this, n) result(node)
312  ! -- dummy
313  class(tvbasetype) :: this
314  integer(I4B), intent(in) :: n !< row index in the period-data arrays
315  ! -- return
316  integer(I4B) :: node
317  ! -- local
318  integer(I4B) :: nodeu
319  !
320  if (this%dis%ndim == 1) then
321  nodeu = this%cellid(1, n)
322  elseif (this%dis%ndim == 2) then
323  nodeu = get_node(this%cellid(1, n), 1, &
324  this%cellid(2, n), &
325  this%dis%mshape(1), 1, &
326  this%dis%mshape(2))
327  else
328  nodeu = get_node(this%cellid(1, n), &
329  this%cellid(2, n), &
330  this%cellid(3, n), &
331  this%dis%mshape(1), &
332  this%dis%mshape(2), &
333  this%dis%mshape(3))
334  end if
335  node = this%dis%get_nodenumber(nodeu, 1)
336  end function tv_get_node
337 
338 end module tvbasemodule
subroutine init()
Definition: GridSorting.f90:25
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 dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:83
This module defines variable data types.
Definition: kind.f90:8
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the base numerical package type.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21
This module contains common time-varying property functionality.
Definition: TvBase.f90:8
subroutine tvbase_source_options(this)
Source common options from the input memory path.
Definition: TvBase.f90:161
integer(i4b) function tv_get_node(this, n)
Return the reduced node number for CELLID row n; caller must bounds-check.
Definition: TvBase.f90:312
subroutine, public tvbase_da(this)
Deallocate package memory.
Definition: TvBase.f90:302
subroutine tvbase_allocate_scalars(this)
Allocate scalar variables.
Definition: TvBase.f90:148
subroutine rp(this)
Read and prepare stress period data for the package.
Definition: TvBase.f90:224
subroutine tvbase_source_package_options(this)
Source package-specific options from the input memory path.
Definition: TvBase.f90:195
subroutine ad(this)
Apply advanced values at each time step.
Definition: TvBase.f90:269
subroutine ar(this, dis)
Allocate and read static data for the package.
Definition: TvBase.f90:204