MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
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  end subroutine src_create
83 
84  !> @brief Deallocate memory
85  !<
86  subroutine src_da(this)
87  ! -- modules
89  ! -- dummy
90  class(gwtsrctype) :: this
91  !
92  ! -- Deallocate parent package
93  call this%BndType%bnd_da()
94  !
95  ! -- scalars
96  end subroutine src_da
97 
98  !> @brief Allocate scalars
99  !!
100  !! Allocate scalars specific to this source loading package
101  !<
102  subroutine src_allocate_scalars(this)
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
113  end subroutine src_allocate_scalars
114 
115  !> @brief Formulate the HCOF and RHS terms
116  !!
117  !! This subroutine:
118  !! - calculates hcof and rhs terms
119  !! - skip if no sources
120  !<
121  subroutine src_cf(this)
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
142  end subroutine src_cf
143 
144  !> @brief Add matrix terms related to specified mass source loading
145  !!
146  !! Copy rhs and hcof into solution rhs and amat
147  !<
148  subroutine src_fc(this, rhs, ia, idxglo, matrix_sln)
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
176  end subroutine src_fc
177 
178  !> @brief Define list labels
179  !!
180  !! Define the list heading that is written to iout when PRINT_INPUT
181  !! option is used.
182  !<
183  subroutine define_listlabel(this)
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
204  end subroutine define_listlabel
205 
206  ! -- Procedures related to observations
207  !> @brief Support function for specified mass source loading observations
208  !!
209  !! This function:
210  !! - returns true because SRC package supports observations.
211  !! - overrides BndType%bnd_obs_supported()
212  !<
213  logical function src_obs_supported(this)
214  implicit none
215  ! -- dummy
216  class(gwtsrctype) :: this
217  !
218  src_obs_supported = .true.
219  end function src_obs_supported
220 
221  !> @brief Define observations
222  !!
223  !! This subroutine:
224  !! - stores observation types supported by SRC package.
225  !! - overrides BndType%bnd_df_obs
226  !<
227  subroutine src_df_obs(this)
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
241  end subroutine src_df_obs
242 
243  !> @brief Procedure related to time series
244  !!
245  !! Assign tsLink%Text appropriately for all time series in use by package.
246  !! In the SRC package only the SMASSRATE variable can be controlled by time
247  !! series.
248  !<
249  subroutine src_rp_ts(this)
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
265  end subroutine src_rp_ts
266 
267 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:103
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:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
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:149
subroutine src_cf(this)
Formulate the HCOF and RHS terms.
Definition: gwt-src.f90:122
character(len=16) text
Definition: gwt-src.f90:18
subroutine src_allocate_scalars(this)
Allocate scalars.
Definition: gwt-src.f90:103
subroutine src_da(this)
Deallocate memory.
Definition: gwt-src.f90:87
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:184
subroutine src_rp_ts(this)
Procedure related to time series.
Definition: gwt-src.f90:250
subroutine src_df_obs(this)
Define observations.
Definition: gwt-src.f90:228
logical function src_obs_supported(this)
Support function for specified mass source loading observations.
Definition: gwt-src.f90:214
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:246
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType