22 integer(I4B),
pointer :: maxhfb => null()
23 integer(I4B),
pointer :: nhfb => null()
24 integer(I4B),
dimension(:),
pointer,
contiguous :: noden => null()
25 integer(I4B),
dimension(:),
pointer,
contiguous :: nodem => null()
26 integer(I4B),
dimension(:),
pointer,
contiguous :: idxloc => null()
27 real(dp),
dimension(:),
pointer,
contiguous :: hydchr => null()
28 real(dp),
dimension(:),
pointer,
contiguous :: csatsav => null()
29 real(dp),
dimension(:),
pointer,
contiguous :: condsav => null()
32 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound => null()
33 integer(I4B),
dimension(:),
pointer,
contiguous :: icelltype => null()
34 integer(I4B),
dimension(:),
pointer,
contiguous :: ihc => null()
35 integer(I4B),
dimension(:),
pointer,
contiguous :: ia => null()
36 integer(I4B),
dimension(:),
pointer,
contiguous :: ja => null()
37 integer(I4B),
dimension(:),
pointer,
contiguous :: jas => null()
38 integer(I4B),
dimension(:),
pointer,
contiguous :: isym => null()
39 real(dp),
dimension(:),
pointer,
contiguous :: condsat => null()
40 real(dp),
dimension(:),
pointer,
contiguous :: top => null()
41 real(dp),
dimension(:),
pointer,
contiguous :: bot => null()
42 real(dp),
dimension(:),
pointer,
contiguous :: hwva => null()
43 real(dp),
dimension(:),
pointer,
contiguous :: hnew => null()
46 integer(I4B),
pointer :: ivsc => null()
70 subroutine hfb_cr(hfbobj, name_model, input_mempath, inunit, iout)
73 character(len=*),
intent(in) :: name_model
74 character(len=*),
intent(in) :: input_mempath
75 integer(I4B),
intent(in) :: inunit
76 integer(I4B),
intent(in) :: iout
82 call hfbobj%set_names(1, name_model,
'HFB',
'HFB', input_mempath)
85 call hfbobj%allocate_scalars()
88 hfbobj%inunit = inunit
94 subroutine hfb_ar(this, ibound, xt3d, dis, invsc, vsc)
100 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound
103 integer(I4B),
pointer :: invsc
106 character(len=*),
parameter :: fmtheader = &
107 "(1x, /1x, 'HFB -- HORIZONTAL FLOW BARRIER PACKAGE, VERSION 8, ', &
108 &'4/24/2015 INPUT READ FROM MEMPATH: ', a, /)"
111 write (this%iout, fmtheader) this%input_mempath
115 this%ibound => ibound
118 call mem_setptr(this%icelltype,
'ICELLTYPE', &
131 call this%source_options()
132 call this%source_dimensions()
133 call this%allocate_arrays()
141 write (this%iout,
'(/1x,a,a)')
'Viscosity active in ', &
142 trim(this%filtyp)//
' Package calculations: '//trim(adjustl(this%packName))
155 integer(I4B),
pointer :: iper
157 character(len=*),
parameter :: fmtlsp = &
158 &
"(1X,/1X,'REUSING ',A,'S FROM LAST STRESS PERIOD')"
160 call mem_setptr(iper,
'IPER', this%input_mempath)
161 if (iper ==
kper)
then
162 call this%condsat_reset()
163 call this%source_data()
164 call this%condsat_modify()
166 write (this%iout, fmtlsp)
'HFB'
178 subroutine hfb_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
183 integer(I4B) :: kiter
185 integer(I4B),
intent(in),
dimension(:) :: idxglo
186 real(DP),
intent(inout),
dimension(:) :: rhs
187 real(DP),
intent(inout),
dimension(:) :: hnew
189 integer(I4B) :: nodes, nja
190 integer(I4B) :: ihfb, n, m
192 integer(I4B) :: idiag, isymcon
193 integer(I4B) :: ixt3d
194 real(DP) :: cond, condhfb, aterm
195 real(DP) :: fawidth, faheight
196 real(DP) :: topn, topm, botn, botm
197 real(DP) :: viscratio
201 nodes = this%dis%nodes
202 nja = this%dis%con%nja
203 if (
associated(this%xt3d%ixt3d))
then
204 ixt3d = this%xt3d%ixt3d
211 do ihfb = 1, this%nhfb
212 n = min(this%noden(ihfb), this%nodem(ihfb))
213 m = max(this%noden(ihfb), this%nodem(ihfb))
215 if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
217 if (this%ivsc /= 0)
then
218 call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
221 if (this%hydchr(ihfb) >
dzero)
then
222 if (this%inewton == 0)
then
223 ipos = this%idxloc(ihfb)
228 if (this%icelltype(n) == 1)
then
229 if (hnew(n) < topn) topn = hnew(n)
231 if (this%icelltype(m) == 1)
then
232 if (hnew(m) < topm) topm = hnew(m)
234 if (this%ihc(this%jas(ipos)) == 2)
then
235 faheight = min(topn, topm) - max(botn, botm)
237 faheight =
dhalf * ((topn - botn) + (topm - botm))
239 fawidth = this%hwva(this%jas(ipos))
240 condhfb = this%hydchr(ihfb) * viscratio * &
243 condhfb = this%hydchr(ihfb) * viscratio
246 condhfb = this%hydchr(ihfb)
249 call this%xt3d%xt3d_fhfb(kiter, nodes, nja, matrix_sln, idxglo, &
250 rhs, hnew, n, m, condhfb)
256 if (this%inewton == 0)
then
257 do ihfb = 1, this%nhfb
258 ipos = this%idxloc(ihfb)
259 aterm = matrix_sln%get_value_pos(idxglo(ipos))
262 if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
264 if (this%ivsc /= 0)
then
265 call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
268 if (this%icelltype(n) == 1 .or. this%icelltype(m) == 1 .or. &
276 if (this%icelltype(n) == 1)
then
277 if (hnew(n) < topn) topn = hnew(n)
279 if (this%icelltype(m) == 1)
then
280 if (hnew(m) < topm) topm = hnew(m)
282 if (this%ihc(this%jas(ipos)) == 2)
then
283 faheight = min(topn, topm) - max(botn, botm)
285 faheight =
dhalf * ((topn - botn) + (topm - botm))
287 if (this%hydchr(ihfb) >
dzero)
then
288 fawidth = this%hwva(this%jas(ipos))
289 condhfb = this%hydchr(ihfb) * viscratio * &
291 cond = aterm * condhfb / (aterm + condhfb)
293 cond = -aterm * this%hydchr(ihfb)
297 this%condsav(ihfb) = cond
301 call matrix_sln%add_value_pos(idxglo(idiag), aterm - cond)
302 call matrix_sln%set_value_pos(idxglo(ipos), cond)
305 isymcon = this%isym(ipos)
307 call matrix_sln%add_value_pos(idxglo(idiag), aterm - cond)
308 call matrix_sln%set_value_pos(idxglo(isymcon), cond)
327 real(DP),
intent(inout),
dimension(:) :: hnew
328 real(DP),
intent(inout),
dimension(:) :: flowja
330 integer(I4B) :: ihfb, n, m
334 integer(I4B) :: ixt3d
336 real(DP) :: fawidth, faheight
337 real(DP) :: topn, topm, botn, botm
338 real(DP) :: viscratio
343 if (
associated(this%xt3d%ixt3d))
then
344 ixt3d = this%xt3d%ixt3d
351 do ihfb = 1, this%nhfb
352 n = min(this%noden(ihfb), this%nodem(ihfb))
353 m = max(this%noden(ihfb), this%nodem(ihfb))
355 if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
357 if (this%ivsc /= 0)
then
358 call this%vsc%get_visc_ratio(n, m, hnew(n), hnew(m), viscratio)
362 if (this%hydchr(ihfb) >
dzero)
then
363 if (this%inewton == 0)
then
364 ipos = this%idxloc(ihfb)
369 if (this%icelltype(n) == 1)
then
370 if (hnew(n) < topn) topn = hnew(n)
372 if (this%icelltype(m) == 1)
then
373 if (hnew(m) < topm) topm = hnew(m)
375 if (this%ihc(this%jas(ipos)) == 2)
then
376 faheight = min(topn, topm) - max(botn, botm)
378 faheight =
dhalf * ((topn - botn) + (topm - botm))
380 fawidth = this%hwva(this%jas(ipos))
381 condhfb = this%hydchr(ihfb) * viscratio * &
384 condhfb = this%hydchr(ihfb)
387 condhfb = this%hydchr(ihfb)
390 call this%xt3d%xt3d_flowjahfb(n, m, hnew, flowja, condhfb)
396 if (this%inewton == 0)
then
397 do ihfb = 1, this%nhfb
400 if (this%ibound(n) == 0 .or. this%ibound(m) == 0) cycle
401 if (this%icelltype(n) == 1 .or. this%icelltype(m) == 1 .or. &
403 ipos = this%dis%con%getjaindex(n, m)
406 cond = this%condsav(ihfb)
407 qnm = cond * (hnew(m) - hnew(n))
409 ipos = this%dis%con%getjaindex(m, n)
433 if (this%inunit > 0)
then
443 call this%NumericalPackageType%da()
447 this%inewton => null()
448 this%ibound => null()
449 this%icelltype => null()
455 this%condsat => null()
471 call this%NumericalPackageType%allocate_scalars()
474 call mem_allocate(this%maxhfb,
'MAXHFB', this%memoryPath)
496 call mem_allocate(this%noden, this%maxhfb,
'NODEN', this%memoryPath)
497 call mem_allocate(this%nodem, this%maxhfb,
'NODEM', this%memoryPath)
498 call mem_allocate(this%hydchr, this%maxhfb,
'HYDCHR', this%memoryPath)
499 call mem_allocate(this%idxloc, this%maxhfb,
'IDXLOC', this%memoryPath)
500 call mem_allocate(this%csatsav, this%maxhfb,
'CSATSAV', this%memoryPath)
501 call mem_allocate(this%condsav, this%maxhfb,
'CONDSAV', this%memoryPath)
504 do ihfb = 1, this%maxhfb
505 this%idxloc(ihfb) = 0
521 call mem_set_value(this%iprpak,
'PRINT_INPUT', this%input_mempath, &
525 write (this%iout,
'(1x,a)')
'PROCESSING HFB OPTIONS'
526 if (found%print_input)
then
527 write (this%iout,
'(4x,a)') &
528 'THE LIST OF HFBS WILL BE PRINTED.'
530 write (this%iout,
'(1x,a)')
'END OF HFB OPTIONS'
545 call mem_set_value(this%maxhfb,
'MAXBOUND', this%input_mempath, &
549 write (this%iout,
'(/1x,a)')
'PROCESSING HFB DIMENSIONS'
550 write (this%iout,
'(4x,a,i7)')
'MAXHFB = ', this%maxhfb
551 write (this%iout,
'(1x,a)')
'END OF HFB DIMENSIONS'
554 if (this%maxhfb <= 0)
then
556 'MAXHFB must be specified with value greater than zero.'
573 integer(I4B),
dimension(:, :),
pointer,
contiguous :: cellids1, cellids2
574 integer(I4B),
dimension(:),
pointer :: cellid1, cellid2
575 real(DP),
dimension(:),
pointer,
contiguous :: hydchr
576 character(len=LINELENGTH) :: nodenstr, nodemstr
577 integer(I4B),
pointer :: nbound
578 integer(I4B) :: n, nodeu1, nodeu2, noder1, noder2
580 character(len=*),
parameter :: fmthfb =
"(i10, 2a10, 1(1pg15.6))"
583 call mem_setptr(nbound,
'NBOUND', this%input_mempath)
584 call mem_setptr(cellids1,
'CELLID1', this%input_mempath)
585 call mem_setptr(cellids2,
'CELLID2', this%input_mempath)
586 call mem_setptr(hydchr,
'HYDCHR', this%input_mempath)
592 write (this%iout,
'(//,1x,a)')
'READING HFB DATA'
593 if (this%iprpak > 0)
then
594 write (this%iout,
'(3a10, 1a15)')
'HFB NUM',
'CELL1',
'CELL2', &
602 cellid1 => cellids1(:, n)
603 cellid2 => cellids2(:, n)
606 if (this%dis%ndim == 1)
then
609 elseif (this%dis%ndim == 2)
then
610 nodeu1 =
get_node(cellid1(1), 1, cellid1(2), &
611 this%dis%mshape(1), 1, &
613 nodeu2 =
get_node(cellid2(1), 1, cellid2(2), &
614 this%dis%mshape(1), 1, &
617 nodeu1 =
get_node(cellid1(1), cellid1(2), cellid1(3), &
618 this%dis%mshape(1), &
619 this%dis%mshape(2), &
621 nodeu2 =
get_node(cellid2(1), cellid2(2), cellid2(3), &
622 this%dis%mshape(1), &
623 this%dis%mshape(2), &
628 noder1 = this%dis%get_nodenumber(nodeu1, 1)
629 noder2 = this%dis%get_nodenumber(nodeu2, 1)
630 if (noder1 <= 0 .or. &
634 this%noden(n) = noder1
635 this%nodem(n) = noder2
638 this%hydchr(n) = hydchr(n)
641 if (this%iprpak /= 0)
then
642 call this%dis%noder_to_string(this%noden(n), nodenstr)
643 call this%dis%noder_to_string(this%nodem(n), nodemstr)
644 write (this%iout, fmthfb) n, trim(adjustl(nodenstr)), &
645 trim(adjustl(nodemstr)), this%hydchr(n)
651 call store_error(
'Errors encountered in HFB input file.')
656 write (this%iout,
'(3x,i0,a,i0)') this%nhfb, &
657 ' HFBs READ FOR STRESS PERIOD ',
kper
658 write (this%iout,
'(1x,a)')
'END READING HFB DATA'
661 call this%check_data()
674 integer(I4B) :: ihfb, n, m
676 character(len=LINELENGTH) :: nodenstr, nodemstr
679 character(len=*),
parameter :: fmterr =
"(1x, 'HFB no. ',i0, &
680 &' is between two unconnected cells: ', a, ' and ', a)"
681 character(len=*),
parameter :: fmtverr =
"(1x, 'HFB no. ',i0, &
682 &' is between two cells not horizontally connected: ', a, ' and ', a)"
684 do ihfb = 1, this%nhfb
688 do ipos = this%ia(n) + 1, this%ia(n + 1) - 1
689 if (m == this%ja(ipos))
then
691 this%idxloc(ihfb) = ipos
697 if (.not. found)
then
698 call this%dis%noder_to_string(n, nodenstr)
699 call this%dis%noder_to_string(m, nodemstr)
700 write (
errmsg, fmterr) ihfb, trim(adjustl(nodenstr)), &
701 trim(adjustl(nodemstr))
706 ipos = this%idxloc(ihfb)
707 if (this%ihc(this%jas(ipos)) == 0)
then
708 call this%dis%noder_to_string(n, nodenstr)
709 call this%dis%noder_to_string(m, nodemstr)
710 write (
errmsg, fmtverr) ihfb, trim(adjustl(nodenstr)), &
711 trim(adjustl(nodemstr))
732 do ihfb = 1, this%nhfb
733 ipos = this%idxloc(ihfb)
734 this%condsat(this%jas(ipos)) = this%csatsav(ihfb)
750 integer(I4B) :: ihfb, n, m
752 real(DP) :: cond, condhfb
753 real(DP) :: fawidth, faheight
754 real(DP) :: topn, topm, botn, botm
756 do ihfb = 1, this%nhfb
757 ipos = this%idxloc(ihfb)
758 cond = this%condsat(this%jas(ipos))
759 this%csatsav(ihfb) = cond
763 if (this%inewton == 1 .or. &
764 (this%icelltype(n) == 0 .and. this%icelltype(m) == 0))
then
771 if (this%ihc(this%jas(ipos)) == 2)
then
772 faheight = min(topn, topm) - max(botn, botm)
774 faheight =
dhalf * ((topn - botn) + (topm - botm))
776 if (this%hydchr(ihfb) >
dzero)
then
777 fawidth = this%hwva(this%jas(ipos))
778 condhfb = this%hydchr(ihfb) * &
780 cond = cond * condhfb / (cond + condhfb)
782 cond = -cond * this%hydchr(ihfb)
784 this%condsat(this%jas(ipos)) = cond
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
real(dp), parameter done
real constant 1
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
subroutine hfb_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
Fill matrix terms.
subroutine check_data(this)
Check for hfb's between two unconnected cells and write a warning.
subroutine hfb_da(this)
Deallocate memory.
subroutine hfb_cq(this, hnew, flowja)
flowja will automatically include the effects of the hfb for confined and newton cases when xt3d is n...
subroutine, public hfb_cr(hfbobj, name_model, input_mempath, inunit, iout)
Create a new hfb object.
subroutine allocate_scalars(this)
Allocate package scalars.
subroutine condsat_modify(this)
Modify condsat.
subroutine condsat_reset(this)
Reset condsat to its value prior to being modified by hfb's.
subroutine source_data(this)
@ brief source hfb period data
subroutine source_dimensions(this)
@ brief Source hfb input options
subroutine hfb_rp(this)
Check for new HFB stress period data.
subroutine allocate_arrays(this)
Allocate package arrays.
subroutine hfb_ar(this, ibound, xt3d, dis, invsc, vsc)
Allocate and read.
subroutine source_options(this)
@ brief Source hfb input options
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 base numerical package type.
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.
This module contains simulation variables.
character(len=maxcharlen) errmsg
error message string
integer(i4b), pointer, public kper
current stress period number