MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
gwt-src.f90
Go to the documentation of this file.
2  !
3  use kindmodule, only: dp, i4b, lgp
5  use simvariablesmodule, only: errmsg
7  use tspfmimodule, only: tspfmitype
8  use bndextmodule, only: bndexttype
11  !
12  implicit none
13  !
14  private
15  public :: src_create
16  !
17  character(len=LENFTYPE) :: ftype = 'SRC'
18  character(len=16) :: text = ' SRC'
19  !
20  type, extends(bndexttype) :: gwtsrctype
21 
22  character(len=LENVARNAME) :: depvartype = '' !< stores string of dependent variable type, depending on model type
23  type(tspfmitype), pointer :: fmi => null() ! pointer to GWE fmi object
24  logical(LGP), pointer :: highest_sat => null()
25  real(dp), dimension(:), pointer, contiguous :: smassrate => null() !< mass source loading rate
26  integer(I4B), dimension(:), pointer, contiguous :: nodesontop => null() ! User provided cell numbers; nodelist is cells where recharge is applied)
27 
28  contains
29 
30  procedure :: allocate_scalars => src_allocate_scalars
31  procedure :: allocate_arrays => src_allocate_arrays
32  procedure :: source_options => src_options
33  procedure :: bnd_rp => src_rp
34  procedure :: bnd_cf => src_cf
35  procedure :: bnd_fc => src_fc
36  procedure :: bnd_da => src_da
37  procedure :: define_listlabel
38  procedure :: set_nodesontop
39  procedure :: bound_value => src_bound_value
40  procedure :: mass_mult
41  ! -- methods for observations
42  procedure, public :: bnd_obs_supported => src_obs_supported
43  procedure, public :: bnd_df_obs => src_df_obs
44 
45  end type gwtsrctype
46 
47 contains
48 
49  !> @brief Create a source loading package
50  !!
51  !! This subroutine points bndobj to the newly created package
52  !<
53  subroutine src_create(packobj, id, ibcnum, inunit, iout, namemodel, &
54  depvartype, pakname, input_mempath, fmi)
55  ! -- modules
56  use bndmodule, only: bndtype
57  ! -- dummy
58  class(bndtype), pointer :: packobj
59  integer(I4B), intent(in) :: id
60  integer(I4B), intent(in) :: ibcnum
61  integer(I4B), intent(in) :: inunit
62  integer(I4B), intent(in) :: iout
63  character(len=*), intent(in) :: namemodel
64  character(len=*), intent(in) :: pakname
65  character(len=LENVARNAME), intent(in) :: depvartype
66  character(len=*), intent(in) :: input_mempath
67  type(tspfmitype), intent(in), target :: fmi
68  ! -- local
69  type(gwtsrctype), pointer :: srcobj
70  !
71  ! -- allocate the object and assign values to object variables
72  allocate (srcobj)
73  packobj => srcobj
74  !
75  ! -- create name and memory path
76  call packobj%set_names(ibcnum, namemodel, pakname, ftype, input_mempath)
77  packobj%text = text
78  !
79  ! -- allocate scalars
80  call srcobj%allocate_scalars()
81  !
82  ! -- initialize package
83  call packobj%pack_initialize()
84 
85  packobj%inunit = inunit
86  packobj%iout = iout
87  packobj%id = id
88  packobj%ibcnum = ibcnum
89  packobj%ncolbnd = 1
90  packobj%iscloc = 1
91  !
92  ! -- Store the appropriate label based on the dependent variable
93  srcobj%depvartype = depvartype
94 
95  srcobj%fmi => fmi
96 
97  end subroutine src_create
98 
99  !> @brief Set additional options specific to the GwtSrcType
100  !<
101  subroutine src_options(this)
102  ! -- modules
106  ! -- dummy
107  class(gwtsrctype), intent(inout) :: this
108  ! -- local
109  type(gwtsrcparamfoundtype) :: found
110  !
111  ! -- source base class options
112  call this%BndExtType%source_options()
113  !
114  ! -- source options from input context
115  call mem_set_value(this%highest_sat, 'HIGHEST_SAT', this%input_mempath, &
116  found%highest_sat)
117 
118  write (this%iout, '(/1x,a)') 'PROCESSING SRC OPTIONS'
119  if (found%highest_sat) then
120  write (this%iout, '(4x,a)') &
121  'Mass source loading rate will be applied to the highest cell at or below &
122  &the specified cellid with a non-zero saturation.'
123  end if
124  write (this%iout, '(1x,a)') 'END OF SRC OPTIONS'
125  end subroutine src_options
126 
127  !> @brief Deallocate memory
128  !<
129  subroutine src_da(this)
130  ! -- modules
132  ! -- dummy
133  class(gwtsrctype) :: this
134  !
135  ! -- Deallocate parent package
136  call this%BndExtType%bnd_da()
137  !
138  ! -- arrays
139  call mem_deallocate(this%smassrate, 'SMASSRATE', this%memoryPath)
140  if (this%highest_sat) then
141  call mem_deallocate(this%nodesontop, "NODESONTOP", this%memoryPath)
142  end if
143  !
144  ! -- scalars
145  call mem_deallocate(this%highest_sat)
146  end subroutine src_da
147 
148  !> @brief Allocate scalars
149  !!
150  !! Allocate scalars specific to this source loading package
151  !<
152  subroutine src_allocate_scalars(this)
154  ! -- dummy
155  class(gwtsrctype) :: this
156  !
157  ! -- base class allocate scalars
158  call this%BndExtType%allocate_scalars()
159  !
160  ! -- allocate the object and assign values to object variables
161  call mem_allocate(this%highest_sat, 'HIGHEST_SAT', this%memoryPath)
162  !
163  ! -- Set values
164  this%highest_sat = .false.
165 
166  end subroutine src_allocate_scalars
167 
168  !> @brief Allocate arrays
169  !!
170  !! Allocate scalars specific to this source loading package
171  !<
172  subroutine src_allocate_arrays(this, nodelist, auxvar)
174  ! -- dummy
175  class(gwtsrctype) :: this
176  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist !< package nodelist
177  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar !< package aux variable array
178  ! local
179  integer(I4B) :: n
180  !
181  ! -- base class allocate arrays
182  call this%BndExtType%allocate_arrays(nodelist, auxvar)
183  !
184  ! -- allocate the object and assign values to object variables
185  if (this%highest_sat) then
186  call mem_allocate(this%nodesontop, this%maxbound, 'NODESONTOP', &
187  this%memoryPath)
188  end if
189 
190  ! -- set input context pointers
191  call mem_setptr(this%smassrate, 'SMASSRATE', this%input_mempath)
192  !
193  ! -- checkin input context pointers
194  call mem_checkin(this%smassrate, 'SMASSRATE', this%memoryPath, &
195  'SMASSRATE', this%input_mempath)
196  !
197  ! -- Set values
198  if (this%highest_sat) then
199  do n = 1, this%maxbound
200  this%nodesontop(n) = 0
201  end do
202  end if
203  end subroutine src_allocate_arrays
204 
205  subroutine src_rp(this)
206  ! -- modules
207  use tdismodule, only: kper
208  ! -- dummy
209  class(gwtsrctype), intent(inout) :: this
210  if (this%iper /= kper) return
211  call this%BndExtType%bnd_rp()
212  if (this%highest_sat) call this%set_nodesontop()
213  end subroutine src_rp
214 
215  !> @brief Store nodelist in nodesontop
216  !<
217  subroutine set_nodesontop(this)
218  implicit none
219  ! -- dummy
220  class(gwtsrctype), intent(inout) :: this
221  ! -- local
222  integer(I4B) :: n
223  ! !
224  ! ! -- allocate if necessary
225  ! if (.not. associated(this%nodesontop)) then
226  ! allocate (this%nodesontop(this%maxbound))
227  ! end if
228  !
229  ! -- copy nodelist into nodesontop
230  do n = 1, this%nbound
231  this%nodesontop(n) = this%nodelist(n)
232  end do
233  end subroutine set_nodesontop
234 
235  !> @brief Formulate the HCOF and RHS terms
236  !!
237  !! This subroutine:
238  !! - calculates hcof and rhs terms
239  !! - skip if no sources
240  !<
241  subroutine src_cf(this)
242  ! -- dummy
243  class(gwtsrctype) :: this
244  ! -- local
245  integer(I4B) :: i, node
246  real(DP) :: q
247  !
248  ! -- Return if no sources
249  if (this%nbound == 0) return
250  !
251  ! -- Calculate hcof and rhs for each source entry
252  do i = 1, this%nbound
253  !
254  ! -- Find the node number
255  if (this%highest_sat) then
256  node = this%nodesontop(i)
257  else
258  node = this%nodelist(i)
259  end if
260  !
261  ! -- reset nodelist to highest active
262  if (this%highest_sat) then
263  if (this%fmi%gwfsat(node) == 0) &
264  call this%dis%highest_saturated(node, this%fmi%gwfsat)
265  this%nodelist(i) = node
266  end if
267 
268  this%hcof(i) = dzero
269  if (this%ibound(node) <= 0) then
270  this%rhs(i) = dzero
271  cycle
272  end if
273  !
274  ! -- set energy loading rate accounting for multiplier
275  q = this%mass_mult(i)
276  !
277  this%rhs(i) = -q
278  end do
279  end subroutine src_cf
280 
281  !> @brief Add matrix terms related to specified mass source loading
282  !!
283  !! Copy rhs and hcof into solution rhs and amat
284  !<
285  subroutine src_fc(this, rhs, ia, idxglo, matrix_sln)
286  ! -- dummy
287  class(gwtsrctype) :: this
288  real(DP), dimension(:), intent(inout) :: rhs
289  integer(I4B), dimension(:), intent(in) :: ia
290  integer(I4B), dimension(:), intent(in) :: idxglo
291  class(matrixbasetype), pointer :: matrix_sln
292  ! -- local
293  integer(I4B) :: i, n, ipos
294  !
295  ! -- pakmvrobj fc
296  if (this%imover == 1) then
297  call this%pakmvrobj%fc()
298  end if
299  !
300  ! -- Copy package rhs and hcof into solution rhs and amat
301  do i = 1, this%nbound
302  n = this%nodelist(i)
303  rhs(n) = rhs(n) + this%rhs(i)
304  ipos = ia(n)
305  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
306  !
307  ! -- If mover is active and mass is being withdrawn,
308  ! store available mass (as positive value).
309  if (this%imover == 1 .and. this%rhs(i) > dzero) then
310  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
311  end if
312  end do
313  end subroutine src_fc
314 
315  !> @brief Define list labels
316  !!
317  !! Define the list heading that is written to iout when PRINT_INPUT
318  !! option is used.
319  !<
320  subroutine define_listlabel(this)
321  ! -- dummy
322  class(gwtsrctype), intent(inout) :: this
323  ! -- local
324  !
325  ! -- create the header list label
326  this%listlabel = trim(this%filtyp)//' NO.'
327  if (this%dis%ndim == 3) then
328  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
329  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
330  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
331  elseif (this%dis%ndim == 2) then
332  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
333  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
334  else
335  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
336  end if
337  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
338  if (this%inamedbound == 1) then
339  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
340  end if
341  end subroutine define_listlabel
342 
343  ! -- Procedures related to observations
344  !> @brief Support function for specified mass source loading observations
345  !!
346  !! This function:
347  !! - returns true because SRC package supports observations.
348  !! - overrides BndType%bnd_obs_supported()
349  !<
350  logical function src_obs_supported(this)
351  implicit none
352  ! -- dummy
353  class(gwtsrctype) :: this
354  !
355  src_obs_supported = .true.
356  end function src_obs_supported
357 
358  !> @brief Define observations
359  !!
360  !! This subroutine:
361  !! - stores observation types supported by SRC package.
362  !! - overrides BndType%bnd_df_obs
363  !<
364  subroutine src_df_obs(this)
365  implicit none
366  ! -- dummy
367  class(gwtsrctype) :: this
368  ! -- local
369  integer(I4B) :: indx
370  !
371  call this%obs%StoreObsType('src', .true., indx)
372  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
373  !
374  ! -- Store obs type and assign procedure pointer
375  ! for to-mvr observation type.
376  call this%obs%StoreObsType('to-mvr', .true., indx)
377  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
378  end subroutine src_df_obs
379 
380  !> @ brief Return a bound value
381  !!
382  !! Return a bound value associated with an ncolbnd index
383  !! and row.
384  !<
385  function src_bound_value(this, col, row) result(bndval)
386  ! -- modules
387  use constantsmodule, only: dzero
388  ! -- dummy variables
389  class(gwtsrctype), intent(inout) :: this
390  integer(I4B), intent(in) :: col
391  integer(I4B), intent(in) :: row
392  ! -- result
393  real(dp) :: bndval
394  !
395  select case (col)
396  case (1)
397  bndval = this%smassrate(row)
398  case default
399  end select
400  end function src_bound_value
401 
402  !> @brief Return a value that applies a multiplier
403  !!
404  !! Apply multiplier to specified mass-source load depending on user-selected
405  !! option
406  !<
407  function mass_mult(this, row) result(ener)
408  ! -- modules
409  use constantsmodule, only: dzero
410  ! -- dummy variables
411  class(gwtsrctype), intent(inout) :: this !< BndExtType object
412  integer(I4B), intent(in) :: row
413  ! -- result
414  real(dp) :: ener
415  !
416  if (this%iauxmultcol > 0) then
417  ener = this%smassrate(row) * this%auxvar(this%iauxmultcol, row)
418  else
419  ener = this%smassrate(row)
420  end if
421  end function mass_mult
422 end module gwtsrcmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:103
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public src_create(packobj, id, ibcnum, inunit, iout, namemodel, depvartype, pakname, input_mempath, fmi)
Create a source loading package.
Definition: gwt-src.f90:55
character(len=lenftype) ftype
Definition: gwt-src.f90:17
subroutine src_fc(this, rhs, ia, idxglo, matrix_sln)
Add matrix terms related to specified mass source loading.
Definition: gwt-src.f90:286
subroutine src_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwt-src.f90:242
character(len=16) text
Definition: gwt-src.f90:18
subroutine src_allocate_arrays(this, nodelist, auxvar)
Allocate arrays.
Definition: gwt-src.f90:173
subroutine src_options(this)
Set additional options specific to the GwtSrcType.
Definition: gwt-src.f90:102
subroutine set_nodesontop(this)
Store nodelist in nodesontop.
Definition: gwt-src.f90:218
subroutine src_allocate_scalars(this)
Allocate scalars.
Definition: gwt-src.f90:153
subroutine src_da(this)
Deallocate memory.
Definition: gwt-src.f90:130
subroutine define_listlabel(this)
Define list labels.
Definition: gwt-src.f90:321
subroutine src_rp(this)
Definition: gwt-src.f90:206
real(dp) function mass_mult(this, row)
Return a value that applies a multiplier.
Definition: gwt-src.f90:408
subroutine src_df_obs(this)
Define observations.
Definition: gwt-src.f90:365
real(dp) function src_bound_value(this, col, row)
@ brief Return a bound value
Definition: gwt-src.f90:386
logical function src_obs_supported(this)
Support function for specified mass source loading observations.
Definition: gwt-src.f90:351
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
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 kper
current stress period number
Definition: tdis.f90:23
@ brief BndType
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23