MODFLOW 6  version 6.6.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_rp (this)
 @ brief SWF read and prepare 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 308 of file swf-flw.f90.

309  ! -- dummy variables
310  class(SwfFlwType), intent(inout) :: this !< SwfFlwType object
311  !
312  ! -- create the header list label
313  this%listlabel = trim(this%filtyp)//' NO.'
314  if (this%dis%ndim == 3) then
315  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
316  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
317  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
318  elseif (this%dis%ndim == 2) then
319  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
320  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
321  else
322  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
323  end if
324  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
325  if (this%inamedbound == 1) then
326  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
327  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 127 of file swf-flw.f90.

128  ! -- modules
130  ! -- dummy
131  class(SwfFlwType) :: this
132  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
133  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
134  ! -- local
135  !
136  ! -- call BndExtType allocate scalars
137  call this%BndExtType%allocate_arrays(nodelist, auxvar)
138  !
139  ! -- set q array input context pointer
140  call mem_setptr(this%q, 'Q', this%input_mempath)
141  !
142  ! -- checkin q array input context pointer
143  call mem_checkin(this%q, 'Q', this%memoryPath, &
144  '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 106 of file swf-flw.f90.

107  ! -- modules
109  ! -- dummy variables
110  class(SwfFlwType) :: this !< SwfFlwType object
111  !
112  ! -- call base type allocate scalars
113  call this%BndExtType%allocate_scalars()
114  !
115  ! -- allocate the object and assign values to object variables
116  ! none for this package
117  !
118  ! -- Set values
119  ! 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 374 of file swf-flw.f90.

375  ! -- dummy variables
376  class(SwfFlwType) :: this !< SwfFlwType object
377  ! -- local variables
378  integer(I4B) :: i
379  integer(I4B) :: n
380  integer(I4B) :: jj
381  real(DP) :: v
382  type(ObserveType), pointer :: obsrv => null()
383  !
384  ! -- clear the observations
385  call this%obs%obs_bd_clear()
386  !
387  ! -- Save simulated values for all of package's observations.
388  do i = 1, this%obs%npakobs
389  obsrv => this%obs%pakobs(i)%obsrv
390  if (obsrv%BndFound) then
391  do n = 1, obsrv%indxbnds_count
392  v = dnodata
393  jj = obsrv%indxbnds(n)
394  select case (obsrv%ObsTypeId)
395  case ('TO-MVR')
396  if (this%imover == 1) then
397  v = this%pakmvrobj%get_qtomvr(jj)
398  if (v > dzero) then
399  v = -v
400  end if
401  end if
402  case ('FLW')
403  v = this%simvals(jj)
404  case default
405  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
406  call store_error(errmsg)
407  end select
408  call this%obs%SaveOneSimval(obsrv, v)
409  end do
410  else
411  call this%obs%SaveOneSimval(obsrv, dnodata)
412  end if
413  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.

Parameters
[in,out]thisBndExtType object

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

466  ! -- modules
467  use constantsmodule, only: dzero
468  ! -- dummy variables
469  class(SwfFlwType), intent(inout) :: this !< BndExtType object
470  integer(I4B), intent(in) :: col
471  integer(I4B), intent(in) :: row
472  ! -- result
473  real(DP) :: bndval
474  !
475  select case (col)
476  case (1)
477  bndval = this%q_mult(row)
478  case default
479  errmsg = 'Programming error. FLW bound value requested column '&
480  &'outside range of ncolbnd (1).'
481  call store_error(errmsg)
482  call store_error_filename(this%input_fname)
483  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 241 of file swf-flw.f90.

242  ! -- dummy variables
243  class(SwfFlwType) :: this !< SwfFlwType object
244  ! -- local variables
245  integer(I4B) :: i, node
246  real(DP) :: q
247  !
248  ! -- Return if no inflows
249  if (this%nbound == 0) return
250  !
251  ! -- Calculate hcof and rhs for each flw entry
252  do i = 1, this%nbound
253  node = this%nodelist(i)
254  this%hcof(i) = dzero
255  if (this%ibound(node) <= 0) then
256  this%rhs(i) = dzero
257  cycle
258  end if
259  q = this%q_mult(i)
260  this%rhs(i) = -q
261  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 63 of file swf-flw.f90.

65  ! -- dummy variables
66  class(BndType), pointer :: packobj !< pointer to default package type
67  integer(I4B), intent(in) :: id !< package id
68  integer(I4B), intent(in) :: ibcnum !< boundary condition number
69  integer(I4B), intent(in) :: inunit !< unit number of FLW package input file
70  integer(I4B), intent(in) :: iout !< unit number of model listing file
71  character(len=*), intent(in) :: namemodel !< model name
72  character(len=*), intent(in) :: pakname !< package name
73  character(len=*), intent(in) :: mempath !< input mempath
74  ! -- local variables
75  type(SwfFlwType), pointer :: flwobj
76  !
77  ! -- allocate the object and assign values to object variables
78  allocate (flwobj)
79  packobj => flwobj
80  !
81  ! -- create name and memory path
82  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
83  packobj%text = text
84  !
85  ! -- allocate scalars
86  call flwobj%allocate_scalars()
87  !
88  ! -- initialize package
89  call packobj%pack_initialize()
90 
91  packobj%inunit = inunit
92  packobj%iout = iout
93  packobj%id = id
94  packobj%ibcnum = ibcnum
95  packobj%ncolbnd = 1
96  packobj%iscloc = 1
97  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 152 of file swf-flw.f90.

153  ! -- modules
155  ! -- dummy variables
156  class(SwfFlwType) :: this !< SwfFlwType object
157  !
158  ! -- Deallocate parent package
159  call this%BndExtType%bnd_da()
160  !
161  ! -- scalars
162  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 353 of file swf-flw.f90.

354  ! -- dummy variables
355  class(SwfFlwType) :: this !< SwfFlwType object
356  ! -- local variables
357  integer(I4B) :: indx
358  !
359  ! -- initialize observations
360  call this%obs%StoreObsType('flw', .true., indx)
361  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
362  !
363  ! -- Store obs type and assign procedure pointer
364  ! for to-mvr observation type.
365  call this%obs%StoreObsType('to-mvr', .true., indx)
366  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 270 of file swf-flw.f90.

271  ! -- dummy variables
272  class(SwfFlwType) :: this !< SwfFlwType object
273  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
274  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
275  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
276  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
277  ! -- local variables
278  integer(I4B) :: i
279  integer(I4B) :: n
280  integer(I4B) :: ipos
281  !
282  ! -- pakmvrobj fc
283  if (this%imover == 1) then
284  call this%pakmvrobj%fc()
285  end if
286  !
287  ! -- Copy package rhs and hcof into solution rhs and amat
288  do i = 1, this%nbound
289  n = this%nodelist(i)
290  rhs(n) = rhs(n) + this%rhs(i)
291  ipos = ia(n)
292  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
293  !
294  ! -- If mover is active and this flw item is discharging,
295  ! store available water (as positive value).
296  if (this%imover == 1 .and. this%rhs(i) > dzero) then
297  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
298  end if
299  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 340 of file swf-flw.f90.

341  ! -- dummy variables
342  class(SwfFlwType) :: this !< SwfFlwType object
343  !
344  ! -- set boolean
345  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 170 of file swf-flw.f90.

171  ! -- modules
172  use inputoutputmodule, only: urword
175  ! -- dummy variables
176  class(SwfFlwType), intent(inout) :: this !< SwfFlwType object
177  ! -- local variables
178  type(SwfFlwParamFoundType) :: found
179  ! -- formats
180  !
181  ! -- source base BndExtType options
182  call this%BndExtType%source_options()
183  !
184  ! -- source options from input context
185  ! none
186  !
187  ! -- log SWF specific options
188  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()

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

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

219  use tdismodule, only: kper
220  ! -- dummy
221  class(SwfFlwType), intent(inout) :: this
222  ! -- local
223  !
224  if (this%iper /= kper) return
225  !
226  ! -- Call the parent class read and prepare
227  call this%BndExtType%bnd_rp()
228  !
229  ! -- Write the list to iout if requested
230  if (this%iprpak /= 0) then
231  call this%write_list()
232  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ 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 424 of file swf-flw.f90.

425  ! -- dummy variables
426  class(SwfFlwType), intent(inout) :: this !< SwfFlwType object
427  ! -- local variables
428  integer(I4B) :: i, nlinks
429  type(TimeSeriesLinkType), pointer :: tslink => null()
430  !
431  ! -- set up the time series links
432  nlinks = this%TsManager%boundtslinks%Count()
433  do i = 1, nlinks
434  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
435  if (associated(tslink)) then
436  if (tslink%JCol == 1) then
437  tslink%Text = 'Q'
438  end if
439  end if
440  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 
)
Parameters
[in,out]thisBndExtType object

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

194  ! -- modules
196  ! -- dummy variables
197  class(SwfFlwType), intent(inout) :: this !< BndExtType object
198  type(SwfFlwParamFoundType), intent(in) :: found
199  ! -- local variables
200  ! -- format
201  !
202  ! -- log found options
203  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
204  //' OPTIONS'
205  !
206  ! if (found%mover) then
207  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
208  ! end if
209  !
210  ! -- close logging block
211  write (this%iout, '(1x,a)') &
212  '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
Parameters
[in,out]thisBndExtType object

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

444  ! -- modules
445  use constantsmodule, only: dzero
446  ! -- dummy variables
447  class(SwfFlwType), intent(inout) :: this !< BndExtType object
448  integer(I4B), intent(in) :: row
449  ! -- result
450  real(DP) :: q
451  !
452  if (this%iauxmultcol > 0) then
453  q = this%q(row) * this%auxvar(this%iauxmultcol, row)
454  else
455  q = this%q(row)
456  end if

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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