MODFLOW 6  version 6.8.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_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 146 of file swf-cdb.f90.

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

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

435  ! -- dummy variables
436  class(SwfCdbType) :: this !< SwfCdbType object
437  ! -- local variables
438  integer(I4B) :: i
439  integer(I4B) :: n
440  integer(I4B) :: jj
441  real(DP) :: v
442  type(ObserveType), pointer :: obsrv => null()
443  !
444  ! -- clear the observations
445  call this%obs%obs_bd_clear()
446  !
447  ! -- Save simulated values for all of package's observations.
448  do i = 1, this%obs%npakobs
449  obsrv => this%obs%pakobs(i)%obsrv
450  if (obsrv%BndFound) then
451  do n = 1, obsrv%indxbnds_count
452  v = dnodata
453  jj = obsrv%indxbnds(n)
454  select case (obsrv%ObsTypeId)
455  case ('TO-MVR')
456  if (this%imover == 1) then
457  v = this%pakmvrobj%get_qtomvr(jj)
458  if (v > dzero) then
459  v = -v
460  end if
461  end if
462  case ('CDB')
463  v = this%simvals(jj)
464  case default
465  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
466  call store_error(errmsg)
467  end select
468  call this%obs%SaveOneSimval(obsrv, v)
469  end do
470  else
471  call this%obs%SaveOneSimval(obsrv, dnodata)
472  end if
473  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.

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

510  ! -- modules
511  use constantsmodule, only: dzero
512  ! -- dummy variables
513  class(SwfCdbType), intent(inout) :: this
514  integer(I4B), intent(in) :: col
515  integer(I4B), intent(in) :: row
516  ! -- result
517  real(DP) :: bndval
518  !
519  select case (col)
520  case (1)
521  bndval = this%idcxs(row)
522  case (2)
523  bndval = this%width(row)
524  case default
525  errmsg = 'Programming error. CDB bound value requested column '&
526  &'outside range of ncolbnd (1).'
527  call store_error(errmsg)
528  call store_error_filename(this%input_fname)
529  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 247 of file swf-cdb.f90.

248  ! modules
250  ! -- dummy variables
251  class(SwfCdbType) :: this !< SwfCdbType object
252  ! -- local variables
253  integer(I4B) :: i, node
254  real(DP) :: q
255  real(DP) :: qeps
256  real(DP) :: depth
257  real(DP) :: derv
258  real(DP) :: eps
259  !
260  ! -- Return if no inflows
261  if (this%nbound == 0) return
262  !
263  ! -- Calculate hcof and rhs for each cdb entry
264  do i = 1, this%nbound
265 
266  node = this%nodelist(i)
267  if (this%ibound(node) <= 0) then
268  this%hcof(i) = dzero
269  this%rhs(i) = dzero
270  cycle
271  end if
272 
273  ! -- calculate terms and add to hcof and rhs
274  depth = this%xnew(node) - this%dis%bot(node)
275 
276  ! -- calculate unperturbed q
277  q = this%qcalc(i, depth)
278 
279  ! -- calculate perturbed q
280  eps = get_perturbation(depth)
281  qeps = this%qcalc(i, depth + eps)
282 
283  ! -- calculate derivative
284  derv = (qeps - q) / eps
285 
286  ! -- add terms to hcof and rhs
287  this%hcof(i) = derv
288  this%rhs(i) = -q + derv * this%xnew(node)
289 
290  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 69 of file swf-cdb.f90.

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

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

414  ! -- dummy variables
415  class(SwfCdbType) :: this !< SwfCdbType object
416  ! -- local variables
417  integer(I4B) :: indx
418  !
419  ! -- initialize observations
420  call this%obs%StoreObsType('cdb', .true., indx)
421  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
422  !
423  ! -- Store obs type and assign procedure pointer
424  ! for to-mvr observation type.
425  call this%obs%StoreObsType('to-mvr', .true., indx)
426  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 330 of file swf-cdb.f90.

331  ! -- dummy variables
332  class(SwfCdbType) :: this !< SwfCdbType object
333  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
334  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
335  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
336  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
337  ! -- local variables
338  integer(I4B) :: i
339  integer(I4B) :: n
340  integer(I4B) :: ipos
341  !
342  ! -- pakmvrobj fc
343  if (this%imover == 1) then
344  call this%pakmvrobj%fc()
345  end if
346  !
347  ! -- Copy package rhs and hcof into solution rhs and amat
348  do i = 1, this%nbound
349  n = this%nodelist(i)
350  rhs(n) = rhs(n) + this%rhs(i)
351  ipos = ia(n)
352  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
353  !
354  ! -- If mover is active and this cdb item is discharging,
355  ! store available water (as positive value).
356  if (this%imover == 1 .and. this%rhs(i) > dzero) then
357  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
358  end if
359  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 400 of file swf-cdb.f90.

401  ! -- dummy variables
402  class(SwfCdbType) :: this !< SwfCdbType object
403  !
404  ! -- set boolean
405  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 196 of file swf-cdb.f90.

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

485  ! -- dummy variables
486  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
487  ! -- local variables
488  integer(I4B) :: i, nlinks
489  type(TimeSeriesLinkType), pointer :: tslink => null()
490  !
491  ! -- set up the time series links
492  nlinks = this%TsManager%boundtslinks%Count()
493  do i = 1, nlinks
494  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
495  if (associated(tslink)) then
496  if (tslink%JCol == 1) then
497  tslink%Text = 'Q'
498  end if
499  end if
500  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 368 of file swf-cdb.f90.

369  ! -- dummy variables
370  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
371  !
372  ! -- create the header list label
373  this%listlabel = trim(this%filtyp)//' NO.'
374  if (this%dis%ndim == 3) then
375  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
376  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
377  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
378  elseif (this%dis%ndim == 2) then
379  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
380  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
381  else
382  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
383  end if
384  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
385  if (this%inamedbound == 1) then
386  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
387  end if

◆ log_cdb_options()

subroutine swfcdbmodule::log_cdb_options ( class(swfcdbtype), intent(inout)  this,
type(swfcdbparamfoundtype), intent(in)  found 
)

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

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

296  ! modules
297  ! dummy
298  class(SwfCdbType) :: this
299  integer(I4B), intent(in) :: i !< boundary number
300  real(DP), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
301  ! return
302  real(DP) :: q
303  ! local
304  integer(I4B) :: idcxs
305  real(DP) :: width
306  real(DP) :: a
307  real(DP) :: r
308 
309  idcxs = this%idcxs(i)
310  width = this%width(i)
311  a = this%cxs%get_area(idcxs, width, depth)
312  r = this%cxs%get_hydraulic_radius(idcxs, width, depth, area=a)
313 
314  q = this%gravconv * a**dtwo * r
315  if (q > dprec) then
316  q = q**dhalf
317  else
318  q = dzero
319  end if
320  q = -q
321 

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