30 character(len=LENFTYPE) ::
ftype =
'CDB'
31 character(len=16) ::
text =
' CDB'
35 integer(I4B),
dimension(:),
pointer,
contiguous :: idcxs => null()
36 real(dp),
dimension(:),
pointer,
contiguous :: width => null()
37 real(dp),
pointer :: gravconv => null()
69 subroutine cdb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
70 mempath, dis, cxs, lengthconv, timeconv)
72 class(
bndtype),
pointer :: packobj
73 integer(I4B),
intent(in) :: id
74 integer(I4B),
intent(in) :: ibcnum
75 integer(I4B),
intent(in) :: inunit
76 integer(I4B),
intent(in) :: iout
77 character(len=*),
intent(in) :: namemodel
78 character(len=*),
intent(in) :: pakname
79 character(len=*),
intent(in) :: mempath
82 real(dp),
intent(in) :: lengthconv
83 real(dp),
intent(in) :: timeconv
92 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
96 call cdbobj%allocate_scalars()
99 call packobj%pack_initialize()
101 packobj%inunit = inunit
104 packobj%ibcnum = ibcnum
116 cdbobj%gravconv =
dgravity * lengthconv * timeconv**2
132 call this%BndExtType%allocate_scalars()
135 call mem_allocate(this%gravconv,
'GRAVCONV', this%memoryPath)
138 this%gravconv =
dzero
151 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
152 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
156 call this%BndExtType%allocate_arrays(nodelist, auxvar)
159 call mem_setptr(this%idcxs,
'IDCXS', this%input_mempath)
160 call mem_setptr(this%width,
'WIDTH', this%input_mempath)
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)
181 call this%BndExtType%bnd_da()
208 call this%BndExtType%source_options()
214 call this%log_cdb_options(found)
229 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
237 write (this%iout,
'(1x,a)') &
238 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
253 integer(I4B) :: i, node
261 if (this%nbound == 0)
return
264 do i = 1, this%nbound
266 node = this%nodelist(i)
267 if (this%ibound(node) <= 0)
then
274 depth = this%xnew(node) - this%dis%bot(node)
277 q = this%qcalc(i, depth)
281 qeps = this%qcalc(i, depth + eps)
284 derv = (qeps - q) / eps
288 this%rhs(i) = -q + derv * this%xnew(node)
295 function qcalc(this, i, depth)
result(q)
299 integer(I4B),
intent(in) :: i
300 real(dp),
intent(in) :: depth
304 integer(I4B) :: idcxs
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)
314 q = this%gravconv * a**
dtwo * r
330 subroutine cdb_fc(this, rhs, ia, idxglo, matrix_sln)
333 real(DP),
dimension(:),
intent(inout) :: rhs
334 integer(I4B),
dimension(:),
intent(in) :: ia
335 integer(I4B),
dimension(:),
intent(in) :: idxglo
343 if (this%imover == 1)
then
344 call this%pakmvrobj%fc()
348 do i = 1, this%nbound
350 rhs(n) = rhs(n) + this%rhs(i)
352 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
356 if (this%imover == 1 .and. this%rhs(i) >
dzero)
then
357 call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
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'
382 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
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'
420 call this%obs%StoreObsType(
'cdb', .true., indx)
425 call this%obs%StoreObsType(
'to-mvr', .true., indx)
445 call this%obs%obs_bd_clear()
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
453 jj = obsrv%indxbnds(n)
454 select case (obsrv%ObsTypeId)
456 if (this%imover == 1)
then
457 v = this%pakmvrobj%get_qtomvr(jj)
465 errmsg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
468 call this%obs%SaveOneSimval(obsrv, v)
471 call this%obs%SaveOneSimval(obsrv,
dnodata)
488 integer(I4B) :: i, nlinks
492 nlinks = this%TsManager%boundtslinks%Count()
495 if (
associated(tslink))
then
496 if (tslink%JCol == 1)
then
514 integer(I4B),
intent(in) :: col
515 integer(I4B),
intent(in) :: row
521 bndval = this%idcxs(row)
523 bndval = this%width(row)
525 errmsg =
'Programming error. CDB bound value requested column '&
526 &
'outside range of ncolbnd (1).'
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
real(dp), parameter dnodata
real no data constant
real(dp), parameter dhalf
real constant 1/2
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
real(dp), parameter dgravity
real constant gravitational acceleration (m/(s s))
real(dp), parameter dzero
real constant zero
real(dp), parameter dprec
real constant machine precision
real(dp), parameter dtwo
real constant 2
This module defines variable data types.
real(dp) function, public get_perturbation(x)
Calculate a numerical perturbation given the value of x.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived types ObserveType and ObsDataType.
This module contains the derived type ObsType.
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
This module contains the CDB package methods.
subroutine cdb_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
subroutine, public cdb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, lengthconv, timeconv)
@ brief Create a new package object
real(dp) function qcalc(this, i, depth)
Calculate critical depth boundary flow.
logical function cdb_obs_supported(this)
Determine if observations are supported.
subroutine log_cdb_options(this, found)
@ brief Log SWF specific package options
subroutine cdb_bd_obs(this)
Save observations for the package.
character(len=16) text
package flow text string
character(len=lenftype) ftype
package ftype
subroutine cdb_rp_ts(this)
Assign time series links for the package.
subroutine cdb_da(this)
@ brief Deallocate package memory
subroutine cdb_allocate_scalars(this)
@ brief Allocate scalars
subroutine cdb_options(this)
@ brief Source additional options for package
subroutine cdb_cf(this)
@ brief Formulate the package hcof and rhs terms.
subroutine cdb_df_obs(this)
Define the observation types available in the package.
subroutine cdb_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
subroutine define_listlabel(this)
@ brief Define the list label for the package
real(dp) function cdb_bound_value(this, col, row)
@ brief Return a bound value
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.