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