25 real(dp),
dimension(:),
pointer :: conc => null()
26 integer(I4B),
dimension(:),
pointer :: icbund => null()
31 integer(I4B),
pointer :: ioutdense => null()
32 integer(I4B),
pointer :: iform => null()
33 integer(I4B),
pointer :: ireadelev => null()
34 integer(I4B),
pointer :: ireadconcbuy => null()
35 integer(I4B),
pointer :: iconcset => null()
36 real(dp),
pointer :: denseref => null()
37 real(dp),
dimension(:),
pointer,
contiguous :: dense => null()
38 real(dp),
dimension(:),
pointer,
contiguous :: concbuy => null()
39 real(dp),
dimension(:),
pointer,
contiguous :: elev => null()
40 integer(I4B),
dimension(:),
pointer :: ibound => null()
42 integer(I4B),
pointer :: nrhospecies => null()
43 real(dp),
dimension(:),
pointer,
contiguous :: drhodc => null()
44 real(dp),
dimension(:),
pointer,
contiguous :: crhoref => null()
45 real(dp),
dimension(:),
pointer,
contiguous :: ctemp => null()
46 character(len=LENMODELNAME),
dimension(:),
allocatable :: cmodelname
47 character(len=LENAUXNAME),
dimension(:),
allocatable :: cauxspeciesname
81 function calcdens(denseref, drhodc, crhoref, conc)
result(dense)
83 real(dp),
intent(in) :: denseref
84 real(dp),
dimension(:),
intent(in) :: drhodc
85 real(dp),
dimension(:),
intent(in) :: crhoref
86 real(dp),
dimension(:),
intent(in) :: conc
90 integer(I4B) :: nrhospec
93 nrhospec =
size(drhodc)
96 dense = dense + drhodc(i) * (conc(i) - crhoref(i))
102 subroutine buy_cr(buyobj, name_model, inunit, iout)
105 character(len=*),
intent(in) :: name_model
106 integer(I4B),
intent(in) :: inunit
107 integer(I4B),
intent(in) :: iout
113 call buyobj%set_names(1, name_model,
'BUY',
'BUY')
116 call buyobj%allocate_scalars()
119 buyobj%inunit = inunit
123 call buyobj%parser%Initialize(buyobj%inunit, buyobj%iout)
135 character(len=*),
parameter :: fmtbuy = &
136 "(1x,/1x,'BUY -- Buoyancy Package, Version 1, 5/16/2018', &
137 &' input read from unit ', i0, //)"
140 if (.not.
present(buy_input))
then
141 write (this%iout, fmtbuy) this%inunit
147 if (.not.
present(buy_input))
then
150 call this%read_options()
153 call this%read_dimensions()
156 call this%set_options(buy_input)
157 this%nrhospecies = buy_input%nrhospecies
161 call this%allocate_arrays(dis%nodes)
163 if (.not.
present(buy_input))
then
166 call this%read_packagedata()
169 call this%set_packagedata(buy_input)
179 integer(I4B),
dimension(:),
pointer :: ibound
183 this%ibound => ibound
186 if (this%npf%ixt3d /= 0)
then
187 call store_error(
'Error in model '//trim(this%name_model)// &
188 '. The XT3D option cannot be used with the BUY Package.')
189 call this%parser%StoreErrorUnit()
193 call this%buy_calcelev()
208 class(
bndtype),
pointer :: packobj
209 real(DP),
intent(in),
dimension(:) :: hnew
212 select case (packobj%filtyp)
216 select type (packobj)
218 call packobj%lak_activate_density()
224 select type (packobj)
226 call packobj%sfr_activate_density()
232 select type (packobj)
234 call packobj%maw_activate_density()
251 character(len=LINELENGTH) :: errmsg
254 character(len=*),
parameter :: fmtc = &
255 "('Buoyancy Package does not have a concentration set &
256 &for species ',i0,'. One or more model names may be specified &
257 &incorrectly in the PACKAGEDATA block or a gwf-gwt exchange may need &
262 do i = 1, this%nrhospecies
263 if (.not.
associated(this%modelconc(i)%conc))
then
264 write (errmsg, fmtc) i
269 call this%parser%StoreErrorUnit()
281 call this%buy_calcdens()
289 integer(I4B) :: kiter
292 if (this%ireadelev == 0)
then
293 if (this%iform == 1 .or. this%iform == 2)
then
294 call this%buy_calcelev()
306 class(
bndtype),
pointer :: packobj
307 real(DP),
intent(in),
dimension(:) :: hnew
310 integer(I4B) :: n, locdense, locelev
311 integer(I4B),
dimension(:),
allocatable :: locconc
315 if (this%iform == 0)
return
320 allocate (locconc(this%nrhospecies))
324 do n = 1, packobj%naux
325 if (packobj%auxname(n) ==
'DENSITY')
then
327 else if (packobj%auxname(n) ==
'ELEVATION')
then
333 do i = 1, this%nrhospecies
335 do j = 1, packobj%naux
336 if (this%cauxspeciesname(i) == packobj%auxname(j))
then
341 if (locconc(i) == 0)
then
349 select case (packobj%filtyp)
353 call buy_cf_ghb(packobj, hnew, this%dense, this%elev, this%denseref, &
354 locelev, locdense, locconc, this%drhodc, this%crhoref, &
355 this%ctemp, this%iform)
359 call buy_cf_riv(packobj, hnew, this%dense, this%elev, this%denseref, &
360 locelev, locdense, locconc, this%drhodc, this%crhoref, &
361 this%ctemp, this%iform)
365 call buy_cf_drn(packobj, hnew, this%dense, this%denseref)
369 call buy_cf_lak(packobj, hnew, this%dense, this%elev, this%denseref, &
370 locdense, locconc, this%drhodc, this%crhoref, &
371 this%ctemp, this%iform)
375 call buy_cf_sfr(packobj, hnew, this%dense, this%elev, this%denseref, &
376 locdense, locconc, this%drhodc, this%crhoref, &
377 this%ctemp, this%iform)
381 call buy_cf_maw(packobj, hnew, this%dense, this%elev, this%denseref, &
382 locdense, locconc, this%drhodc, this%crhoref, &
383 this%ctemp, this%iform)
400 ctemp, auxvar)
result(densebnd)
402 integer(I4B),
intent(in) :: n
403 integer(I4B),
intent(in) :: locdense
404 integer(I4B),
dimension(:),
intent(in) :: locconc
405 real(dp),
intent(in) :: denseref
406 real(dp),
dimension(:),
intent(in) :: drhodc
407 real(dp),
dimension(:),
intent(in) :: crhoref
408 real(dp),
dimension(:),
intent(inout) :: ctemp
409 real(dp),
dimension(:, :),
intent(in) :: auxvar
416 if (locdense > 0)
then
418 densebnd = auxvar(locdense, n)
419 else if (locconc(1) > 0)
then
421 do i = 1,
size(locconc)
423 if (locconc(i) > 0)
then
424 ctemp(i) = auxvar(locconc(i), n)
427 densebnd =
calcdens(denseref, drhodc, crhoref, ctemp)
436 subroutine buy_cf_ghb(packobj, hnew, dense, elev, denseref, locelev, &
437 locdense, locconc, drhodc, crhoref, ctemp, &
442 class(
bndtype),
pointer :: packobj
444 real(DP),
intent(in),
dimension(:) :: hnew
445 real(DP),
intent(in),
dimension(:) :: dense
446 real(DP),
intent(in),
dimension(:) :: elev
447 real(DP),
intent(in) :: denseref
448 integer(I4B),
intent(in) :: locelev
449 integer(I4B),
intent(in) :: locdense
450 integer(I4B),
dimension(:),
intent(in) :: locconc
451 real(DP),
dimension(:),
intent(in) :: drhodc
452 real(DP),
dimension(:),
intent(in) :: crhoref
453 real(DP),
dimension(:),
intent(inout) :: ctemp
454 integer(I4B),
intent(in) :: iform
462 real(DP) :: hcofterm, rhsterm
465 select type (packobj)
467 do n = 1, packobj%nbound
468 node = packobj%nodelist(n)
469 if (packobj%ibound(node) <= 0) cycle
473 drhodc, crhoref, ctemp, packobj%auxvar)
477 if (locelev > 0) elevghb = packobj%auxvar(locelev, n)
480 hghb = packobj%bhead(n)
481 cond = packobj%cond(n)
485 elevghb, elev(node), hghb, hnew(node), &
486 cond, iform, rhsterm, hcofterm)
487 packobj%hcof(n) = packobj%hcof(n) + hcofterm
488 packobj%rhs(n) = packobj%rhs(n) - rhsterm
497 elevghb, elevnode, hghb, hnode, &
498 cond, iform, rhsterm, hcofterm)
500 real(DP),
intent(in) :: denseref
501 real(DP),
intent(in) :: denseghb
502 real(DP),
intent(in) :: densenode
503 real(DP),
intent(in) :: elevghb
504 real(DP),
intent(in) :: elevnode
505 real(DP),
intent(in) :: hghb
506 real(DP),
intent(in) :: hnode
507 real(DP),
intent(in) :: cond
508 integer(I4B),
intent(in) :: iform
509 real(DP),
intent(inout) :: rhsterm
510 real(DP),
intent(inout) :: hcofterm
513 real(DP) :: avgdense, avgelev
516 avgdense =
dhalf * denseghb +
dhalf * densenode
518 t1 = avgdense / denseref -
done
519 t2 = (denseghb - densenode) / denseref
522 hcofterm = -cond * t1
526 hcofterm = hcofterm +
dhalf * cond * t2
530 rhsterm = cond * t1 * hghb
531 rhsterm = rhsterm - cond * t2 * avgelev
532 rhsterm = rhsterm +
dhalf * cond * t2 * hghb
536 rhsterm = rhsterm +
dhalf * cond * t2 * hnode
542 subroutine buy_cf_riv(packobj, hnew, dense, elev, denseref, locelev, &
543 locdense, locconc, drhodc, crhoref, ctemp, &
548 class(
bndtype),
pointer :: packobj
550 real(DP),
intent(in),
dimension(:) :: hnew
551 real(DP),
intent(in),
dimension(:) :: dense
552 real(DP),
intent(in),
dimension(:) :: elev
553 real(DP),
intent(in) :: denseref
554 integer(I4B),
intent(in) :: locelev
555 integer(I4B),
intent(in) :: locdense
556 integer(I4B),
dimension(:),
intent(in) :: locconc
557 real(DP),
dimension(:),
intent(in) :: drhodc
558 real(DP),
dimension(:),
intent(in) :: crhoref
559 real(DP),
dimension(:),
intent(inout) :: ctemp
560 integer(I4B),
intent(in) :: iform
573 select type (packobj)
575 do n = 1, packobj%nbound
576 node = packobj%nodelist(n)
577 if (packobj%ibound(node) <= 0) cycle
581 drhodc, crhoref, ctemp, packobj%auxvar)
585 if (locelev > 0) elevriv = packobj%auxvar(locelev, n)
588 hriv = packobj%stage(n)
589 cond = packobj%cond(n)
590 rbot = packobj%rbot(n)
593 if (hnew(node) > rbot)
then
597 elevriv, elev(node), hriv, hnew(node), &
598 cond, iform, rhsterm, hcofterm)
601 rhsterm = cond * (denseriv / denseref -
done) * (hriv - rbot)
605 packobj%hcof(n) = packobj%hcof(n) + hcofterm
606 packobj%rhs(n) = packobj%rhs(n) - rhsterm
617 class(
bndtype),
pointer :: packobj
619 real(DP),
intent(in),
dimension(:) :: hnew
620 real(DP),
intent(in),
dimension(:) :: dense
621 real(DP),
intent(in) :: denseref
632 select type (packobj)
634 do n = 1, packobj%nbound
635 node = packobj%nodelist(n)
636 if (packobj%ibound(node) <= 0) cycle
638 hbnd = packobj%elev(n)
639 cond = packobj%cond(n)
640 if (hnew(node) > hbnd)
then
641 hcofterm = -cond * (rho / denseref -
done)
642 rhsterm = hcofterm * hbnd
643 packobj%hcof(n) = packobj%hcof(n) + hcofterm
644 packobj%rhs(n) = packobj%rhs(n) + rhsterm
654 subroutine buy_cf_lak(packobj, hnew, dense, elev, denseref, locdense, &
655 locconc, drhodc, crhoref, ctemp, iform)
659 class(
bndtype),
pointer :: packobj
661 real(DP),
intent(in),
dimension(:) :: hnew
662 real(DP),
intent(in),
dimension(:) :: dense
663 real(DP),
intent(in),
dimension(:) :: elev
664 real(DP),
intent(in) :: denseref
665 integer(I4B),
intent(in) :: locdense
666 integer(I4B),
dimension(:),
intent(in) :: locconc
667 real(DP),
dimension(:),
intent(in) :: drhodc
668 real(DP),
dimension(:),
intent(in) :: crhoref
669 real(DP),
dimension(:),
intent(inout) :: ctemp
670 integer(I4B),
intent(in) :: iform
678 select type (packobj)
680 do n = 1, packobj%nbound
683 node = packobj%nodelist(n)
684 if (packobj%ibound(node) <= 0) cycle
688 drhodc, crhoref, ctemp, packobj%auxvar)
691 packobj%denseterms(1, n) = denselak / denseref
694 packobj%denseterms(2, n) = dense(node) / denseref
697 packobj%denseterms(3, n) = elev(node)
707 subroutine buy_cf_sfr(packobj, hnew, dense, elev, denseref, locdense, &
708 locconc, drhodc, crhoref, ctemp, iform)
712 class(
bndtype),
pointer :: packobj
714 real(DP),
intent(in),
dimension(:) :: hnew
715 real(DP),
intent(in),
dimension(:) :: dense
716 real(DP),
intent(in),
dimension(:) :: elev
717 real(DP),
intent(in) :: denseref
718 integer(I4B),
intent(in) :: locdense
719 integer(I4B),
dimension(:),
intent(in) :: locconc
720 real(DP),
dimension(:),
intent(in) :: drhodc
721 real(DP),
dimension(:),
intent(in) :: crhoref
722 real(DP),
dimension(:),
intent(inout) :: ctemp
723 integer(I4B),
intent(in) :: iform
731 select type (packobj)
733 do n = 1, packobj%nbound
736 node = packobj%nodelist(n)
737 if (packobj%ibound(node) <= 0) cycle
741 drhodc, crhoref, ctemp, packobj%auxvar)
744 packobj%denseterms(1, n) = densesfr / denseref
747 packobj%denseterms(2, n) = dense(node) / denseref
750 packobj%denseterms(3, n) = elev(node)
760 subroutine buy_cf_maw(packobj, hnew, dense, elev, denseref, locdense, &
761 locconc, drhodc, crhoref, ctemp, iform)
765 class(
bndtype),
pointer :: packobj
767 real(DP),
intent(in),
dimension(:) :: hnew
768 real(DP),
intent(in),
dimension(:) :: dense
769 real(DP),
intent(in),
dimension(:) :: elev
770 real(DP),
intent(in) :: denseref
771 integer(I4B),
intent(in) :: locdense
772 integer(I4B),
dimension(:),
intent(in) :: locconc
773 real(DP),
dimension(:),
intent(in) :: drhodc
774 real(DP),
dimension(:),
intent(in) :: crhoref
775 real(DP),
dimension(:),
intent(inout) :: ctemp
776 integer(I4B),
intent(in) :: iform
784 select type (packobj)
786 do n = 1, packobj%nbound
789 node = packobj%nodelist(n)
790 if (packobj%ibound(node) <= 0) cycle
794 drhodc, crhoref, ctemp, packobj%auxvar)
797 packobj%denseterms(1, n) = densemaw / denseref
800 packobj%denseterms(2, n) = dense(node) / denseref
803 packobj%denseterms(3, n) = elev(node)
811 subroutine buy_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
814 integer(I4B) :: kiter
816 integer(I4B),
intent(in),
dimension(:) :: idxglo
817 real(DP),
dimension(:),
intent(inout) :: rhs
818 real(DP),
intent(inout),
dimension(:) :: hnew
820 integer(I4B) :: n, m, ipos, idiag
821 real(DP) :: rhsterm, amatnn, amatnm
828 do n = 1, this%dis%nodes
829 if (this%ibound(n) == 0) cycle
830 idiag = this%dis%con%ia(n)
831 do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
832 m = this%dis%con%ja(ipos)
833 if (this%ibound(m) == 0) cycle
834 if (this%iform == 0)
then
835 call this%calcbuy(n, m, ipos, hnew(n), hnew(m), rhsterm)
837 call this%calchhterms(n, m, ipos, hnew(n), hnew(m), rhsterm, &
842 rhs(n) = rhs(n) - rhsterm
843 call matrix_sln%add_value_pos(idxglo(idiag), -amatnn)
844 call matrix_sln%add_value_pos(idxglo(ipos), amatnm)
854 integer(I4B),
intent(in) :: idvfl
856 character(len=1) :: cdatafmp =
' ', editdesc =
' '
857 integer(I4B) :: ibinun
858 integer(I4B) :: iprint
859 integer(I4B) :: nvaluesp
860 integer(I4B) :: nwidthp
864 if (this%ioutdense /= 0)
then
869 if (idvfl == 0) ibinun = 0
872 if (ibinun /= 0)
then
877 if (this%ioutdense /= 0)
then
878 ibinun = this%ioutdense
879 call this%dis%record_array(this%dense, this%iout, iprint, ibinun, &
880 ' DENSITY', cdatafmp, nvaluesp, &
881 nwidthp, editdesc, dinact)
891 real(DP),
intent(in),
dimension(:) :: hnew
892 real(DP),
intent(inout),
dimension(:) :: flowja
893 integer(I4B) :: n, m, ipos
895 real(DP) :: rhsterm, amatnn, amatnm
898 do n = 1, this%dis%nodes
899 do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
900 m = this%dis%con%ja(ipos)
902 if (this%iform == 0)
then
904 call this%calcbuy(n, m, ipos, hnew(n), hnew(m), deltaq)
907 call this%calchhterms(n, m, ipos, hnew(n), hnew(m), rhsterm, &
909 deltaq = amatnm * hnew(m) - amatnn * hnew(n) + rhsterm
911 flowja(ipos) = flowja(ipos) + deltaq
912 flowja(this%dis%con%isym(ipos)) = flowja(this%dis%con%isym(ipos)) - &
925 if (this%inunit > 0)
then
932 deallocate (this%cmodelname)
933 deallocate (this%cauxspeciesname)
934 deallocate (this%modelconc)
948 call this%NumericalPackageType%da()
957 character(len=LINELENGTH) :: errmsg, keyword
959 logical :: isfound, endOfBlock
963 call this%parser%GetBlock(
'DIMENSIONS', isfound, ierr, &
964 supportopenclose=.true.)
968 write (this%iout,
'(/1x,a)')
'Processing BUY DIMENSIONS block'
970 call this%parser%GetNextLine(endofblock)
972 call this%parser%GetStringCaps(keyword)
973 select case (keyword)
975 this%nrhospecies = this%parser%GetInteger()
976 write (this%iout,
'(4x,a,i0)')
'NRHOSPECIES = ', this%nrhospecies
978 write (errmsg,
'(a,a)') &
979 'Unknown BUY dimension: ', trim(keyword)
981 call this%parser%StoreErrorUnit()
984 write (this%iout,
'(1x,a)')
'End of BUY DIMENSIONS block'
986 call store_error(
'Required BUY DIMENSIONS block not found.')
987 call this%parser%StoreErrorUnit()
991 if (this%nrhospecies < 1)
then
992 call store_error(
'NRHOSPECIES must be greater than zero.')
993 call this%parser%StoreErrorUnit()
1003 character(len=LINELENGTH) :: errmsg
1004 character(len=LINELENGTH) :: line
1005 integer(I4B) :: ierr
1006 integer(I4B) :: irhospec
1007 logical :: isfound, endOfBlock
1008 logical :: blockrequired
1009 integer(I4B),
dimension(:),
allocatable :: itemp
1010 character(len=10) :: c10
1011 character(len=16) :: c16
1013 character(len=*),
parameter :: fmterr = &
1014 "('Invalid value for IRHOSPEC (',i0,') detected in BUY Package. &
1015 &IRHOSPEC must be > 0 and <= NRHOSPECIES, and duplicate values &
1016 &are not allowed.')"
1019 allocate (itemp(this%nrhospecies))
1023 blockrequired = .true.
1024 call this%parser%GetBlock(
'PACKAGEDATA', isfound, ierr, &
1025 blockrequired=blockrequired, &
1026 supportopenclose=.true.)
1030 write (this%iout,
'(1x,a)')
'Processing BUY PACKAGEDATA block'
1032 call this%parser%GetNextLine(endofblock)
1033 if (endofblock)
exit
1034 irhospec = this%parser%GetInteger()
1035 if (irhospec < 1 .or. irhospec > this%nrhospecies)
then
1036 write (errmsg, fmterr) irhospec
1039 if (itemp(irhospec) /= 0)
then
1040 write (errmsg, fmterr) irhospec
1044 this%drhodc(irhospec) = this%parser%GetDouble()
1045 this%crhoref(irhospec) = this%parser%GetDouble()
1046 call this%parser%GetStringCaps(this%cmodelname(irhospec))
1047 call this%parser%GetStringCaps(this%cauxspeciesname(irhospec))
1049 write (this%iout,
'(1x,a)')
'End of BUY PACKAGEDATA block'
1051 call store_error(
'Required BUY PACKAGEDATA block not found.')
1052 call this%parser%StoreErrorUnit()
1057 call this%parser%StoreErrorUnit()
1061 write (this%iout,
'(/,a)')
'Summary of species information in BUY Package'
1062 write (this%iout,
'(1a11, 4a17)') &
1063 'SPECIES',
'DRHODC',
'CRHOREF',
'MODEL', &
1065 do irhospec = 1, this%nrhospecies
1066 write (c10,
'(i0)') irhospec
1067 line =
' '//adjustr(c10)
1068 write (c16,
'(g15.6)') this%drhodc(irhospec)
1069 line = trim(line)//
' '//adjustr(c16)
1070 write (c16,
'(g15.6)') this%crhoref(irhospec)
1071 line = trim(line)//
' '//adjustr(c16)
1072 write (c16,
'(a)') this%cmodelname(irhospec)
1073 line = trim(line)//
' '//adjustr(c16)
1074 write (c16,
'(a)') this%cauxspeciesname(irhospec)
1075 line = trim(line)//
' '//adjustr(c16)
1076 write (this%iout,
'(a)') trim(line)
1090 integer(I4B) :: ispec
1092 do ispec = 1, this%nrhospecies
1093 this%drhodc(ispec) = input_data%drhodc(ispec)
1094 this%crhoref(ispec) = input_data%crhoref(ispec)
1095 this%cmodelname(ispec) = input_data%cmodelname(ispec)
1096 this%cauxspeciesname(ispec) = input_data%cauxspeciesname(ispec)
1107 integer(I4B),
intent(in) :: n
1108 integer(I4B),
intent(in) :: m
1109 integer(I4B),
intent(in) :: icon
1110 real(DP),
intent(in) :: hn
1111 real(DP),
intent(in) :: hm
1112 real(DP),
intent(inout) :: buy
1115 real(DP) :: densen, densem, cl1, cl2, avgdense, wt, elevn, elevm, &
1121 densen = this%dense(n)
1122 densem = this%dense(m)
1124 cl1 = this%dis%con%cl1(this%dis%con%jas(icon))
1125 cl2 = this%dis%con%cl2(this%dis%con%jas(icon))
1127 cl1 = this%dis%con%cl2(this%dis%con%jas(icon))
1128 cl2 = this%dis%con%cl1(this%dis%con%jas(icon))
1130 wt = cl1 / (cl1 + cl2)
1131 avgdense = wt * densen + (
done - wt) * densem
1134 if (this%ireadelev == 0)
then
1135 tp = this%dis%top(n)
1136 bt = this%dis%bot(n)
1137 elevn = bt +
dhalf * this%npf%sat(n) * (tp - bt)
1138 tp = this%dis%top(m)
1139 bt = this%dis%bot(m)
1140 elevm = bt +
dhalf * this%npf%sat(m) * (tp - bt)
1142 elevn = this%elev(n)
1143 elevm = this%elev(m)
1146 ihc = this%dis%con%ihc(this%dis%con%jas(icon))
1147 hyn = this%npf%hy_eff(n, m, ihc, ipos=icon)
1148 hym = this%npf%hy_eff(m, n, ihc, ipos=icon)
1151 if (this%dis%con%ihc(this%dis%con%jas(icon)) == 0)
then
1152 cond =
vcond(this%ibound(n), this%ibound(m), &
1153 this%npf%icelltype(n), this%npf%icelltype(m), &
1155 this%npf%ivarcv, this%npf%idewatcv, &
1156 this%npf%condsat(this%dis%con%jas(icon)), hn, hm, &
1158 this%npf%sat(n), this%npf%sat(m), &
1159 this%dis%top(n), this%dis%top(m), &
1160 this%dis%bot(n), this%dis%bot(m), &
1161 this%dis%con%hwva(this%dis%con%jas(icon)))
1163 cond =
hcond(this%ibound(n), this%ibound(m), &
1164 this%npf%icelltype(n), this%npf%icelltype(m), &
1166 this%dis%con%ihc(this%dis%con%jas(icon)), &
1167 this%npf%icellavg, &
1168 this%npf%condsat(this%dis%con%jas(icon)), &
1169 hn, hm, this%npf%sat(n), this%npf%sat(m), &
1171 this%dis%top(n), this%dis%top(m), &
1172 this%dis%bot(n), this%dis%bot(m), &
1173 this%dis%con%cl1(this%dis%con%jas(icon)), &
1174 this%dis%con%cl2(this%dis%con%jas(icon)), &
1175 this%dis%con%hwva(this%dis%con%jas(icon)))
1179 buy = cond * (avgdense - this%denseref) / this%denseref * (elevm - elevn)
1184 subroutine calchhterms(this, n, m, icon, hn, hm, rhsterm, amatnn, amatnm)
1189 integer(I4B),
intent(in) :: n
1190 integer(I4B),
intent(in) :: m
1191 integer(I4B),
intent(in) :: icon
1192 real(DP),
intent(in) :: hn
1193 real(DP),
intent(in) :: hm
1194 real(DP),
intent(inout) :: rhsterm
1195 real(DP),
intent(inout) :: amatnn
1196 real(DP),
intent(inout) :: amatnm
1199 real(DP) :: densen, densem, cl1, cl2, avgdense, wt, elevn, elevm, cond
1200 real(DP) :: rhonormn, rhonormm
1208 densen = this%dense(n)
1209 densem = this%dense(m)
1211 cl1 = this%dis%con%cl1(this%dis%con%jas(icon))
1212 cl2 = this%dis%con%cl2(this%dis%con%jas(icon))
1214 cl1 = this%dis%con%cl2(this%dis%con%jas(icon))
1215 cl2 = this%dis%con%cl1(this%dis%con%jas(icon))
1217 wt = cl1 / (cl1 + cl2)
1218 avgdense = wt * densen + (1.0 - wt) * densem
1221 elevn = this%elev(n)
1222 elevm = this%elev(m)
1223 elevnm = (
done - wt) * elevn + wt * elevm
1225 ihc = this%dis%con%ihc(this%dis%con%jas(icon))
1226 hyn = this%npf%hy_eff(n, m, ihc, ipos=icon)
1227 hym = this%npf%hy_eff(m, n, ihc, ipos=icon)
1231 cond =
vcond(this%ibound(n), this%ibound(m), &
1232 this%npf%icelltype(n), this%npf%icelltype(m), &
1234 this%npf%ivarcv, this%npf%idewatcv, &
1235 this%npf%condsat(this%dis%con%jas(icon)), hn, hm, &
1237 this%npf%sat(n), this%npf%sat(m), &
1238 this%dis%top(n), this%dis%top(m), &
1239 this%dis%bot(n), this%dis%bot(m), &
1240 this%dis%con%hwva(this%dis%con%jas(icon)))
1242 cond =
hcond(this%ibound(n), this%ibound(m), &
1243 this%npf%icelltype(n), this%npf%icelltype(m), &
1245 this%dis%con%ihc(this%dis%con%jas(icon)), &
1246 this%npf%icellavg, &
1247 this%npf%condsat(this%dis%con%jas(icon)), &
1248 hn, hm, this%npf%sat(n), this%npf%sat(m), &
1250 this%dis%top(n), this%dis%top(m), &
1251 this%dis%bot(n), this%dis%bot(m), &
1252 this%dis%con%cl1(this%dis%con%jas(icon)), &
1253 this%dis%con%cl2(this%dis%con%jas(icon)), &
1254 this%dis%con%hwva(this%dis%con%jas(icon)))
1258 rhonormn = densen / this%denseref
1259 rhonormm = densem / this%denseref
1260 rhoterm = wt * rhonormn + (
done - wt) * rhonormm
1261 amatnn = cond * (rhoterm -
done)
1263 rhsterm = -cond * (rhonormm - rhonormn) * elevnm
1264 if (this%iform == 1)
then
1266 hphi = (
done - wt) * hn + wt * hm
1267 rhsterm = rhsterm + cond * hphi * (rhonormm - rhonormn)
1270 amatnn = amatnn - cond * (
done - wt) * (rhonormm - rhonormn)
1271 amatnm = amatnm + cond * wt * (rhonormm - rhonormn)
1286 do n = 1, this%dis%nodes
1287 do i = 1, this%nrhospecies
1288 if (this%modelconc(i)%icbund(n) == 0)
then
1291 this%ctemp(i) = this%modelconc(i)%conc(n)
1294 this%dense(n) =
calcdens(this%denseref, this%drhodc, this%crhoref, &
1306 real(DP) :: tp, bt, frac
1309 do n = 1, this%dis%nodes
1310 tp = this%dis%top(n)
1311 bt = this%dis%bot(n)
1312 frac = this%npf%sat(n)
1313 this%elev(n) = bt +
dhalf * frac * (tp - bt)
1327 call this%NumericalPackageType%allocate_scalars()
1330 call mem_allocate(this%ioutdense,
'IOUTDENSE', this%memoryPath)
1331 call mem_allocate(this%iform,
'IFORM', this%memoryPath)
1332 call mem_allocate(this%ireadelev,
'IREADELEV', this%memoryPath)
1333 call mem_allocate(this%ireadconcbuy,
'IREADCONCBUY', this%memoryPath)
1334 call mem_allocate(this%iconcset,
'ICONCSET', this%memoryPath)
1335 call mem_allocate(this%denseref,
'DENSEREF', this%memoryPath)
1337 call mem_allocate(this%nrhospecies,
'NRHOSPECIES', this%memoryPath)
1343 this%ireadconcbuy = 0
1344 this%denseref = 1000.d0
1346 this%nrhospecies = 0
1358 integer(I4B),
intent(in) :: nodes
1363 call mem_allocate(this%dense, nodes,
'DENSE', this%memoryPath)
1364 call mem_allocate(this%concbuy, 0,
'CONCBUY', this%memoryPath)
1365 call mem_allocate(this%elev, nodes,
'ELEV', this%memoryPath)
1366 call mem_allocate(this%drhodc, this%nrhospecies,
'DRHODC', this%memoryPath)
1367 call mem_allocate(this%crhoref, this%nrhospecies,
'CRHOREF', this%memoryPath)
1368 call mem_allocate(this%ctemp, this%nrhospecies,
'CTEMP', this%memoryPath)
1369 allocate (this%cmodelname(this%nrhospecies))
1370 allocate (this%cauxspeciesname(this%nrhospecies))
1371 allocate (this%modelconc(this%nrhospecies))
1375 this%dense(i) = this%denseref
1376 this%elev(i) =
dzero
1380 do i = 1, this%nrhospecies
1381 this%drhodc(i) =
dzero
1382 this%crhoref(i) =
dzero
1383 this%ctemp(i) =
dzero
1384 this%cmodelname(i) =
''
1385 this%cauxspeciesname(i) =
''
1398 character(len=LINELENGTH) :: errmsg, keyword
1399 character(len=MAXCHARLEN) :: fname
1400 integer(I4B) :: ierr
1401 logical :: isfound, endOfBlock
1403 character(len=*),
parameter :: fmtfileout = &
1404 "(4x, 'BUY ', 1x, a, 1x, ' will be saved to file: ', &
1405 &a, /4x, 'opened on unit: ', I7)"
1408 call this%parser%GetBlock(
'OPTIONS', isfound, ierr, &
1409 supportopenclose=.true., blockrequired=.false.)
1413 write (this%iout,
'(1x,a)')
'Processing BUY OPTIONS block'
1415 call this%parser%GetNextLine(endofblock)
1416 if (endofblock)
exit
1417 call this%parser%GetStringCaps(keyword)
1418 select case (keyword)
1419 case (
'HHFORMULATION_RHS')
1422 write (this%iout,
'(4x,a)') &
1423 'Hydraulic head formulation set to right-hand side'
1425 this%denseref = this%parser%GetDouble()
1426 write (this%iout,
'(4x,a,1pg15.6)') &
1427 'Reference density has been set to: ', &
1429 case (
'DEV_EFH_FORMULATION')
1430 call this%parser%DevOpt()
1433 write (this%iout,
'(4x,a)') &
1434 'Formulation set to equivalent freshwater head'
1436 call this%parser%GetStringCaps(keyword)
1437 if (keyword ==
'FILEOUT')
then
1438 call this%parser%GetString(fname)
1440 call openfile(this%ioutdense, this%iout, fname,
'DATA(BINARY)', &
1442 write (this%iout, fmtfileout) &
1443 'DENSITY', fname, this%ioutdense
1445 errmsg =
'Optional density keyword must be '// &
1446 'followed by FILEOUT'
1450 write (errmsg,
'(a,a)')
'Unknown BUY option: ', &
1453 call this%parser%StoreErrorUnit()
1456 write (this%iout,
'(1x,a)')
'End of BUY OPTIONS block'
1467 this%iform = input_data%iform
1468 this%denseref = input_data%denseref
1472 if (this%iform == 0 .or. this%iform == 1)
then
1486 character(len=LENMODELNAME),
intent(in) :: modelname
1487 real(DP),
dimension(:),
pointer :: conc
1488 integer(I4B),
dimension(:),
pointer :: icbund
1495 do i = 1, this%nrhospecies
1496 if (this%cmodelname(i) == modelname)
then
1497 this%modelconc(i)%conc => conc
1498 this%modelconc(i)%icbund => icbund
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 lenmodelname
maximum length of the model name
real(dp), parameter dhnoflo
real no flow constant
real(dp), parameter dhalf
real constant 1/2
integer(i4b), parameter lenauxname
maximum length of a aux variable
real(dp), parameter dzero
real constant zero
integer(i4b), parameter maxcharlen
maximum length of char string
real(dp), parameter done
real constant 1
real(dp) function calcdens(denseref, drhodc, crhoref, conc)
Generic function to calculate fluid density from concentration.
subroutine buy_cf(this, kiter)
Fill coefficients.
subroutine buy_cf_bnd(this, packobj, hnew)
Fill coefficients.
subroutine read_options(this)
Read package options.
subroutine buy_cf_lak(packobj, hnew, dense, elev, denseref, locdense, locconc, drhodc, crhoref, ctemp, iform)
Pass density information into lak package; density terms are calculated in the lake package as part o...
subroutine set_options(this, input_data)
Sets options as opposed to reading them from a file.
subroutine calchhterms(this, n, m, icon, hn, hm, rhsterm, amatnn, amatnm)
Calculate hydraulic head term for this connection.
subroutine buy_rp(this)
Check for new buy period data.
subroutine, public buy_cr(buyobj, name_model, inunit, iout)
Create a new BUY object.
subroutine buy_cf_maw(packobj, hnew, dense, elev, denseref, locdense, locconc, drhodc, crhoref, ctemp, iform)
Pass density information into maw package; density terms are calculated in the maw package as part of...
subroutine buy_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
Fill coefficients.
subroutine calc_ghb_hcof_rhs_terms(denseref, denseghb, densenode, elevghb, elevnode, hghb, hnode, cond, iform, rhsterm, hcofterm)
Calculate density hcof and rhs terms for ghb conditions.
subroutine allocate_scalars(this)
Allocate scalars used by the package.
subroutine buy_cq(this, hnew, flowja)
Add buy term to flowja.
subroutine buy_calcelev(this)
Calculate cell elevations to use in density flow equations.
subroutine set_concentration_pointer(this, modelname, conc, icbund)
Pass in a gwt model name, concentration array and ibound, and store a pointer to these in the BUY pac...
subroutine buy_cf_riv(packobj, hnew, dense, elev, denseref, locelev, locdense, locconc, drhodc, crhoref, ctemp, iform)
Fill riv coefficients.
subroutine set_packagedata(this, input_data)
Sets package data instead of reading from file.
subroutine buy_df(this, dis, buy_input)
Read options and package data, or set from argument.
subroutine read_packagedata(this)
Read PACKAGEDATA block.
subroutine buy_cf_sfr(packobj, hnew, dense, elev, denseref, locdense, locconc, drhodc, crhoref, ctemp, iform)
Pass density information into sfr package; density terms are calculated in the sfr package as part of...
subroutine buy_ot_dv(this, idvfl)
Save density array to binary file.
real(dp) function get_bnd_density(n, locdense, locconc, denseref, drhodc, crhoref, ctemp, auxvar)
Return the density of the boundary package using one of several different options in the following or...
subroutine buy_cf_ghb(packobj, hnew, dense, elev, denseref, locelev, locdense, locconc, drhodc, crhoref, ctemp, iform)
Fill ghb coefficients.
subroutine read_dimensions(this)
Read the dimensions for this package.
subroutine buy_ar(this, npf, ibound)
Allocate and Read.
subroutine buy_da(this)
Deallocate.
subroutine buy_cf_drn(packobj, hnew, dense, denseref)
Fill drn coefficients.
subroutine buy_ad(this)
Advance.
subroutine allocate_arrays(this, nodes)
Allocate arrays used by the package.
subroutine buy_ar_bnd(this, packobj, hnew)
Buoyancy ar_bnd routine to activate density in packages.
subroutine calcbuy(this, n, m, icon, hn, hm, buy)
Calculate buyancy term for this connection.
subroutine buy_calcdens(this)
calculate fluid density from concentration
This module contains stateless conductance functions.
real(dp) function, public hcond(ibdn, ibdm, ictn, ictm, iupstream, ihc, icellavg, condsat, hn, hm, satn, satm, hkn, hkm, topn, topm, botn, botm, cln, clm, fawidth)
Horizontal conductance between two cells.
real(dp) function, public vcond(ibdn, ibdm, ictn, ictm, inewton, ivarcv, idewatcv, condsat, hn, hm, vkn, vkm, satn, satm, topn, topm, botn, botm, flowarea)
Vertical conductance between two cells.
This module defines variable data types.
This module contains the base numerical package type.
This module contains the SFR package methods.
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.
integer(i4b), pointer, public kstp
current time step number
integer(i4b), pointer, public kper
current stress period number