MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
gwtsrcmodule Module Reference

Data Types

type  gwtsrctype
 

Functions/Subroutines

subroutine, public src_create (packobj, id, ibcnum, inunit, iout, namemodel, depvartype, pakname, input_mempath, fmi)
 Create a source loading package. More...
 
subroutine src_options (this)
 Set additional options specific to the GwtSrcType. More...
 
subroutine src_da (this)
 Deallocate memory. More...
 
subroutine src_allocate_scalars (this)
 Allocate scalars. More...
 
subroutine src_allocate_arrays (this, nodelist, auxvar)
 Allocate arrays. More...
 
subroutine src_rp (this)
 
subroutine set_nodesontop (this)
 Store nodelist in nodesontop. More...
 
subroutine src_cf (this)
 Formulate the HCOF and RHS terms. More...
 
subroutine src_fc (this, rhs, ia, idxglo, matrix_sln)
 Add matrix terms related to specified mass source loading. More...
 
subroutine define_listlabel (this)
 Define list labels. More...
 
logical function src_obs_supported (this)
 Support function for specified mass source loading observations. More...
 
subroutine src_df_obs (this)
 Define observations. More...
 
real(dp) function src_bound_value (this, col, row)
 @ brief Return a bound value More...
 
real(dp) function mass_mult (this, row)
 Return a value that applies a multiplier. More...
 

Variables

character(len=lenftype) ftype = 'SRC'
 
character(len=16) text = ' SRC'
 

Function/Subroutine Documentation

◆ define_listlabel()

subroutine gwtsrcmodule::define_listlabel ( class(gwtsrctype), intent(inout)  this)
private

Define the list heading that is written to iout when PRINT_INPUT option is used.

Definition at line 320 of file gwt-src.f90.

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

◆ mass_mult()

real(dp) function gwtsrcmodule::mass_mult ( class(gwtsrctype), intent(inout)  this,
integer(i4b), intent(in)  row 
)

Apply multiplier to specified mass-source load depending on user-selected option

Parameters
[in,out]thisBndExtType object

Definition at line 407 of file gwt-src.f90.

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
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ set_nodesontop()

subroutine gwtsrcmodule::set_nodesontop ( class(gwtsrctype), intent(inout)  this)

Definition at line 217 of file gwt-src.f90.

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

◆ src_allocate_arrays()

subroutine gwtsrcmodule::src_allocate_arrays ( class(gwtsrctype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)

Allocate scalars specific to this source loading package

Parameters
nodelistpackage nodelist
auxvarpackage aux variable array

Definition at line 172 of file gwt-src.f90.

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

◆ src_allocate_scalars()

subroutine gwtsrcmodule::src_allocate_scalars ( class(gwtsrctype this)

Allocate scalars specific to this source loading package

Definition at line 152 of file gwt-src.f90.

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 

◆ src_bound_value()

real(dp) function gwtsrcmodule::src_bound_value ( class(gwtsrctype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private

Return a bound value associated with an ncolbnd index and row.

Definition at line 385 of file gwt-src.f90.

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

◆ src_cf()

subroutine gwtsrcmodule::src_cf ( class(gwtsrctype this)
private

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

Definition at line 241 of file gwt-src.f90.

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

◆ src_create()

subroutine, public gwtsrcmodule::src_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=lenvarname), intent(in)  depvartype,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  input_mempath,
type(tspfmitype), intent(in), target  fmi 
)

This subroutine points bndobj to the newly created package

Definition at line 53 of file gwt-src.f90.

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 
This module contains the base boundary package.
@ brief BndType
Here is the caller graph for this function:

◆ src_da()

subroutine gwtsrcmodule::src_da ( class(gwtsrctype this)

Definition at line 129 of file gwt-src.f90.

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)

◆ src_df_obs()

subroutine gwtsrcmodule::src_df_obs ( class(gwtsrctype this)
private

This subroutine:

  • stores observation types supported by SRC package.
  • overrides BndTypebnd_df_obs

Definition at line 364 of file gwt-src.f90.

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
Here is the call graph for this function:

◆ src_fc()

subroutine gwtsrcmodule::src_fc ( class(gwtsrctype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)
private

Copy rhs and hcof into solution rhs and amat

Definition at line 285 of file gwt-src.f90.

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

◆ src_obs_supported()

logical function gwtsrcmodule::src_obs_supported ( class(gwtsrctype this)
private

This function:

  • returns true because SRC package supports observations.
  • overrides BndTypebnd_obs_supported()

Definition at line 350 of file gwt-src.f90.

351  implicit none
352  ! -- dummy
353  class(GwtSrcType) :: this
354  !
355  src_obs_supported = .true.

◆ src_options()

subroutine gwtsrcmodule::src_options ( class(gwtsrctype), intent(inout)  this)

Definition at line 101 of file gwt-src.f90.

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'
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23

◆ src_rp()

subroutine gwtsrcmodule::src_rp ( class(gwtsrctype), intent(inout)  this)

Definition at line 205 of file gwt-src.f90.

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()
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

Variable Documentation

◆ ftype

character(len=lenftype) gwtsrcmodule::ftype = 'SRC'
private

Definition at line 17 of file gwt-src.f90.

17  character(len=LENFTYPE) :: ftype = 'SRC'

◆ text

character(len=16) gwtsrcmodule::text = ' SRC'
private

Definition at line 18 of file gwt-src.f90.

18  character(len=16) :: text = ' SRC'