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()
73 subroutine zdg_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
74 mempath, dis, cxs, unitconv)
76 class(
bndtype),
pointer :: packobj
77 integer(I4B),
intent(in) :: id
78 integer(I4B),
intent(in) :: ibcnum
79 integer(I4B),
intent(in) :: inunit
80 integer(I4B),
intent(in) :: iout
81 character(len=*),
intent(in) :: namemodel
82 character(len=*),
intent(in) :: pakname
83 character(len=*),
intent(in) :: mempath
86 real(dp),
intent(in) :: unitconv
95 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
99 call zdgobj%allocate_scalars()
102 call packobj%pack_initialize()
104 packobj%inunit = inunit
107 packobj%ibcnum = ibcnum
122 zdgobj%unitconv = unitconv
138 call this%BndExtType%allocate_scalars()
141 call mem_allocate(this%unitconv,
'UNITCONV', this%memoryPath)
144 this%unitconv =
dzero
157 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
158 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
162 call this%BndExtType%allocate_arrays(nodelist, auxvar)
165 call mem_setptr(this%idcxs,
'IDCXS', this%input_mempath)
166 call mem_setptr(this%width,
'WIDTH', this%input_mempath)
167 call mem_setptr(this%slope,
'SLOPE', this%input_mempath)
168 call mem_setptr(this%rough,
'ROUGH', this%input_mempath)
171 call mem_checkin(this%idcxs,
'IDCXS', this%memoryPath, &
172 'IDCXS', this%input_mempath)
173 call mem_checkin(this%width,
'WIDTH', this%memoryPath, &
174 'WIDTH', this%input_mempath)
175 call mem_checkin(this%slope,
'SLOPE', this%memoryPath, &
176 'SLOPE', this%input_mempath)
177 call mem_checkin(this%rough,
'ROUGH', this%memoryPath, &
178 'ROUGH', this%input_mempath)
193 call this%BndExtType%bnd_da()
222 call this%BndExtType%source_options()
228 call this%log_zdg_options(found)
243 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
251 write (this%iout,
'(1x,a)') &
252 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
264 if (this%iper /=
kper)
return
267 call this%BndExtType%bnd_rp()
270 if (this%iprpak /= 0)
then
271 call this%write_list()
288 integer(I4B) :: i, node
291 real(DP) :: absdhdxsq
295 real(DP) :: range = 1.d-6
297 real(DP) :: smooth_factor
300 if (this%nbound == 0)
return
303 do i = 1, this%nbound
305 node = this%nodelist(i)
306 if (this%ibound(node) <= 0)
then
313 absdhdxsq = this%slope(i)**
dhalf
314 depth = this%xnew(node) - this%dis%bot(node)
317 call squadratic(depth, range, dydx, smooth_factor)
318 depth = depth * smooth_factor
321 q = -this%qcalc(i, depth, this%unitconv)
325 qeps = -this%qcalc(i, depth + eps, this%unitconv)
328 derv = (qeps - q) / eps
332 this%rhs(i) = -q + derv * this%xnew(node)
342 function qcalc(this, i, depth, unitconv)
result(q)
345 integer(I4B),
intent(in) :: i
346 real(dp),
intent(in) :: depth
347 real(dp),
intent(in) :: unitconv
351 integer(I4B) :: idcxs
355 real(dp) :: conveyance
357 idcxs = this%idcxs(i)
358 width = this%width(i)
359 rough = this%rough(i)
360 slope = this%slope(i)
361 conveyance = this%cxs%get_conveyance(idcxs, width, depth, rough)
362 q = conveyance * slope**
dhalf * unitconv
372 subroutine zdg_fc(this, rhs, ia, idxglo, matrix_sln)
375 real(DP),
dimension(:),
intent(inout) :: rhs
376 integer(I4B),
dimension(:),
intent(in) :: ia
377 integer(I4B),
dimension(:),
intent(in) :: idxglo
385 if (this%imover == 1)
then
386 call this%pakmvrobj%fc()
390 do i = 1, this%nbound
392 rhs(n) = rhs(n) + this%rhs(i)
394 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
398 if (this%imover == 1 .and. this%rhs(i) >
dzero)
then
399 call this%pakmvrobj%accumulate_qformvr(i, this%rhs(i))
415 this%listlabel = trim(this%filtyp)//
' NO.'
416 if (this%dis%ndim == 3)
then
417 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
418 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'ROW'
419 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'COL'
420 elseif (this%dis%ndim == 2)
then
421 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
422 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'CELL2D'
424 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
426 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'FLOW RATE'
427 if (this%inamedbound == 1)
then
428 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
462 call this%obs%StoreObsType(
'zdg', .true., indx)
467 call this%obs%StoreObsType(
'to-mvr', .true., indx)
487 call this%obs%obs_bd_clear()
490 do i = 1, this%obs%npakobs
491 obsrv => this%obs%pakobs(i)%obsrv
492 if (obsrv%BndFound)
then
493 do n = 1, obsrv%indxbnds_count
495 jj = obsrv%indxbnds(n)
496 select case (obsrv%ObsTypeId)
498 if (this%imover == 1)
then
499 v = this%pakmvrobj%get_qtomvr(jj)
507 errmsg =
'Unrecognized observation type: '//trim(obsrv%ObsTypeId)
510 call this%obs%SaveOneSimval(obsrv, v)
513 call this%obs%SaveOneSimval(obsrv,
dnodata)
530 integer(I4B) :: i, nlinks
534 nlinks = this%TsManager%boundtslinks%Count()
537 if (
associated(tslink))
then
538 if (tslink%JCol == 1)
then
556 integer(I4B),
intent(in) :: col
557 integer(I4B),
intent(in) :: row
563 bndval = this%idcxs(row)
565 bndval = this%width(row)
567 bndval = this%slope(row)
569 bndval = this%rough(row)
571 errmsg =
'Programming error. ZDG bound value requested column '&
572 &
'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
subroutine zdg_rp(this)
@ brief SWF read and prepare
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.
integer(i4b), pointer, public kper
current stress period number
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.