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

Data Types

type  gwtsrctype
 

Functions/Subroutines

subroutine, public src_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype)
 Create a source loading package. More...
 
subroutine src_da (this)
 Deallocate memory. More...
 
subroutine src_allocate_scalars (this)
 Allocate scalars. 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...
 
subroutine src_rp_ts (this)
 Procedure related to time series. 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 183 of file gwt-src.f90.

184  ! -- dummy
185  class(GwtSrcType), intent(inout) :: this
186  ! -- local
187  !
188  ! -- create the header list label
189  this%listlabel = trim(this%filtyp)//' NO.'
190  if (this%dis%ndim == 3) then
191  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
192  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
193  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
194  elseif (this%dis%ndim == 2) then
195  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
196  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
197  else
198  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
199  end if
200  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
201  if (this%inamedbound == 1) then
202  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
203  end if

◆ src_allocate_scalars()

subroutine gwtsrcmodule::src_allocate_scalars ( class(gwtsrctype this)

Allocate scalars specific to this source loading package

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

104  ! -- dummy
105  class(GwtSrcType) :: this
106  !
107  ! -- call standard BndType allocate scalars
108  call this%BndType%allocate_scalars()
109  !
110  ! -- allocate the object and assign values to object variables
111  !
112  ! -- Set values

◆ src_cf()

subroutine gwtsrcmodule::src_cf ( class(gwtsrctype this)

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

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

122  ! -- dummy
123  class(GwtSrcType) :: this
124  ! -- local
125  integer(I4B) :: i, node
126  real(DP) :: q
127  !
128  ! -- Return if no sources
129  if (this%nbound == 0) return
130  !
131  ! -- Calculate hcof and rhs for each source entry
132  do i = 1, this%nbound
133  node = this%nodelist(i)
134  this%hcof(i) = dzero
135  if (this%ibound(node) <= 0) then
136  this%rhs(i) = dzero
137  cycle
138  end if
139  q = this%bound(1, i)
140  this%rhs(i) = -q
141  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=*), intent(in)  pakname,
character(len=lenvarname), intent(in)  depvartype 
)

This subroutine points bndobj to the newly created package

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

47  ! -- dummy
48  class(BndType), pointer :: packobj
49  integer(I4B), intent(in) :: id
50  integer(I4B), intent(in) :: ibcnum
51  integer(I4B), intent(in) :: inunit
52  integer(I4B), intent(in) :: iout
53  character(len=*), intent(in) :: namemodel
54  character(len=*), intent(in) :: pakname
55  character(len=LENVARNAME), intent(in) :: depvartype
56  ! -- local
57  type(GwtSrcType), pointer :: srcobj
58  !
59  ! -- allocate the object and assign values to object variables
60  allocate (srcobj)
61  packobj => srcobj
62  !
63  ! -- create name and memory path
64  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
65  packobj%text = text
66  !
67  ! -- allocate scalars
68  call srcobj%allocate_scalars()
69  !
70  ! -- initialize package
71  call packobj%pack_initialize()
72 
73  packobj%inunit = inunit
74  packobj%iout = iout
75  packobj%id = id
76  packobj%ibcnum = ibcnum
77  packobj%ncolbnd = 1
78  packobj%iscloc = 1
79  !
80  ! -- Store the appropriate label based on the dependent variable
81  srcobj%depvartype = depvartype
Here is the caller graph for this function:

◆ src_da()

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

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

87  ! -- modules
89  ! -- dummy
90  class(GwtSrcType) :: this
91  !
92  ! -- Deallocate parent package
93  call this%BndType%bnd_da()
94  !
95  ! -- scalars

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

228  implicit none
229  ! -- dummy
230  class(GwtSrcType) :: this
231  ! -- local
232  integer(I4B) :: indx
233  !
234  call this%obs%StoreObsType('src', .true., indx)
235  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
236  !
237  ! -- Store obs type and assign procedure pointer
238  ! for to-mvr observation type.
239  call this%obs%StoreObsType('to-mvr', .true., indx)
240  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 148 of file gwt-src.f90.

149  ! -- dummy
150  class(GwtSrcType) :: this
151  real(DP), dimension(:), intent(inout) :: rhs
152  integer(I4B), dimension(:), intent(in) :: ia
153  integer(I4B), dimension(:), intent(in) :: idxglo
154  class(MatrixBaseType), pointer :: matrix_sln
155  ! -- local
156  integer(I4B) :: i, n, ipos
157  !
158  ! -- pakmvrobj fc
159  if (this%imover == 1) then
160  call this%pakmvrobj%fc()
161  end if
162  !
163  ! -- Copy package rhs and hcof into solution rhs and amat
164  do i = 1, this%nbound
165  n = this%nodelist(i)
166  rhs(n) = rhs(n) + this%rhs(i)
167  ipos = ia(n)
168  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
169  !
170  ! -- If mover is active and mass is being withdrawn,
171  ! store available mass (as positive value).
172  if (this%imover == 1 .and. this%rhs(i) > dzero) then
173  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
174  end if
175  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 213 of file gwt-src.f90.

214  implicit none
215  ! -- dummy
216  class(GwtSrcType) :: this
217  !
218  src_obs_supported = .true.

◆ src_rp_ts()

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

Assign tsLinkText appropriately for all time series in use by package. In the SRC package only the SMASSRATE variable can be controlled by time series.

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

250  ! -- dummy
251  class(GwtSrcType), intent(inout) :: this
252  ! -- local
253  integer(I4B) :: i, nlinks
254  type(TimeSeriesLinkType), pointer :: tslink => null()
255  !
256  nlinks = this%TsManager%boundtslinks%Count()
257  do i = 1, nlinks
258  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
259  if (associated(tslink)) then
260  if (tslink%JCol == 1) then
261  tslink%Text = 'SMASSRATE'
262  end if
263  end if
264  end do
Here is the call graph for this function:

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'