MODFLOW 6  version 6.7.0.dev2
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...
 

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 315 of file gwt-src.f90.

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

◆ set_nodesontop()

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

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

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

◆ 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 169 of file gwt-src.f90.

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 

◆ src_allocate_scalars()

subroutine gwtsrcmodule::src_allocate_scalars ( class(gwtsrctype this)

Allocate scalars specific to this source loading package

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

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 

◆ 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 380 of file gwt-src.f90.

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

◆ src_cf()

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

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

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

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

◆ 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 50 of file gwt-src.f90.

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 
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 126 of file gwt-src.f90.

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)

◆ 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 359 of file gwt-src.f90.

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
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 280 of file gwt-src.f90.

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

◆ 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 345 of file gwt-src.f90.

346  implicit none
347  ! -- dummy
348  class(GwtSrcType) :: this
349  !
350  src_obs_supported = .true.

◆ src_options()

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

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

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'
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 203 of file gwt-src.f90.

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()
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 15 of file gwt-src.f90.

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

◆ text

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

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

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