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