MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
swf-flw.f90
Go to the documentation of this file.
1 !> @brief This module contains the FLW package methods
2 !!
3 !! This module can be used to represent inflow to streams. It is
4 !! designed similarly to the GWF WEL package.
5 !!
6 !<
8  ! -- modules
9  use kindmodule, only: dp, i4b
11  use simvariablesmodule, only: errmsg
13  use bndmodule, only: bndtype
14  use bndextmodule, only: bndexttype
16  use observemodule, only: observetype
20  !
21  implicit none
22  !
23  private
24  public :: flw_create
25  !
26  character(len=LENFTYPE) :: ftype = 'FLW' !< package ftype
27  character(len=16) :: text = ' FLW' !< package flow text string
28  !
29  type, extends(bndexttype) :: swfflwtype
30  real(dp), dimension(:), pointer, contiguous :: q => null() !< volumetric rate
31  contains
32  procedure :: allocate_scalars => flw_allocate_scalars
33  procedure :: allocate_arrays => flw_allocate_arrays
34  procedure :: source_options => flw_options
35  procedure :: log_flw_options
36  procedure :: bnd_cf => flw_cf
37  procedure :: bnd_fc => flw_fc
38  procedure :: bnd_da => flw_da
39  procedure :: define_listlabel
40  procedure :: bound_value => flw_bound_value
41  procedure :: q_mult
42  ! -- methods for observations
43  procedure, public :: bnd_obs_supported => flw_obs_supported
44  procedure, public :: bnd_df_obs => flw_df_obs
45  procedure, public :: bnd_bd_obs => flw_bd_obs
46  ! -- methods for time series
47  procedure, public :: bnd_rp_ts => flw_rp_ts
48  end type swfflwtype
49 
50 contains
51 
52  !> @ brief Create a new package object
53  !!
54  !! Create a new FLW Package object
55  !!
56  !<
57  subroutine flw_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
58  mempath)
59  ! -- dummy variables
60  class(bndtype), pointer :: packobj !< pointer to default package type
61  integer(I4B), intent(in) :: id !< package id
62  integer(I4B), intent(in) :: ibcnum !< boundary condition number
63  integer(I4B), intent(in) :: inunit !< unit number of FLW package input file
64  integer(I4B), intent(in) :: iout !< unit number of model listing file
65  character(len=*), intent(in) :: namemodel !< model name
66  character(len=*), intent(in) :: pakname !< package name
67  character(len=*), intent(in) :: mempath !< input mempath
68  ! -- local variables
69  type(swfflwtype), pointer :: flwobj
70  !
71  ! -- allocate the object and assign values to object variables
72  allocate (flwobj)
73  packobj => flwobj
74  !
75  ! -- create name and memory path
76  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
77  packobj%text = text
78  !
79  ! -- allocate scalars
80  call flwobj%allocate_scalars()
81  !
82  ! -- initialize package
83  call packobj%pack_initialize()
84 
85  packobj%inunit = inunit
86  packobj%iout = iout
87  packobj%id = id
88  packobj%ibcnum = ibcnum
89  packobj%ncolbnd = 1
90  packobj%iscloc = 1
91  packobj%ictMemPath = ''
92  end subroutine flw_create
93 
94  !> @ brief Allocate scalars
95  !!
96  !! Allocate and initialize scalars for the FLW package. The base model
97  !! allocate scalars method is also called.
98  !!
99  !<
100  subroutine flw_allocate_scalars(this)
101  ! -- modules
103  ! -- dummy variables
104  class(swfflwtype) :: this !< SwfFlwType object
105  !
106  ! -- call base type allocate scalars
107  call this%BndExtType%allocate_scalars()
108  !
109  ! -- allocate the object and assign values to object variables
110  ! none for this package
111  !
112  ! -- Set values
113  ! none for this package
114  end subroutine flw_allocate_scalars
115 
116  !> @ brief Allocate arrays
117  !!
118  !! Allocate and initialize arrays for the SWF package
119  !!
120  !<
121  subroutine flw_allocate_arrays(this, nodelist, auxvar)
122  ! -- modules
124  ! -- dummy
125  class(swfflwtype) :: this
126  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
127  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
128  ! -- local
129  !
130  ! -- call BndExtType allocate scalars
131  call this%BndExtType%allocate_arrays(nodelist, auxvar)
132  !
133  ! -- set q array input context pointer
134  call mem_setptr(this%q, 'Q', this%input_mempath)
135  !
136  ! -- checkin q array input context pointer
137  call mem_checkin(this%q, 'Q', this%memoryPath, &
138  'Q', this%input_mempath)
139  end subroutine flw_allocate_arrays
140 
141  !> @ brief Deallocate package memory
142  !!
143  !! Deallocate SWF package scalars and arrays.
144  !!
145  !<
146  subroutine flw_da(this)
147  ! -- modules
149  ! -- dummy variables
150  class(swfflwtype) :: this !< SwfFlwType object
151  !
152  ! -- Deallocate parent package
153  call this%BndExtType%bnd_da()
154  !
155  ! -- scalars
156  call mem_deallocate(this%q, 'Q', this%memoryPath)
157  end subroutine flw_da
158 
159  !> @ brief Source additional options for package
160  !!
161  !! Source additional options for SWF package.
162  !!
163  !<
164  subroutine flw_options(this)
165  ! -- modules
166  use inputoutputmodule, only: urword
169  ! -- dummy variables
170  class(swfflwtype), intent(inout) :: this !< SwfFlwType object
171  ! -- local variables
172  type(swfflwparamfoundtype) :: found
173  ! -- formats
174  !
175  ! -- source base BndExtType options
176  call this%BndExtType%source_options()
177  !
178  ! -- source options from input context
179  ! none
180  !
181  ! -- log SWF specific options
182  call this%log_flw_options(found)
183  end subroutine flw_options
184 
185  !> @ brief Log SWF specific package options
186  !<
187  subroutine log_flw_options(this, found)
188  ! -- modules
190  ! -- dummy variables
191  class(swfflwtype), intent(inout) :: this
192  type(swfflwparamfoundtype), intent(in) :: found
193  ! -- local variables
194  ! -- format
195  !
196  ! -- log found options
197  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
198  //' OPTIONS'
199  !
200  ! if (found%mover) then
201  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
202  ! end if
203  !
204  ! -- close logging block
205  write (this%iout, '(1x,a)') &
206  'END OF '//trim(adjustl(this%text))//' OPTIONS'
207  end subroutine log_flw_options
208 
209  !> @ brief Formulate the package hcof and rhs terms.
210  !!
211  !! Formulate the hcof and rhs terms for the FLW package that will be
212  !! added to the coefficient matrix and right-hand side vector.
213  !!
214  !<
215  subroutine flw_cf(this)
216  ! -- dummy variables
217  class(swfflwtype) :: this !< SwfFlwType object
218  ! -- local variables
219  integer(I4B) :: i, node
220  real(DP) :: q
221  !
222  ! -- Return if no inflows
223  if (this%nbound == 0) return
224  !
225  ! -- Calculate hcof and rhs for each flw entry
226  do i = 1, this%nbound
227  node = this%nodelist(i)
228  this%hcof(i) = dzero
229  if (this%ibound(node) <= 0) then
230  this%rhs(i) = dzero
231  cycle
232  end if
233  q = this%q_mult(i)
234  this%rhs(i) = -q
235  end do
236  end subroutine flw_cf
237 
238  !> @ brief Copy hcof and rhs terms into solution.
239  !!
240  !! Add the hcof and rhs terms for the FLW package to the
241  !! coefficient matrix and right-hand side vector.
242  !!
243  !<
244  subroutine flw_fc(this, rhs, ia, idxglo, matrix_sln)
245  ! -- dummy variables
246  class(swfflwtype) :: this !< SwfFlwType object
247  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
248  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
249  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
250  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
251  ! -- local variables
252  integer(I4B) :: i
253  integer(I4B) :: n
254  integer(I4B) :: ipos
255  !
256  ! -- pakmvrobj fc
257  if (this%imover == 1) then
258  call this%pakmvrobj%fc()
259  end if
260  !
261  ! -- Copy package rhs and hcof into solution rhs and amat
262  do i = 1, this%nbound
263  n = this%nodelist(i)
264  rhs(n) = rhs(n) + this%rhs(i)
265  ipos = ia(n)
266  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
267  !
268  ! -- If mover is active and this flw item is discharging,
269  ! store available water (as positive value).
270  if (this%imover == 1 .and. this%rhs(i) > dzero) then
271  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
272  end if
273  end do
274  end subroutine flw_fc
275 
276  !> @ brief Define the list label for the package
277  !!
278  !! Method defined the list label for the FLW package. The list label is
279  !! the heading that is written to iout when PRINT_INPUT option is used.
280  !!
281  !<
282  subroutine define_listlabel(this)
283  ! -- dummy variables
284  class(swfflwtype), intent(inout) :: this !< SwfFlwType object
285  !
286  ! -- create the header list label
287  this%listlabel = trim(this%filtyp)//' NO.'
288  if (this%dis%ndim == 3) then
289  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
290  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
291  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
292  elseif (this%dis%ndim == 2) then
293  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
294  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
295  else
296  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
297  end if
298  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
299  if (this%inamedbound == 1) then
300  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
301  end if
302  end subroutine define_listlabel
303 
304  ! -- Procedures related to observations
305 
306  !> @brief Determine if observations are supported.
307  !!
308  !! Function to determine if observations are supported by the FLW package.
309  !! Observations are supported by the FLW package.
310  !!
311  !! @return flw_obs_supported boolean indicating if observations are supported
312  !!
313  !<
314  logical function flw_obs_supported(this)
315  ! -- dummy variables
316  class(swfflwtype) :: this !< SwfFlwType object
317  !
318  ! -- set boolean
319  flw_obs_supported = .true.
320  end function flw_obs_supported
321 
322  !> @brief Define the observation types available in the package
323  !!
324  !! Method to define the observation types available in the FLW package.
325  !!
326  !<
327  subroutine flw_df_obs(this)
328  ! -- dummy variables
329  class(swfflwtype) :: this !< SwfFlwType object
330  ! -- local variables
331  integer(I4B) :: indx
332  !
333  ! -- initialize observations
334  call this%obs%StoreObsType('flw', .true., indx)
335  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
336  !
337  ! -- Store obs type and assign procedure pointer
338  ! for to-mvr observation type.
339  call this%obs%StoreObsType('to-mvr', .true., indx)
340  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
341  end subroutine flw_df_obs
342 
343  !> @brief Save observations for the package
344  !!
345  !! Method to save simulated values for the FLW package.
346  !!
347  !<
348  subroutine flw_bd_obs(this)
349  ! -- dummy variables
350  class(swfflwtype) :: this !< SwfFlwType object
351  ! -- local variables
352  integer(I4B) :: i
353  integer(I4B) :: n
354  integer(I4B) :: jj
355  real(DP) :: v
356  type(observetype), pointer :: obsrv => null()
357  !
358  ! -- clear the observations
359  call this%obs%obs_bd_clear()
360  !
361  ! -- Save simulated values for all of package's observations.
362  do i = 1, this%obs%npakobs
363  obsrv => this%obs%pakobs(i)%obsrv
364  if (obsrv%BndFound) then
365  do n = 1, obsrv%indxbnds_count
366  v = dnodata
367  jj = obsrv%indxbnds(n)
368  select case (obsrv%ObsTypeId)
369  case ('TO-MVR')
370  if (this%imover == 1) then
371  v = this%pakmvrobj%get_qtomvr(jj)
372  if (v > dzero) then
373  v = -v
374  end if
375  end if
376  case ('FLW')
377  v = this%simvals(jj)
378  case default
379  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
380  call store_error(errmsg)
381  end select
382  call this%obs%SaveOneSimval(obsrv, v)
383  end do
384  else
385  call this%obs%SaveOneSimval(obsrv, dnodata)
386  end if
387  end do
388  end subroutine flw_bd_obs
389 
390  ! -- Procedure related to time series
391 
392  !> @brief Assign time series links for the package
393  !!
394  !! Assign the time series links for the FLW package. Only
395  !! the Q variable can be defined with time series.
396  !!
397  !<
398  subroutine flw_rp_ts(this)
399  ! -- dummy variables
400  class(swfflwtype), intent(inout) :: this !< SwfFlwType object
401  ! -- local variables
402  integer(I4B) :: i, nlinks
403  type(timeserieslinktype), pointer :: tslink => null()
404  !
405  ! -- set up the time series links
406  nlinks = this%TsManager%boundtslinks%Count()
407  do i = 1, nlinks
408  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
409  if (associated(tslink)) then
410  if (tslink%JCol == 1) then
411  tslink%Text = 'Q'
412  end if
413  end if
414  end do
415  end subroutine flw_rp_ts
416 
417  function q_mult(this, row) result(q)
418  ! -- modules
419  use constantsmodule, only: dzero
420  ! -- dummy variables
421  class(swfflwtype), intent(inout) :: this
422  integer(I4B), intent(in) :: row
423  ! -- result
424  real(dp) :: q
425  !
426  if (this%iauxmultcol > 0) then
427  q = this%q(row) * this%auxvar(this%iauxmultcol, row)
428  else
429  q = this%q(row)
430  end if
431  end function q_mult
432 
433  !> @ brief Return a bound value
434  !!
435  !! Return a bound value associated with an ncolbnd index
436  !! and row.
437  !!
438  !<
439  function flw_bound_value(this, col, row) result(bndval)
440  ! -- modules
441  use constantsmodule, only: dzero
442  ! -- dummy variables
443  class(swfflwtype), intent(inout) :: this
444  integer(I4B), intent(in) :: col
445  integer(I4B), intent(in) :: row
446  ! -- result
447  real(dp) :: bndval
448  !
449  select case (col)
450  case (1)
451  bndval = this%q_mult(row)
452  case default
453  errmsg = 'Programming error. FLW bound value requested column '&
454  &'outside range of ncolbnd (1).'
455  call store_error(errmsg)
456  call store_error_filename(this%input_fname)
457  end select
458  end function flw_bound_value
459 
460 end module swfflwmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
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
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module defines variable data types.
Definition: kind.f90:8
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
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
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
This module contains the FLW package methods.
Definition: swf-flw.f90:7
subroutine flw_df_obs(this)
Define the observation types available in the package.
Definition: swf-flw.f90:328
character(len=lenftype) ftype
package ftype
Definition: swf-flw.f90:26
subroutine flw_cf(this)
@ brief Formulate the package hcof and rhs terms.
Definition: swf-flw.f90:216
subroutine, public flw_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
Definition: swf-flw.f90:59
subroutine flw_da(this)
@ brief Deallocate package memory
Definition: swf-flw.f90:147
real(dp) function flw_bound_value(this, col, row)
@ brief Return a bound value
Definition: swf-flw.f90:440
logical function flw_obs_supported(this)
Determine if observations are supported.
Definition: swf-flw.f90:315
subroutine flw_allocate_scalars(this)
@ brief Allocate scalars
Definition: swf-flw.f90:101
character(len=16) text
package flow text string
Definition: swf-flw.f90:27
subroutine flw_bd_obs(this)
Save observations for the package.
Definition: swf-flw.f90:349
subroutine define_listlabel(this)
@ brief Define the list label for the package
Definition: swf-flw.f90:283
subroutine log_flw_options(this, found)
@ brief Log SWF specific package options
Definition: swf-flw.f90:188
subroutine flw_options(this)
@ brief Source additional options for package
Definition: swf-flw.f90:165
real(dp) function q_mult(this, row)
Definition: swf-flw.f90:418
subroutine flw_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
Definition: swf-flw.f90:122
subroutine flw_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
Definition: swf-flw.f90:245
subroutine flw_rp_ts(this)
Assign time series links for the package.
Definition: swf-flw.f90:399
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
@ brief BndType