18 character(len=LENFTYPE) ::
ftype =
'GHB'
19 character(len=LENPACKAGENAME) ::
text =
' GHB'
22 real(dp),
dimension(:),
pointer,
contiguous :: bhead => null()
23 real(dp),
dimension(:),
pointer,
contiguous :: cond => null()
46 subroutine ghb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
49 class(
bndtype),
pointer :: packobj
50 integer(I4B),
intent(in) :: id
51 integer(I4B),
intent(in) :: ibcnum
52 integer(I4B),
intent(in) :: inunit
53 integer(I4B),
intent(in) :: iout
54 character(len=*),
intent(in) :: namemodel
55 character(len=*),
intent(in) :: pakname
56 character(len=*),
intent(in) :: mempath
58 type(
ghbtype),
pointer :: ghbobj
65 call packobj%set_names(ibcnum, namemodel, pakname,
ftype, mempath)
69 call packobj%allocate_scalars()
72 call packobj%pack_initialize()
74 packobj%inunit = inunit
77 packobj%ibcnum = ibcnum
90 call this%BndExtType%bnd_da()
105 class(
ghbtype),
intent(inout) :: this
110 call this%BndExtType%source_options()
113 call mem_set_value(this%imover,
'MOVER', this%input_mempath, found%mover)
116 call this%log_ghb_options(found)
125 class(
ghbtype),
intent(inout) :: this
129 write (this%iout,
'(/1x,a)')
'PROCESSING '//trim(adjustl(this%text)) &
132 if (found%mover)
then
133 write (this%iout,
'(4x,A)')
'MOVER OPTION ENABLED'
137 write (this%iout,
'(1x,a)') &
138 'END OF '//trim(adjustl(this%text))//
' OPTIONS'
148 integer(I4B),
dimension(:),
pointer,
contiguous,
optional :: nodelist
149 real(DP),
dimension(:, :),
pointer,
contiguous,
optional :: auxvar
152 call this%BndExtType%allocate_arrays(nodelist, auxvar)
155 call mem_setptr(this%bhead,
'BHEAD', this%input_mempath)
156 call mem_setptr(this%cond,
'COND', this%input_mempath)
159 call mem_checkin(this%bhead,
'BHEAD', this%memoryPath, &
160 'BHEAD', this%input_mempath)
161 call mem_checkin(this%cond,
'COND', this%memoryPath, &
162 'COND', this%input_mempath)
171 class(
ghbtype),
intent(inout) :: this
173 if (this%iper /=
kper)
return
176 call this%BndExtType%bnd_rp()
179 if (this%ivsc == 1)
then
180 call this%ghb_store_user_cond()
184 if (this%iprpak /= 0)
then
185 call this%write_list()
196 class(
ghbtype),
intent(inout) :: this
198 character(len=LINELENGTH) :: errmsg
203 character(len=*),
parameter :: fmtghberr = &
204 "('GHB BOUNDARY (',i0,') HEAD (',f10.3,') IS LESS THAN CELL &
205 &BOTTOM (',f10.3,')')"
206 character(len=*),
parameter :: fmtcondmulterr = &
207 "('GHB BOUNDARY (',i0,') CONDUCTANCE MULTIPLIER (',g10.3,') IS &
209 character(len=*),
parameter :: fmtconderr = &
210 "('GHB BOUNDARY (',i0,') CONDUCTANCE (',g10.3,') IS LESS THAN &
214 do i = 1, this%nbound
215 node = this%nodelist(i)
216 bt = this%dis%bot(node)
218 if (this%bhead(i) < bt .and. this%icelltype(node) /= 0)
then
219 write (errmsg, fmt=fmtghberr) i, this%bhead(i), bt
222 if (this%iauxmultcol > 0)
then
223 if (this%auxvar(this%iauxmultcol, i) < dzero)
then
224 write (errmsg, fmt=fmtcondmulterr) &
225 i, this%auxvar(this%iauxmultcol, i)
229 if (this%cond(i) < dzero)
then
230 write (errmsg, fmt=fmtconderr) i, this%cond(i)
249 integer(I4B) :: i, node
252 if (this%nbound .eq. 0)
return
255 do i = 1, this%nbound
256 node = this%nodelist(i)
257 if (this%ibound(node) .le. 0)
then
262 this%hcof(i) = -this%cond_mult(i)
263 this%rhs(i) = -this%cond_mult(i) * this%bhead(i)
269 subroutine ghb_fc(this, rhs, ia, idxglo, matrix_sln)
272 real(DP),
dimension(:),
intent(inout) :: rhs
273 integer(I4B),
dimension(:),
intent(in) :: ia
274 integer(I4B),
dimension(:),
intent(in) :: idxglo
277 integer(I4B) :: i, n, ipos
278 real(DP) :: cond, bhead, qghb
281 if (this%imover == 1)
then
282 call this%pakmvrobj%fc()
286 do i = 1, this%nbound
288 rhs(n) = rhs(n) + this%rhs(i)
290 call matrix_sln%add_value_pos(idxglo(ipos), this%hcof(i))
294 bhead = this%bhead(i)
295 if (this%imover == 1 .and. this%xnew(n) > bhead)
then
296 cond = this%cond_mult(i)
297 qghb = cond * (this%xnew(n) - bhead)
298 call this%pakmvrobj%accumulate_qformvr(i, qghb)
308 class(
ghbtype),
intent(inout) :: this
311 this%listlabel = trim(this%filtyp)//
' NO.'
312 if (this%dis%ndim == 3)
then
313 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
314 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'ROW'
315 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'COL'
316 elseif (this%dis%ndim == 2)
then
317 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'LAYER'
318 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'CELL2D'
320 write (this%listlabel,
'(a, a7)') trim(this%listlabel),
'NODE'
322 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'STAGE'
323 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'CONDUCTANCE'
324 if (this%inamedbound == 1)
then
325 write (this%listlabel,
'(a, a16)') trim(this%listlabel),
'BOUNDARY NAME'
354 call this%obs%StoreObsType(
'ghb', .true., indx)
359 call this%obs%StoreObsType(
'to-mvr', .true., indx)
367 class(
ghbtype),
intent(inout) :: this
372 do n = 1, this%nbound
373 this%condinput(n) = this%cond_mult(n)
383 class(
ghbtype),
intent(inout) :: this
384 integer(I4B),
intent(in) :: row
388 if (this%iauxmultcol > 0)
then
389 cond = this%cond(row) * this%auxvar(this%iauxmultcol, row)
391 cond = this%cond(row)
401 class(
ghbtype),
intent(inout) :: this
402 integer(I4B),
intent(in) :: col
403 integer(I4B),
intent(in) :: row
409 bndval = this%bhead(row)
411 bndval = this%cond_mult(row)
413 errmsg =
'Programming error. GHB bound value requested column '&
414 &
'outside range of ncolbnd (2).'
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
integer(i4b), parameter lenpackagename
maximum length of the package name
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
real(dp), parameter dzero
real constant zero
subroutine ghb_cf(this)
Formulate the HCOF and RHS terms.
subroutine ghb_allocate_arrays(this, nodelist, auxvar)
Allocate arrays.
subroutine log_ghb_options(this, found)
Log options specific to GhbType.
subroutine ghb_ck(this)
Check ghb boundary condition data.
subroutine ghb_da(this)
Deallocate memory.
subroutine ghb_fc(this, rhs, ia, idxglo, matrix_sln)
Copy rhs and hcof into solution rhs and amat.
subroutine ghb_df_obs(this)
Store observation type supported by GHB package.
character(len=lenpackagename) text
real(dp) function ghb_bound_value(this, col, row)
Return requested boundary value.
real(dp) function cond_mult(this, row)
Apply multiplier to GHB conductance if option AUXMULTCOL is used.
character(len=lenftype) ftype
subroutine ghb_store_user_cond(this)
Store user-specified conductance for GHB boundary type.
logical function ghb_obs_supported(this)
Return true because GHB package supports observations.
subroutine ghb_rp(this)
Read and prepare.
subroutine, public ghb_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
Create a New Ghb Package and point bndobj to the new package.
subroutine ghb_options(this)
Set options specific to GhbType.
subroutine define_listlabel(this)
Define the list heading that is written to iout when PRINT_INPUT option is used.
This module defines variable data types.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
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.
integer(i4b) function, public count_errors()
Return number of errors.
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kper
current stress period number
This class is used to store a single deferred-length character string. It was designed to work in an ...