MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwt-src.f90
Go to the documentation of this file.
2  !
3  use kindmodule, only: dp, i4b
5  use bndmodule, only: bndtype
11  !
12  implicit none
13  !
14  private
15  public :: src_create
16  !
17  character(len=LENFTYPE) :: ftype = 'SRC'
18  character(len=16) :: text = ' SRC'
19  !
20  type, extends(bndtype) :: gwtsrctype
21 
22  character(len=LENVARNAME) :: depvartype = '' !< stores string of dependent variable type, depending on model type
23 
24  contains
25 
26  procedure :: allocate_scalars => src_allocate_scalars
27  procedure :: bnd_cf => src_cf
28  procedure :: bnd_fc => src_fc
29  procedure :: bnd_da => src_da
30  procedure :: define_listlabel
31  ! -- methods for observations
32  procedure, public :: bnd_obs_supported => src_obs_supported
33  procedure, public :: bnd_df_obs => src_df_obs
34  ! -- methods for time series
35  procedure, public :: bnd_rp_ts => src_rp_ts
36 
37  end type gwtsrctype
38 
39 contains
40 
41  !> @brief Create a source loading package
42  !!
43  !! This subroutine points bndobj to the newly created package
44  !<
45  subroutine src_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
46  depvartype)
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
85  end subroutine src_create
86 
87  !> @brief Deallocate memory
88  !<
89  subroutine src_da(this)
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
102  end subroutine src_da
103 
104  !> @brief Allocate scalars
105  !!
106  !! Allocate scalars specific to this source loading package
107  !<
108  subroutine src_allocate_scalars(this)
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
122  end subroutine src_allocate_scalars
123 
124  !> @brief Formulate the HCOF and RHS terms
125  !!
126  !! This subroutine:
127  !! - calculates hcof and rhs terms
128  !! - skip if no sources
129  !<
130  subroutine src_cf(this)
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
154  end subroutine src_cf
155 
156  !> @brief Add matrix terms related to specified mass source loading
157  !!
158  !! Copy rhs and hcof into solution rhs and amat
159  !<
160  subroutine src_fc(this, rhs, ia, idxglo, matrix_sln)
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
191  end subroutine src_fc
192 
193  !> @brief Define list labels
194  !!
195  !! Define the list heading that is written to iout when PRINT_INPUT
196  !! option is used.
197  !<
198  subroutine define_listlabel(this)
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
222  end subroutine define_listlabel
223 
224  ! -- Procedures related to observations
225  !> @brief Support function for specified mass source loading observations
226  !!
227  !! This function:
228  !! - returns true because SRC package supports observations.
229  !! - overrides BndType%bnd_obs_supported()
230  !<
231  logical function src_obs_supported(this)
232  implicit none
233  ! -- dummy
234  class(gwtsrctype) :: this
235  !
236  src_obs_supported = .true.
237  !
238  ! -- Return
239  return
240  end function src_obs_supported
241 
242  !> @brief Define observations
243  !!
244  !! This subroutine:
245  !! - stores observation types supported by SRC package.
246  !! - overrides BndType%bnd_df_obs
247  !<
248  subroutine src_df_obs(this)
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
265  end subroutine src_df_obs
266 
267  !> @brief Procedure related to time series
268  !!
269  !! Assign tsLink%Text appropriately for all time series in use by package.
270  !! In the SRC package only the SMASSRATE variable can be controlled by time
271  !! series.
272  !<
273  subroutine src_rp_ts(this)
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
292  end subroutine src_rp_ts
293 
294 end module gwtsrcmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
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:102
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:38
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
character(len=lenftype) ftype
Definition: gwt-src.f90:17
subroutine src_fc(this, rhs, ia, idxglo, matrix_sln)
Add matrix terms related to specified mass source loading.
Definition: gwt-src.f90:161
subroutine src_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwt-src.f90:131
character(len=16) text
Definition: gwt-src.f90:18
subroutine src_allocate_scalars(this)
Allocate scalars.
Definition: gwt-src.f90:109
subroutine src_da(this)
Deallocate memory.
Definition: gwt-src.f90:90
subroutine, public src_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype)
Create a source loading package.
Definition: gwt-src.f90:47
subroutine define_listlabel(this)
Define list labels.
Definition: gwt-src.f90:199
subroutine src_rp_ts(this)
Procedure related to time series.
Definition: gwt-src.f90:274
subroutine src_df_obs(this)
Define observations.
Definition: gwt-src.f90:249
logical function src_obs_supported(this)
Support function for specified mass source loading observations.
Definition: gwt-src.f90:232
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:249
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType