MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 198 of file gwt-src.f90.

199  ! -- dummy
200  class(GwtSrcType), intent(inout) :: this
201  ! -- local
202  !
203  ! -- create the header list label
204  this%listlabel = trim(this%filtyp)//' NO.'
205  if (this%dis%ndim == 3) then
206  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
207  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
208  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
209  elseif (this%dis%ndim == 2) then
210  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
211  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
212  else
213  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
214  end if
215  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
216  if (this%inamedbound == 1) then
217  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
218  end if
219  !
220  ! -- Return
221  return

◆ src_allocate_scalars()

subroutine gwtsrcmodule::src_allocate_scalars ( class(gwtsrctype this)

Allocate scalars specific to this source loading package

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

110  ! -- dummy
111  class(GwtSrcType) :: this
112  !
113  ! -- call standard BndType allocate scalars
114  call this%BndType%allocate_scalars()
115  !
116  ! -- allocate the object and assign values to object variables
117  !
118  ! -- Set values
119  !
120  ! -- Return
121  return

◆ src_cf()

subroutine gwtsrcmodule::src_cf ( class(gwtsrctype this)

This subroutine:

  • calculates hcof and rhs terms
  • skip if no sources

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

131  ! -- dummy
132  class(GwtSrcType) :: this
133  ! -- local
134  integer(I4B) :: i, node
135  real(DP) :: q
136  !
137  ! -- Return if no sources
138  if (this%nbound == 0) return
139  !
140  ! -- Calculate hcof and rhs for each source entry
141  do i = 1, this%nbound
142  node = this%nodelist(i)
143  this%hcof(i) = dzero
144  if (this%ibound(node) <= 0) then
145  this%rhs(i) = dzero
146  cycle
147  end if
148  q = this%bound(1, i)
149  this%rhs(i) = -q
150  end do
151  !
152  ! -- Return
153  return

◆ 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
82  !
83  ! -- return
84  return
Here is the caller graph for this function:

◆ src_da()

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

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

90  ! -- modules
92  ! -- dummy
93  class(GwtSrcType) :: this
94  !
95  ! -- Deallocate parent package
96  call this%BndType%bnd_da()
97  !
98  ! -- scalars
99  !
100  ! -- Return
101  return

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

249  implicit none
250  ! -- dummy
251  class(GwtSrcType) :: this
252  ! -- local
253  integer(I4B) :: indx
254  !
255  call this%obs%StoreObsType('src', .true., indx)
256  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
257  !
258  ! -- Store obs type and assign procedure pointer
259  ! for to-mvr observation type.
260  call this%obs%StoreObsType('to-mvr', .true., indx)
261  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
262  !
263  ! -- Return
264  return
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 160 of file gwt-src.f90.

161  ! -- dummy
162  class(GwtSrcType) :: this
163  real(DP), dimension(:), intent(inout) :: rhs
164  integer(I4B), dimension(:), intent(in) :: ia
165  integer(I4B), dimension(:), intent(in) :: idxglo
166  class(MatrixBaseType), pointer :: matrix_sln
167  ! -- local
168  integer(I4B) :: i, n, ipos
169  !
170  ! -- pakmvrobj fc
171  if (this%imover == 1) then
172  call this%pakmvrobj%fc()
173  end if
174  !
175  ! -- Copy package rhs and hcof into solution rhs and amat
176  do i = 1, this%nbound
177  n = this%nodelist(i)
178  rhs(n) = rhs(n) + this%rhs(i)
179  ipos = ia(n)
180  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
181  !
182  ! -- If mover is active and mass is being withdrawn,
183  ! store available mass (as positive value).
184  if (this%imover == 1 .and. this%rhs(i) > dzero) then
185  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
186  end if
187  end do
188  !
189  ! -- Return
190  return

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

232  implicit none
233  ! -- dummy
234  class(GwtSrcType) :: this
235  !
236  src_obs_supported = .true.
237  !
238  ! -- Return
239  return

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

274  ! -- dummy
275  class(GwtSrcType), intent(inout) :: this
276  ! -- local
277  integer(I4B) :: i, nlinks
278  type(TimeSeriesLinkType), pointer :: tslink => null()
279  !
280  nlinks = this%TsManager%boundtslinks%Count()
281  do i = 1, nlinks
282  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
283  if (associated(tslink)) then
284  if (tslink%JCol == 1) then
285  tslink%Text = 'SMASSRATE'
286  end if
287  end if
288  end do
289  !
290  ! -- Return
291  return
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'