MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
swfcdbmodule Module Reference

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

Data Types

type  swfcdbtype
 

Functions/Subroutines

subroutine, public cdb_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, lengthconv, timeconv)
 @ brief Create a new package object More...
 
subroutine cdb_allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine cdb_allocate_arrays (this, nodelist, auxvar)
 @ brief Allocate arrays More...
 
subroutine cdb_da (this)
 @ brief Deallocate package memory More...
 
subroutine cdb_options (this)
 @ brief Source additional options for package More...
 
subroutine log_cdb_options (this, found)
 @ brief Log SWF specific package options More...
 
subroutine cdb_rp (this)
 @ brief SWF read and prepare More...
 
subroutine cdb_cf (this)
 @ brief Formulate the package hcof and rhs terms. More...
 
real(dp) function qcalc (this, i, depth)
 Calculate critical depth boundary flow. More...
 
subroutine cdb_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 cdb_obs_supported (this)
 Determine if observations are supported. More...
 
subroutine cdb_df_obs (this)
 Define the observation types available in the package. More...
 
subroutine cdb_bd_obs (this)
 Save observations for the package. More...
 
subroutine cdb_rp_ts (this)
 Assign time series links for the package. More...
 
real(dp) function cdb_bound_value (this, col, row)
 @ brief Return a bound value More...
 

Variables

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

Detailed Description

This module can be used to represent outflow from streams using a critical depth boundary.

Function/Subroutine Documentation

◆ cdb_allocate_arrays()

subroutine swfcdbmodule::cdb_allocate_arrays ( class(swfcdbtype 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 147 of file swf-cdb.f90.

148  ! -- modules
150  ! -- dummy
151  class(SwfCdbType) :: this
152  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
153  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
154  ! -- local
155  !
156  ! -- call BndExtType allocate scalars
157  call this%BndExtType%allocate_arrays(nodelist, auxvar)
158  !
159  ! -- set array input context pointer
160  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
161  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
162  !
163  ! -- checkin array input context pointer
164  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
165  'IDCXS', this%input_mempath)
166  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
167  'WIDTH', this%input_mempath)

◆ cdb_allocate_scalars()

subroutine swfcdbmodule::cdb_allocate_scalars ( class(swfcdbtype this)
private

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

Parameters
thisSwfcdbType object

Definition at line 126 of file swf-cdb.f90.

127  ! -- modules
129  ! -- dummy variables
130  class(SwfCdbType) :: this !< SwfcdbType object
131  !
132  ! -- call base type allocate scalars
133  call this%BndExtType%allocate_scalars()
134  !
135  ! -- allocate the object and assign values to object variables
136  call mem_allocate(this%gravconv, 'GRAVCONV', this%memoryPath)
137  !
138  ! -- Set values
139  this%gravconv = dzero

◆ cdb_bd_obs()

subroutine swfcdbmodule::cdb_bd_obs ( class(swfcdbtype this)
private

Method to save simulated values for the CDB package.

Parameters
thisSwfCdbType object

Definition at line 455 of file swf-cdb.f90.

456  ! -- dummy variables
457  class(SwfCdbType) :: this !< SwfCdbType object
458  ! -- local variables
459  integer(I4B) :: i
460  integer(I4B) :: n
461  integer(I4B) :: jj
462  real(DP) :: v
463  type(ObserveType), pointer :: obsrv => null()
464  !
465  ! -- clear the observations
466  call this%obs%obs_bd_clear()
467  !
468  ! -- Save simulated values for all of package's observations.
469  do i = 1, this%obs%npakobs
470  obsrv => this%obs%pakobs(i)%obsrv
471  if (obsrv%BndFound) then
472  do n = 1, obsrv%indxbnds_count
473  v = dnodata
474  jj = obsrv%indxbnds(n)
475  select case (obsrv%ObsTypeId)
476  case ('TO-MVR')
477  if (this%imover == 1) then
478  v = this%pakmvrobj%get_qtomvr(jj)
479  if (v > dzero) then
480  v = -v
481  end if
482  end if
483  case ('CDB')
484  v = this%simvals(jj)
485  case default
486  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
487  call store_error(errmsg)
488  end select
489  call this%obs%SaveOneSimval(obsrv, v)
490  end do
491  else
492  call this%obs%SaveOneSimval(obsrv, dnodata)
493  end if
494  end do
Here is the call graph for this function:

◆ cdb_bound_value()

real(dp) function swfcdbmodule::cdb_bound_value ( class(swfcdbtype), intent(inout)  this,
integer(i4b), intent(in)  col,
integer(i4b), intent(in)  row 
)
private

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

Parameters
[in,out]thisBndExtType object

Definition at line 530 of file swf-cdb.f90.

531  ! -- modules
532  use constantsmodule, only: dzero
533  ! -- dummy variables
534  class(SwfCdbType), intent(inout) :: this !< BndExtType object
535  integer(I4B), intent(in) :: col
536  integer(I4B), intent(in) :: row
537  ! -- result
538  real(DP) :: bndval
539  !
540  select case (col)
541  case (1)
542  bndval = this%idcxs(row)
543  case (2)
544  bndval = this%width(row)
545  case default
546  errmsg = 'Programming error. CDB bound value requested column '&
547  &'outside range of ncolbnd (1).'
548  call store_error(errmsg)
549  call store_error_filename(this%input_fname)
550  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:

◆ cdb_cf()

subroutine swfcdbmodule::cdb_cf ( class(swfcdbtype this)

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

Parameters
thisSwfCdbType object

Definition at line 268 of file swf-cdb.f90.

269  ! modules
271  ! -- dummy variables
272  class(SwfCdbType) :: this !< SwfCdbType object
273  ! -- local variables
274  integer(I4B) :: i, node
275  real(DP) :: q
276  real(DP) :: qeps
277  real(DP) :: depth
278  real(DP) :: derv
279  real(DP) :: eps
280  !
281  ! -- Return if no inflows
282  if (this%nbound == 0) return
283  !
284  ! -- Calculate hcof and rhs for each cdb entry
285  do i = 1, this%nbound
286 
287  node = this%nodelist(i)
288  if (this%ibound(node) <= 0) then
289  this%hcof(i) = dzero
290  this%rhs(i) = dzero
291  cycle
292  end if
293 
294  ! -- calculate terms and add to hcof and rhs
295  depth = this%xnew(node) - this%dis%bot(node)
296 
297  ! -- calculate unperturbed q
298  q = this%qcalc(i, depth)
299 
300  ! -- calculate perturbed q
301  eps = get_perturbation(depth)
302  qeps = this%qcalc(i, depth + eps)
303 
304  ! -- calculate derivative
305  derv = (qeps - q) / eps
306 
307  ! -- add terms to hcof and rhs
308  this%hcof(i) = derv
309  this%rhs(i) = -q + derv * this%xnew(node)
310 
311  end do
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:372
Here is the call graph for this function:

◆ cdb_create()

subroutine, public swfcdbmodule::cdb_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,
class(disbasetype), intent(inout), pointer  dis,
type(swfcxstype), intent(in), pointer  cxs,
real(dp), intent(in)  lengthconv,
real(dp), intent(in)  timeconv 
)

Create a new CDB Package object

Parameters
packobjpointer to default package type
[in]idpackage id
[in]ibcnumboundary condition number
[in]inunitunit number of CDB package input file
[in]ioutunit number of model listing file
[in]namemodelmodel name
[in]paknamepackage name
[in]mempathinput mempath
[in,out]disthe pointer to the discretization
[in]cxsthe pointer to the cxs package
[in]lengthconvconversion factor from model length to meters
[in]timeconvconversion factor from model time units to seconds

Definition at line 70 of file swf-cdb.f90.

72  ! -- dummy variables
73  class(BndType), pointer :: packobj !< pointer to default package type
74  integer(I4B), intent(in) :: id !< package id
75  integer(I4B), intent(in) :: ibcnum !< boundary condition number
76  integer(I4B), intent(in) :: inunit !< unit number of CDB package input file
77  integer(I4B), intent(in) :: iout !< unit number of model listing file
78  character(len=*), intent(in) :: namemodel !< model name
79  character(len=*), intent(in) :: pakname !< package name
80  character(len=*), intent(in) :: mempath !< input mempath
81  class(DisBaseType), pointer, intent(inout) :: dis !< the pointer to the discretization
82  type(SwfCxsType), pointer, intent(in) :: cxs !< the pointer to the cxs package
83  real(DP), intent(in) :: lengthconv !< conversion factor from model length to meters
84  real(DP), intent(in) :: timeconv !< conversion factor from model time units to seconds
85  ! -- local variables
86  type(SwfCdbType), pointer :: cdbobj
87  !
88  ! -- allocate the object and assign values to object variables
89  allocate (cdbobj)
90  packobj => cdbobj
91  !
92  ! -- create name and memory path
93  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
94  packobj%text = text
95  !
96  ! -- allocate scalars
97  call cdbobj%allocate_scalars()
98  !
99  ! -- initialize package
100  call packobj%pack_initialize()
101 
102  packobj%inunit = inunit
103  packobj%iout = iout
104  packobj%id = id
105  packobj%ibcnum = ibcnum
106  packobj%ncolbnd = 1
107  packobj%iscloc = 1
108  packobj%ictMemPath = create_mem_path(namemodel, 'DFW')
109 
110  ! -- store pointer to dis
111  cdbobj%dis => dis
112 
113  ! -- store pointer to cxs
114  cdbobj%cxs => cxs
115  !
116  ! -- store unit conversion
117  cdbobj%gravconv = dgravity * lengthconv * timeconv**2
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cdb_da()

subroutine swfcdbmodule::cdb_da ( class(swfcdbtype this)

Deallocate SWF package scalars and arrays.

Parameters
thisSwfcdbType object

Definition at line 175 of file swf-cdb.f90.

176  ! -- modules
178  ! -- dummy variables
179  class(SwfCdbType) :: this !< SwfcdbType object
180  !
181  ! -- Deallocate parent package
182  call this%BndExtType%bnd_da()
183  !
184  ! -- arrays
185  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
186  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
187  !
188  ! -- scalars
189  call mem_deallocate(this%gravconv)

◆ cdb_df_obs()

subroutine swfcdbmodule::cdb_df_obs ( class(swfcdbtype this)
private

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

Parameters
thisSwfCdbType object

Definition at line 434 of file swf-cdb.f90.

435  ! -- dummy variables
436  class(SwfCdbType) :: this !< SwfCdbType object
437  ! -- local variables
438  integer(I4B) :: indx
439  !
440  ! -- initialize observations
441  call this%obs%StoreObsType('cdb', .true., indx)
442  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
443  !
444  ! -- Store obs type and assign procedure pointer
445  ! for to-mvr observation type.
446  call this%obs%StoreObsType('to-mvr', .true., indx)
447  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ cdb_fc()

subroutine swfcdbmodule::cdb_fc ( class(swfcdbtype 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 CDB package to the coefficient matrix and right-hand side vector.

Parameters
thisSwfCdbType 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 351 of file swf-cdb.f90.

352  ! -- dummy variables
353  class(SwfCdbType) :: this !< SwfCdbType object
354  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
355  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
356  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
357  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
358  ! -- local variables
359  integer(I4B) :: i
360  integer(I4B) :: n
361  integer(I4B) :: ipos
362  !
363  ! -- pakmvrobj fc
364  if (this%imover == 1) then
365  call this%pakmvrobj%fc()
366  end if
367  !
368  ! -- Copy package rhs and hcof into solution rhs and amat
369  do i = 1, this%nbound
370  n = this%nodelist(i)
371  rhs(n) = rhs(n) + this%rhs(i)
372  ipos = ia(n)
373  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
374  !
375  ! -- If mover is active and this cdb item is discharging,
376  ! store available water (as positive value).
377  if (this%imover == 1 .and. this%rhs(i) > dzero) then
378  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
379  end if
380  end do

◆ cdb_obs_supported()

logical function swfcdbmodule::cdb_obs_supported ( class(swfcdbtype this)
private

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

Returns
cdb_obs_supported boolean indicating if observations are supported
Parameters
thisSwfCdbType object

Definition at line 421 of file swf-cdb.f90.

422  ! -- dummy variables
423  class(SwfCdbType) :: this !< SwfCdbType object
424  !
425  ! -- set boolean
426  cdb_obs_supported = .true.

◆ cdb_options()

subroutine swfcdbmodule::cdb_options ( class(swfcdbtype), intent(inout)  this)

Source additional options for SWF package.

Parameters
[in,out]thisSwfCdbType object

Definition at line 197 of file swf-cdb.f90.

198  ! -- modules
199  use inputoutputmodule, only: urword
202  ! -- dummy variables
203  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
204  ! -- local variables
205  type(SwfCdbParamFoundType) :: found
206  ! -- formats
207  !
208  ! -- source base BndExtType options
209  call this%BndExtType%source_options()
210  !
211  ! -- source options from input context
212  ! none
213  !
214  ! -- log SWF specific options
215  call this%log_cdb_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:

◆ cdb_rp()

subroutine swfcdbmodule::cdb_rp ( class(swfcdbtype), intent(inout)  this)

Definition at line 245 of file swf-cdb.f90.

246  use tdismodule, only: kper
247  ! -- dummy
248  class(SwfCdbType), intent(inout) :: this
249  ! -- local
250  !
251  if (this%iper /= kper) return
252  !
253  ! -- Call the parent class read and prepare
254  call this%BndExtType%bnd_rp()
255  !
256  ! -- Write the list to iout if requested
257  if (this%iprpak /= 0) then
258  call this%write_list()
259  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ cdb_rp_ts()

subroutine swfcdbmodule::cdb_rp_ts ( class(swfcdbtype), intent(inout)  this)
private

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

Parameters
[in,out]thisSwfCdbType object

Definition at line 505 of file swf-cdb.f90.

506  ! -- dummy variables
507  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
508  ! -- local variables
509  integer(I4B) :: i, nlinks
510  type(TimeSeriesLinkType), pointer :: tslink => null()
511  !
512  ! -- set up the time series links
513  nlinks = this%TsManager%boundtslinks%Count()
514  do i = 1, nlinks
515  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
516  if (associated(tslink)) then
517  if (tslink%JCol == 1) then
518  tslink%Text = 'Q'
519  end if
520  end if
521  end do
Here is the call graph for this function:

◆ define_listlabel()

subroutine swfcdbmodule::define_listlabel ( class(swfcdbtype), intent(inout)  this)
private

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

Parameters
[in,out]thisSwfCdbType object

Definition at line 389 of file swf-cdb.f90.

390  ! -- dummy variables
391  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
392  !
393  ! -- create the header list label
394  this%listlabel = trim(this%filtyp)//' NO.'
395  if (this%dis%ndim == 3) then
396  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
397  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
398  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
399  elseif (this%dis%ndim == 2) then
400  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
401  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
402  else
403  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
404  end if
405  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
406  if (this%inamedbound == 1) then
407  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
408  end if

◆ log_cdb_options()

subroutine swfcdbmodule::log_cdb_options ( class(swfcdbtype), intent(inout)  this,
type(swfcdbparamfoundtype), intent(in)  found 
)
Parameters
[in,out]thisBndExtType object

Definition at line 220 of file swf-cdb.f90.

221  ! -- modules
223  ! -- dummy variables
224  class(SwfCdbType), intent(inout) :: this !< BndExtType object
225  type(SwfCdbParamFoundType), intent(in) :: found
226  ! -- local variables
227  ! -- format
228  !
229  ! -- log found options
230  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
231  //' OPTIONS'
232  !
233  ! if (found%mover) then
234  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
235  ! end if
236  !
237  ! -- close logging block
238  write (this%iout, '(1x,a)') &
239  'END OF '//trim(adjustl(this%text))//' OPTIONS'

◆ qcalc()

real(dp) function swfcdbmodule::qcalc ( class(swfcdbtype this,
integer(i4b), intent(in)  i,
real(dp), intent(in)  depth 
)
Parameters
[in]iboundary number
[in]depthsimulated depth (stage - elevation) in reach n for this iteration

Definition at line 316 of file swf-cdb.f90.

317  ! modules
318  ! dummy
319  class(SwfCdbType) :: this
320  integer(I4B), intent(in) :: i !< boundary number
321  real(DP), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
322  ! return
323  real(DP) :: q
324  ! local
325  integer(I4B) :: idcxs
326  real(DP) :: width
327  real(DP) :: a
328  real(DP) :: r
329 
330  idcxs = this%idcxs(i)
331  width = this%width(i)
332  a = this%cxs%get_area(idcxs, width, depth)
333  r = this%cxs%get_hydraulic_radius(idcxs, width, depth, area=a)
334 
335  q = this%gravconv * a**dtwo * r
336  if (q > dprec) then
337  q = q**dhalf
338  else
339  q = dzero
340  end if
341  q = -q
342 

Variable Documentation

◆ ftype

character(len=lenftype) swfcdbmodule::ftype = 'CDB'
private

Definition at line 30 of file swf-cdb.f90.

30  character(len=LENFTYPE) :: ftype = 'CDB' !< package ftype

◆ text

character(len=16) swfcdbmodule::text = ' CDB'
private

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

31  character(len=16) :: text = ' CDB' !< package flow text string