MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
swfflwmodule Module Reference

This module contains the FLW package methods. More...

Data Types

type  swfflwtype
 

Functions/Subroutines

subroutine, public flw_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
 @ brief Create a new package object More...
 
subroutine flw_allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine flw_allocate_arrays (this, nodelist, auxvar)
 @ brief Allocate arrays More...
 
subroutine flw_da (this)
 @ brief Deallocate package memory More...
 
subroutine flw_options (this)
 @ brief Source additional options for package More...
 
subroutine log_flw_options (this, found)
 @ brief Log SWF specific package options More...
 
subroutine flw_cf (this)
 @ brief Formulate the package hcof and rhs terms. More...
 
subroutine flw_fc (this, rhs, ia, idxglo, matrix_sln)
 @ brief Copy hcof and rhs terms into solution. More...
 
subroutine define_listlabel (this)
 @ brief Define the list label for the package More...
 
logical function flw_obs_supported (this)
 Determine if observations are supported. More...
 
subroutine flw_df_obs (this)
 Define the observation types available in the package. More...
 
subroutine flw_bd_obs (this)
 Save observations for the package. More...
 
subroutine flw_rp_ts (this)
 Assign time series links for the package. More...
 
real(dp) function q_mult (this, row)
 
real(dp) function flw_bound_value (this, col, row)
 @ brief Return a bound value More...
 

Variables

character(len=lenftype) ftype = 'FLW'
 package ftype More...
 
character(len=16) text = ' FLW'
 package flow text string More...
 

Detailed Description

This module can be used to represent inflow to streams. It is designed similarly to the GWF WEL package.

Function/Subroutine Documentation

◆ define_listlabel()

subroutine swfflwmodule::define_listlabel ( class(swfflwtype), intent(inout)  this)
private

Method defined the list label for the FLW package. The list label is the heading that is written to iout when PRINT_INPUT option is used.

Parameters
[in,out]thisSwfFlwType object

Definition at line 282 of file swf-flw.f90.

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

◆ flw_allocate_arrays()

subroutine swfflwmodule::flw_allocate_arrays ( class(swfflwtype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)

Allocate and initialize arrays for the SWF package

Definition at line 121 of file swf-flw.f90.

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)

◆ flw_allocate_scalars()

subroutine swfflwmodule::flw_allocate_scalars ( class(swfflwtype this)
private

Allocate and initialize scalars for the FLW package. The base model allocate scalars method is also called.

Parameters
thisSwfFlwType object

Definition at line 100 of file swf-flw.f90.

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

◆ flw_bd_obs()

subroutine swfflwmodule::flw_bd_obs ( class(swfflwtype this)
private

Method to save simulated values for the FLW package.

Parameters
thisSwfFlwType object

Definition at line 348 of file swf-flw.f90.

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
Here is the call graph for this function:

◆ flw_bound_value()

real(dp) function swfflwmodule::flw_bound_value ( class(swfflwtype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)

Return a bound value associated with an ncolbnd index and row.

Definition at line 439 of file swf-flw.f90.

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
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
Here is the call graph for this function:

◆ flw_cf()

subroutine swfflwmodule::flw_cf ( class(swfflwtype this)

Formulate the hcof and rhs terms for the FLW package that will be added to the coefficient matrix and right-hand side vector.

Parameters
thisSwfFlwType object

Definition at line 215 of file swf-flw.f90.

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

◆ flw_create()

subroutine, public swfflwmodule::flw_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=*), intent(in)  mempath 
)

Create a new FLW Package object

Parameters
packobjpointer to default package type
[in]idpackage id
[in]ibcnumboundary condition number
[in]inunitunit number of FLW package input file
[in]ioutunit number of model listing file
[in]namemodelmodel name
[in]paknamepackage name
[in]mempathinput mempath

Definition at line 57 of file swf-flw.f90.

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 = ''
Here is the caller graph for this function:

◆ flw_da()

subroutine swfflwmodule::flw_da ( class(swfflwtype this)

Deallocate SWF package scalars and arrays.

Parameters
thisSwfFlwType object

Definition at line 146 of file swf-flw.f90.

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)

◆ flw_df_obs()

subroutine swfflwmodule::flw_df_obs ( class(swfflwtype this)
private

Method to define the observation types available in the FLW package.

Parameters
thisSwfFlwType object

Definition at line 327 of file swf-flw.f90.

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
Here is the call graph for this function:

◆ flw_fc()

subroutine swfflwmodule::flw_fc ( class(swfflwtype 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

Add the hcof and rhs terms for the FLW package to the coefficient matrix and right-hand side vector.

Parameters
thisSwfFlwType object
[in,out]rhsright-hand side vector for model
[in]iasolution CRS row pointers
[in]idxglomapping vector for model (local) to solution (global)
matrix_slnsolution coefficient matrix

Definition at line 244 of file swf-flw.f90.

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

◆ flw_obs_supported()

logical function swfflwmodule::flw_obs_supported ( class(swfflwtype this)
private

Function to determine if observations are supported by the FLW package. Observations are supported by the FLW package.

Returns
flw_obs_supported boolean indicating if observations are supported
Parameters
thisSwfFlwType object

Definition at line 314 of file swf-flw.f90.

315  ! -- dummy variables
316  class(SwfFlwType) :: this !< SwfFlwType object
317  !
318  ! -- set boolean
319  flw_obs_supported = .true.

◆ flw_options()

subroutine swfflwmodule::flw_options ( class(swfflwtype), intent(inout)  this)

Source additional options for SWF package.

Parameters
[in,out]thisSwfFlwType object

Definition at line 164 of file swf-flw.f90.

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)
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
Here is the call graph for this function:

◆ flw_rp_ts()

subroutine swfflwmodule::flw_rp_ts ( class(swfflwtype), intent(inout)  this)
private

Assign the time series links for the FLW package. Only the Q variable can be defined with time series.

Parameters
[in,out]thisSwfFlwType object

Definition at line 398 of file swf-flw.f90.

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
Here is the call graph for this function:

◆ log_flw_options()

subroutine swfflwmodule::log_flw_options ( class(swfflwtype), intent(inout)  this,
type(swfflwparamfoundtype), intent(in)  found 
)

Definition at line 187 of file swf-flw.f90.

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'

◆ q_mult()

real(dp) function swfflwmodule::q_mult ( class(swfflwtype), intent(inout)  this,
integer(i4b), intent(in)  row 
)
private

Definition at line 417 of file swf-flw.f90.

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

Variable Documentation

◆ ftype

character(len=lenftype) swfflwmodule::ftype = 'FLW'
private

Definition at line 26 of file swf-flw.f90.

26  character(len=LENFTYPE) :: ftype = 'FLW' !< package ftype

◆ text

character(len=16) swfflwmodule::text = ' FLW'
private

Definition at line 27 of file swf-flw.f90.

27  character(len=16) :: text = ' FLW' !< package flow text string