MODFLOW 6  version 6.7.0.dev1
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 303 of file swf-flw.f90.

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

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

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)

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

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

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

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

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

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

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

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 = ''
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 147 of file swf-flw.f90.

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)

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

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

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

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

336  ! -- dummy variables
337  class(SwfFlwType) :: this !< SwfFlwType object
338  !
339  ! -- set boolean
340  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 165 of file swf-flw.f90.

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)
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 213 of file swf-flw.f90.

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

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

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'

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

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

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