MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 156 of file swf-cdb.f90.

157  ! -- modules
159  ! -- dummy
160  class(SwfCdbType) :: this
161  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
162  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
163  ! -- local
164  !
165  ! -- call BndExtType allocate scalars
166  call this%BndExtType%allocate_arrays(nodelist, auxvar)
167  !
168  ! -- set array input context pointer
169  call mem_setptr(this%idcxs, 'IDCXS', this%input_mempath)
170  call mem_setptr(this%width, 'WIDTH', this%input_mempath)
171  !
172  ! -- checkin array input context pointer
173  call mem_checkin(this%idcxs, 'IDCXS', this%memoryPath, &
174  'IDCXS', this%input_mempath)
175  call mem_checkin(this%width, 'WIDTH', this%memoryPath, &
176  'WIDTH', this%input_mempath)
177  !
178  ! -- return
179  return

◆ 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 132 of file swf-cdb.f90.

133  ! -- modules
135  ! -- dummy variables
136  class(SwfCdbType) :: this !< SwfcdbType object
137  !
138  ! -- call base type allocate scalars
139  call this%BndExtType%allocate_scalars()
140  !
141  ! -- allocate the object and assign values to object variables
142  call mem_allocate(this%gravconv, 'GRAVCONV', this%memoryPath)
143  !
144  ! -- Set values
145  this%gravconv = dzero
146  !
147  ! -- return
148  return

◆ 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 493 of file swf-cdb.f90.

494  ! -- dummy variables
495  class(SwfCdbType) :: this !< SwfCdbType object
496  ! -- local variables
497  integer(I4B) :: i
498  integer(I4B) :: n
499  integer(I4B) :: jj
500  real(DP) :: v
501  type(ObserveType), pointer :: obsrv => null()
502  !
503  ! -- clear the observations
504  call this%obs%obs_bd_clear()
505  !
506  ! -- Save simulated values for all of package's observations.
507  do i = 1, this%obs%npakobs
508  obsrv => this%obs%pakobs(i)%obsrv
509  if (obsrv%BndFound) then
510  do n = 1, obsrv%indxbnds_count
511  v = dnodata
512  jj = obsrv%indxbnds(n)
513  select case (obsrv%ObsTypeId)
514  case ('TO-MVR')
515  if (this%imover == 1) then
516  v = this%pakmvrobj%get_qtomvr(jj)
517  if (v > dzero) then
518  v = -v
519  end if
520  end if
521  case ('CDB')
522  v = this%simvals(jj)
523  case default
524  errmsg = 'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
525  call store_error(errmsg)
526  end select
527  call this%obs%SaveOneSimval(obsrv, v)
528  end do
529  else
530  call this%obs%SaveOneSimval(obsrv, dnodata)
531  end if
532  end do
533  !
534  ! -- return
535  return
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 574 of file swf-cdb.f90.

575  ! -- modules
576  use constantsmodule, only: dzero
577  ! -- dummy variables
578  class(SwfCdbType), intent(inout) :: this !< BndExtType object
579  integer(I4B), intent(in) :: col
580  integer(I4B), intent(in) :: row
581  ! -- result
582  real(DP) :: bndval
583  !
584  select case (col)
585  case (1)
586  bndval = this%idcxs(row)
587  case (2)
588  bndval = this%width(row)
589  case default
590  errmsg = 'Programming error. CDB bound value requested column '&
591  &'outside range of ncolbnd (1).'
592  call store_error(errmsg)
593  call store_error_filename(this%input_fname)
594  end select
595  !
596  ! -- return
597  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
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 292 of file swf-cdb.f90.

293  ! modules
295  ! -- dummy variables
296  class(SwfCdbType) :: this !< SwfCdbType object
297  ! -- local variables
298  integer(I4B) :: i, node
299  real(DP) :: q
300  real(DP) :: qeps
301  real(DP) :: depth
302  real(DP) :: derv
303  real(DP) :: eps
304  !
305  ! -- Return if no inflows
306  if (this%nbound == 0) return
307  !
308  ! -- Calculate hcof and rhs for each cdb entry
309  do i = 1, this%nbound
310 
311  node = this%nodelist(i)
312  if (this%ibound(node) <= 0) then
313  this%hcof(i) = dzero
314  this%rhs(i) = dzero
315  cycle
316  end if
317 
318  ! -- calculate terms and add to hcof and rhs
319  depth = this%xnew(node) - this%dis%bot(node)
320 
321  ! -- calculate unperturbed q
322  q = this%qcalc(i, depth)
323 
324  ! -- calculate perturbed q
325  eps = get_perturbation(depth)
326  qeps = this%qcalc(i, depth + eps)
327 
328  ! -- calculate derivative
329  derv = (qeps - q) / eps
330 
331  ! -- add terms to hcof and rhs
332  this%hcof(i) = derv
333  this%rhs(i) = -q + derv * this%xnew(node)
334 
335  end do
336  !
337  return
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
Definition: MathUtil.f90:446
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
121  !
122  ! -- return
123  return
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 187 of file swf-cdb.f90.

188  ! -- modules
190  ! -- dummy variables
191  class(SwfCdbType) :: this !< SwfcdbType object
192  !
193  ! -- Deallocate parent package
194  call this%BndExtType%bnd_da()
195  !
196  ! -- arrays
197  call mem_deallocate(this%idcxs, 'IDCXS', this%memoryPath)
198  call mem_deallocate(this%width, 'WIDTH', this%memoryPath)
199  !
200  ! -- scalars
201  call mem_deallocate(this%gravconv)
202  !
203  ! -- return
204  return

◆ 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 469 of file swf-cdb.f90.

470  ! -- dummy variables
471  class(SwfCdbType) :: this !< SwfCdbType object
472  ! -- local variables
473  integer(I4B) :: indx
474  !
475  ! -- initialize observations
476  call this%obs%StoreObsType('cdb', .true., indx)
477  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
478  !
479  ! -- Store obs type and assign procedure pointer
480  ! for to-mvr observation type.
481  call this%obs%StoreObsType('to-mvr', .true., indx)
482  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
483  !
484  ! -- return
485  return
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 377 of file swf-cdb.f90.

378  ! -- dummy variables
379  class(SwfCdbType) :: this !< SwfCdbType object
380  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
381  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
382  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
383  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
384  ! -- local variables
385  integer(I4B) :: i
386  integer(I4B) :: n
387  integer(I4B) :: ipos
388  !
389  ! -- pakmvrobj fc
390  if (this%imover == 1) then
391  call this%pakmvrobj%fc()
392  end if
393  !
394  ! -- Copy package rhs and hcof into solution rhs and amat
395  do i = 1, this%nbound
396  n = this%nodelist(i)
397  rhs(n) = rhs(n) + this%rhs(i)
398  ipos = ia(n)
399  call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
400  !
401  ! -- If mover is active and this cdb item is discharging,
402  ! store available water (as positive value).
403  if (this%imover == 1 .and. this%rhs(i) > dzero) then
404  call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
405  end if
406  end do
407  !
408  ! -- return
409  return

◆ 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 453 of file swf-cdb.f90.

454  ! -- dummy variables
455  class(SwfCdbType) :: this !< SwfCdbType object
456  !
457  ! -- set boolean
458  cdb_obs_supported = .true.
459  !
460  ! -- return
461  return

◆ 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 212 of file swf-cdb.f90.

213  ! -- modules
214  use inputoutputmodule, only: urword
217  ! -- dummy variables
218  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
219  ! -- local variables
220  type(SwfCdbParamFoundType) :: found
221  ! -- formats
222  !
223  ! -- source base BndExtType options
224  call this%BndExtType%source_options()
225  !
226  ! -- source options from input context
227  ! none
228  !
229  ! -- log SWF specific options
230  call this%log_cdb_options(found)
231  !
232  ! -- return
233  return
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 266 of file swf-cdb.f90.

267  use tdismodule, only: kper
268  ! -- dummy
269  class(SwfCdbType), intent(inout) :: this
270  ! -- local
271  !
272  if (this%iper /= kper) return
273  !
274  ! -- Call the parent class read and prepare
275  call this%BndExtType%bnd_rp()
276  !
277  ! -- Write the list to iout if requested
278  if (this%iprpak /= 0) then
279  call this%write_list()
280  end if
281  !
282  ! -- return
283  return
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 546 of file swf-cdb.f90.

547  ! -- dummy variables
548  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
549  ! -- local variables
550  integer(I4B) :: i, nlinks
551  type(TimeSeriesLinkType), pointer :: tslink => null()
552  !
553  ! -- set up the time series links
554  nlinks = this%TsManager%boundtslinks%Count()
555  do i = 1, nlinks
556  tslink => gettimeserieslinkfromlist(this%TsManager%boundtslinks, i)
557  if (associated(tslink)) then
558  if (tslink%JCol == 1) then
559  tslink%Text = 'Q'
560  end if
561  end if
562  end do
563  !
564  ! -- return
565  return
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 418 of file swf-cdb.f90.

419  ! -- dummy variables
420  class(SwfCdbType), intent(inout) :: this !< SwfCdbType object
421  !
422  ! -- create the header list label
423  this%listlabel = trim(this%filtyp)//' NO.'
424  if (this%dis%ndim == 3) then
425  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
426  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
427  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
428  elseif (this%dis%ndim == 2) then
429  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
430  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
431  else
432  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
433  end if
434  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'FLOW RATE'
435  if (this%inamedbound == 1) then
436  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
437  end if
438  !
439  ! -- return
440  return

◆ 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 238 of file swf-cdb.f90.

239  ! -- modules
241  ! -- dummy variables
242  class(SwfCdbType), intent(inout) :: this !< BndExtType object
243  type(SwfCdbParamFoundType), intent(in) :: found
244  ! -- local variables
245  ! -- format
246  !
247  ! -- log found options
248  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%text)) &
249  //' OPTIONS'
250  !
251  ! if (found%mover) then
252  ! write (this%iout, '(4x,A)') 'MOVER OPTION ENABLED'
253  ! end if
254  !
255  ! -- close logging block
256  write (this%iout, '(1x,a)') &
257  'END OF '//trim(adjustl(this%text))//' OPTIONS'
258  !
259  ! -- return
260  return

◆ 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 342 of file swf-cdb.f90.

343  ! modules
344  ! dummy
345  class(SwfCdbType) :: this
346  integer(I4B), intent(in) :: i !< boundary number
347  real(DP), intent(in) :: depth !< simulated depth (stage - elevation) in reach n for this iteration
348  ! return
349  real(DP) :: q
350  ! local
351  integer(I4B) :: idcxs
352  real(DP) :: width
353  real(DP) :: a
354  real(DP) :: r
355 
356  idcxs = this%idcxs(i)
357  width = this%width(i)
358  a = this%cxs%get_area(idcxs, width, depth)
359  r = this%cxs%get_hydraulic_radius(idcxs, width, depth, area=a)
360 
361  q = this%gravconv * a**dtwo * r
362  if (q > dprec) then
363  q = q**dhalf
364  else
365  q = dzero
366  end if
367  q = -q
368 

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