30 character(len=LENFTYPE) ::
ftype =
'ZDG'
31 character(len=16) ::
text =
' ZDG'
35 integer(I4B),
dimension(:),
pointer,
contiguous :: idcxs => null()
36 real(dp),
dimension(:),
pointer,
contiguous :: width => null()
37 real(dp),
dimension(:),
pointer,
contiguous :: slope => null()
38 real(dp),
dimension(:),
pointer,
contiguous :: rough => null()
39 real(dp),
pointer :: unitconv => null()
72 subroutine zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
73 mempath, dis, cxs, unitconv)
75 class(
bndtype),
pointer :: packobj
76 integer(I4B),
intent(in) :: id
77 integer(I4B),
intent(in) :: ibcnum
78 integer(I4B),
intent(in) :: inunit
79 integer(I4B),
intent(in) :: iout
80 character(len=*),
intent(in) :: namemodel
81 character(len=*),
intent(in) :: pakname
82 character(len=*),
intent(in) :: mempath
85 real(dp),
intent(in) :: unitconv
94 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
98 call zdgobj%allocate_scalars()
101 call packobj%pack_initialize()
103 packobj%inunit = inunit
106 packobj%ibcnum = ibcnum
121 zdgobj%unitconv = unitconv
137 call this%BndExtType%allocate_scalars()
140 call mem_allocate(this%unitconv,
'UNITCONV', this%memoryPath)
143 this%unitconv =
dzero
156 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
157 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
161 call this%BndExtType%allocate_arrays(nodelist, auxvar)
164 call mem_setptr(this%idcxs,
'IDCXS', this%input_mempath)
165 call mem_setptr(this%width,
'WIDTH', this%input_mempath)
166 call mem_setptr(this%slope,
'SLOPE', this%input_mempath)
167 call mem_setptr(this%rough,
'ROUGH', this%input_mempath)
170 call mem_checkin(this%idcxs,
'IDCXS', this%memoryPath, &
171 'IDCXS', this%input_mempath)
172 call mem_checkin(this%width,
'WIDTH', this%memoryPath, &
173 'WIDTH', this%input_mempath)
174 call mem_checkin(this%slope,
'SLOPE', this%memoryPath, &
175 'SLOPE', this%input_mempath)
176 call mem_checkin(this%rough,
'ROUGH', this%memoryPath, &
177 'ROUGH', this%input_mempath)
192 call this%BndExtType%bnd_da()
221 call this%BndExtType%source_options()
227 call this%log_zdg_options(found)
242 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
250 write (this%iout,
'(1x,a)') &
251 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
267 integer(I4B) :: i, node
270 real(DP) :: absdhdxsq
274 real(DP) :: range = 1.d-6
276 real(DP) :: smooth_factor
279 if (this%nbound == 0)
return
282 do i = 1, this%nbound
284 node = this%nodelist(i)
285 if (this%ibound(node) <= 0)
then
292 absdhdxsq = this%slope(i)**
dhalf
293 depth = this%xnew(node) - this%dis%bot(node)
296 call squadratic(depth, range, dydx, smooth_factor)
297 depth = depth * smooth_factor
300 q = -this%qcalc(i, depth, this%unitconv)
304 qeps = -this%qcalc(i, depth + eps, this%unitconv)
307 derv = (qeps - q) / eps
311 this%rhs(i) = -q + derv * this%xnew(node)
321 function qcalc(this, i, depth, unitconv)
result(q)
324 integer(I4B),
intent(in) :: i
325 real(dp),
intent(in) :: depth
326 real(dp),
intent(in) :: unitconv
330 integer(I4B) :: idcxs
334 real(dp) :: conveyance
336 idcxs = this%idcxs(i)
337 width = this%width(i)
338 rough = this%rough(i)
339 slope = this%slope(i)
340 conveyance = this%cxs%get_conveyance(idcxs, width, depth, rough)
341 q = conveyance * slope**
dhalf * unitconv
351 subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
354 real(DP),
dimension(:),
intent(inout) :: rhs
355 integer(I4B),
dimension(:),
intent(in) :: ia
356 integer(I4B),
dimension(:),
intent(in) :: idxglo
364 if (this%imover == 1)
then
365 call this%pakmvrobj%fc()
369 do i = 1, this%nbound
371 rhs(n) = rhs(n) + this%rhs(i)
373 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
377 if (this%imover == 1 .and. this%rhs(i) >
dzero)
then
378 call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
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'
403 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
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'
441 call this%obs%StoreObsType(
'zdg', .true., indx)
446 call this%obs%StoreObsType(
'to-mvr', .true., indx)
466 call this%obs%obs_bd_clear()
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
474 jj = obsrv%indxbnds(n)
475 select case (obsrv%ObsTypeId)
477 if (this%imover == 1)
then
478 v = this%pakmvrobj%get_qtomvr(jj)
486 errmsg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
489 call this%obs%SaveOneSimval(obsrv, v)
492 call this%obs%SaveOneSimval(obsrv,
dnodata)
509 integer(I4B) :: i, nlinks
513 nlinks = this%TsManager%boundtslinks%Count()
516 if (
associated(tslink))
then
517 if (tslink%JCol == 1)
then
535 integer(I4B),
intent(in) :: col
536 integer(I4B),
intent(in) :: row
542 bndval = this%idcxs(row)
544 bndval = this%width(row)
546 bndval = this%slope(row)
548 bndval = this%rough(row)
550 errmsg =
'Programming error. ZDG bound value requested column '&
551 &
'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 dzero
real constant zero
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
subroutine squadratic(x, range, dydx, y)
@ brief sQuadratic
This module contains the ZDG package methods.
subroutine zdg_allocate_scalars(this)
@ brief Allocate scalars
subroutine define_listlabel(this)
@ brief Define the list label for the package
logical function zdg_obs_supported(this)
Determine if observations are supported.
subroutine zdg_cf(this)
@ brief Formulate the package hcof and rhs terms.
subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Copy hcof and rhs terms into solution.
real(dp) function zdg_bound_value(this, col, row)
@ brief Return a bound value
subroutine log_zdg_options(this, found)
@ brief Log SWF specific package options
subroutine, public zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, dis, cxs, unitconv)
@ brief Create a new package object
subroutine zdg_rp_ts(this)
Assign time series links for the package.
character(len=16) text
package flow text string
subroutine zdg_da(this)
@ brief Deallocate package memory
subroutine zdg_options(this)
@ brief Source additional options for package
character(len=lenftype) ftype
package ftype
subroutine zdg_allocate_arrays(this, nodelist, auxvar)
@ brief Allocate arrays
subroutine zdg_bd_obs(this)
Save observations for the package.
real(dp) function qcalc(this, i, depth, unitconv)
subroutine zdg_df_obs(this)
Define the observation types available in the package.
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.