MODFLOW 6  version 6.6.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 150 of file swf-cdb.f90.

151  ! -- modules
153  ! -- dummy
154  class(SwfCdbType) :: this
155  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
156  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
157  ! -- local
158  !
159  ! -- call BndExtType allocate scalars
160  call this%BndExtType%allocate_arrays(nodelist, auxvar)
161  !
162  ! -- set array input context pointer
163  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
164  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
165  !
166  ! -- checkin array input context pointer
167  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
168  'IDCXS', this%input_mempath)
169  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
170  '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 129 of file swf-cdb.f90.

130  ! -- modules
132  ! -- dummy variables
133  class(SwfCdbType) :: this !< SwfcdbType object
134  !
135  ! -- call base type allocate scalars
136  call this%BndExtType%allocate_scalars()
137  !
138  ! -- allocate the object and assign values to object variables
139  call mem_allocate(this%gravconv, 'GRAVCONV', this%memoryPath)
140  !
141  ! -- Set values
142  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 458 of file swf-cdb.f90.

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

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

272  ! modules
274  ! -- dummy variables
275  class(SwfCdbType) :: this !< SwfCdbType object
276  ! -- local variables
277  integer(I4B) :: i, node
278  real(DP) :: q
279  real(DP) :: qeps
280  real(DP) :: depth
281  real(DP) :: derv
282  real(DP) :: eps
283  !
284  ! -- Return if no inflows
285  if (this%nbound == 0) return
286  !
287  ! -- Calculate hcof and rhs for each cdb entry
288  do i = 1, this%nbound
289 
290  node = this%nodelist(i)
291  if (this%ibound(node) <= 0) then
292  this%hcof(i) = dzero
293  this%rhs(i) = dzero
294  cycle
295  end if
296 
297  ! -- calculate terms and add to hcof and rhs
298  depth = this%xnew(node) - this%dis%bot(node)
299 
300  ! -- calculate unperturbed q
301  q = this%qcalc(i, depth)
302 
303  ! -- calculate perturbed q
304  eps = get_perturbation(depth)
305  qeps = this%qcalc(i, depth + eps)
306 
307  ! -- calculate derivative
308  derv = (qeps - q) / eps
309 
310  ! -- add terms to hcof and rhs
311  this%hcof(i) = derv
312  this%rhs(i) = -q + derv * this%xnew(node)
313 
314  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 73 of file swf-cdb.f90.

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

179  ! -- modules
181  ! -- dummy variables
182  class(SwfCdbType) :: this !< SwfcdbType object
183  !
184  ! -- Deallocate parent package
185  call this%BndExtType%bnd_da()
186  !
187  ! -- arrays
188  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
189  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
190  !
191  ! -- scalars
192  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 437 of file swf-cdb.f90.

438  ! -- dummy variables
439  class(SwfCdbType) :: this !< SwfCdbType object
440  ! -- local variables
441  integer(I4B) :: indx
442  !
443  ! -- initialize observations
444  call this%obs%StoreObsType('cdb', .true., indx)
445  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
446  !
447  ! -- Store obs type and assign procedure pointer
448  ! for to-mvr observation type.
449  call this%obs%StoreObsType('to-mvr', .true., indx)
450  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 354 of file swf-cdb.f90.

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

425  ! -- dummy variables
426  class(SwfCdbType) :: this !< SwfCdbType object
427  !
428  ! -- set boolean
429  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 200 of file swf-cdb.f90.

201  ! -- modules
202  use inputoutputmodule, only: urword
205  ! -- dummy variables
206  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
207  ! -- local variables
208  type(SwfCdbParamFoundType) :: found
209  ! -- formats
210  !
211  ! -- source base BndExtType options
212  call this%BndExtType%source_options()
213  !
214  ! -- source options from input context
215  ! none
216  !
217  ! -- log SWF specific options
218  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 248 of file swf-cdb.f90.

249  use tdismodule, only: kper
250  ! -- dummy
251  class(SwfCdbType), intent(inout) :: this
252  ! -- local
253  !
254  if (this%iper /= kper) return
255  !
256  ! -- Call the parent class read and prepare
257  call this%BndExtType%bnd_rp()
258  !
259  ! -- Write the list to iout if requested
260  if (this%iprpak /= 0) then
261  call this%write_list()
262  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 508 of file swf-cdb.f90.

509  ! -- dummy variables
510  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
511  ! -- local variables
512  integer(I4B) :: i, nlinks
513  type(TimeSeriesLinkType), pointer :: tslink => null()
514  !
515  ! -- set up the time series links
516  nlinks = this%TsManager%boundtslinks%Count()
517  do i = 1, nlinks
518  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
519  if (associated(tslink)) then
520  if (tslink%JCol == 1) then
521  tslink%Text = 'Q'
522  end if
523  end if
524  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 392 of file swf-cdb.f90.

393  ! -- dummy variables
394  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
395  !
396  ! -- create the header list label
397  this%listlabel = trim(this%filtyp)//' NO.'
398  if (this%dis%ndim == 3) then
399  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
400  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
401  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
402  elseif (this%dis%ndim == 2) then
403  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
404  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
405  else
406  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
407  end if
408  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
409  if (this%inamedbound == 1) then
410  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
411  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 223 of file swf-cdb.f90.

224  ! -- modules
226  ! -- dummy variables
227  class(SwfCdbType), intent(inout) :: this !< BndExtType object
228  type(SwfCdbParamFoundType), intent(in) :: found
229  ! -- local variables
230  ! -- format
231  !
232  ! -- log found options
233  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
234  //' OPTIONS'
235  !
236  ! if (found%mover) then
237  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
238  ! end if
239  !
240  ! -- close logging block
241  write (this%iout, '(1x,a)') &
242  '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 319 of file swf-cdb.f90.

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

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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