28 character(len=LENBUDTXT),
dimension(NBDITEMS) ::
budtxt
29 data budtxt/
' STORAGE-AQUEOUS',
' DECAY-AQUEOUS', &
30 ' STORAGE-SORBED',
' DECAY-SORBED'/
53 real(dp),
dimension(:),
pointer,
contiguous :: porosity => null()
54 real(dp),
dimension(:),
pointer,
contiguous :: thetam => null()
55 real(dp),
dimension(:),
pointer,
contiguous :: volfracim => null()
56 real(dp),
dimension(:),
pointer,
contiguous :: ratesto => null()
59 integer(I4B),
pointer :: idcy => null()
60 real(dp),
dimension(:),
pointer,
contiguous :: decay => null()
61 real(dp),
dimension(:),
pointer,
contiguous :: decay_sorbed => null()
62 real(dp),
dimension(:),
pointer,
contiguous :: ratedcy => null()
63 real(dp),
dimension(:),
pointer,
contiguous :: decaylast => null()
64 real(dp),
dimension(:),
pointer,
contiguous :: decayslast => null()
67 integer(I4B),
pointer :: isrb => null()
68 integer(I4B),
pointer :: ioutsorbate => null()
69 real(dp),
dimension(:),
pointer,
contiguous :: bulk_density => null()
70 real(dp),
dimension(:),
pointer,
contiguous :: distcoef => null()
71 real(dp),
dimension(:),
pointer,
contiguous :: sp2 => null()
72 real(dp),
dimension(:),
pointer,
contiguous :: ratesrb => null()
73 real(dp),
dimension(:),
pointer,
contiguous :: ratedcys => null()
74 real(dp),
dimension(:),
pointer,
contiguous :: csrb => null()
77 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound => null()
115 subroutine mst_cr(mstobj, name_model, input_mempath, inunit, iout, fmi)
118 character(len=*),
intent(in) :: name_model
119 character(len=*),
intent(in) :: input_mempath
120 integer(I4B),
intent(in) :: inunit
121 integer(I4B),
intent(in) :: iout
128 call mstobj%set_names(1, name_model,
'MST',
'MST', input_mempath)
131 call mstobj%allocate_scalars()
134 mstobj%inunit = inunit
148 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound
151 character(len=*),
parameter :: fmtmst = &
152 "(1x,/1x,'MST -- MOBILE STORAGE AND TRANSFER PACKAGE, VERSION 1, &
153 &7/29/2020 INPUT READ FROM MEMPATH: ', A, //)"
156 write (this%iout, fmtmst) this%input_mempath
159 call this%source_options()
163 this%ibound => ibound
166 call this%allocate_arrays(dis%nodes)
169 call this%source_data()
176 subroutine mst_fc(this, nodes, cold, nja, matrix_sln, idxglo, cnew, &
181 integer,
intent(in) :: nodes
182 real(DP),
intent(in),
dimension(nodes) :: cold
183 integer(I4B),
intent(in) :: nja
185 integer(I4B),
intent(in),
dimension(nja) :: idxglo
186 real(DP),
intent(inout),
dimension(nodes) :: rhs
187 real(DP),
intent(in),
dimension(nodes) :: cnew
188 integer(I4B),
intent(in) :: kiter
191 call this%mst_fc_sto(nodes, cold, nja, matrix_sln, idxglo, rhs)
195 call this%mst_fc_dcy(nodes, cold, cnew, nja, matrix_sln, idxglo, &
201 call this%mst_fc_srb(nodes, cold, nja, matrix_sln, idxglo, rhs, cnew)
206 call this%mst_fc_dcy_srb(nodes, cold, nja, matrix_sln, idxglo, rhs, &
215 subroutine mst_fc_sto(this, nodes, cold, nja, matrix_sln, idxglo, rhs)
220 integer,
intent(in) :: nodes
221 real(DP),
intent(in),
dimension(nodes) :: cold
222 integer(I4B),
intent(in) :: nja
224 integer(I4B),
intent(in),
dimension(nja) :: idxglo
225 real(DP),
intent(inout),
dimension(nodes) :: rhs
227 integer(I4B) :: n, idiag
229 real(DP) :: hhcof, rrhs
230 real(DP) :: vnew, vold
236 do n = 1, this%dis%nodes
239 if (this%ibound(n) <= 0) cycle
242 vnew = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n)) * &
243 this%fmi%gwfsat(n) * this%thetam(n)
245 if (this%fmi%igwfstrgss /= 0) vold = vold + this%fmi%gwfstrgss(n) *
delt
246 if (this%fmi%igwfstrgsy /= 0) vold = vold + this%fmi%gwfstrgsy(n) *
delt
250 rrhs = -vold * tled * cold(n)
251 idiag = this%dis%con%ia(n)
252 call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
253 rhs(n) = rhs(n) + rrhs
261 subroutine mst_fc_dcy(this, nodes, cold, cnew, nja, matrix_sln, &
267 integer,
intent(in) :: nodes
268 real(DP),
intent(in),
dimension(nodes) :: cold
269 real(DP),
intent(in),
dimension(nodes) :: cnew
270 integer(I4B),
intent(in) :: nja
272 integer(I4B),
intent(in),
dimension(nja) :: idxglo
273 real(DP),
intent(inout),
dimension(nodes) :: rhs
274 integer(I4B),
intent(in) :: kiter
276 integer(I4B) :: n, idiag
277 real(DP) :: hhcof, rrhs
280 real(DP) :: decay_rate
283 do n = 1, this%dis%nodes
286 if (this%ibound(n) <= 0) cycle
289 vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
290 swtpdt = this%fmi%gwfsat(n)
293 idiag = this%dis%con%ia(n)
294 select case (this%idcy)
299 if (cnew(n) >
dzero)
then
300 hhcof = -this%decay(n) * vcell * swtpdt * this%thetam(n)
301 call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
308 kiter, cold(n), cnew(n),
delt)
309 this%decaylast(n) = decay_rate
310 rrhs = decay_rate * vcell * swtpdt * this%thetam(n)
311 rhs(n) = rhs(n) + rrhs
321 subroutine mst_fc_srb(this, nodes, cold, nja, matrix_sln, idxglo, rhs, &
327 integer,
intent(in) :: nodes
328 real(DP),
intent(in),
dimension(nodes) :: cold
329 integer(I4B),
intent(in) :: nja
331 integer(I4B),
intent(in),
dimension(nja) :: idxglo
332 real(DP),
intent(inout),
dimension(nodes) :: rhs
333 real(DP),
intent(in),
dimension(nodes) :: cnew
335 integer(I4B) :: n, idiag
337 real(DP) :: hhcof, rrhs
338 real(DP) :: swt, swtpdt
349 do n = 1, this%dis%nodes
352 if (this%ibound(n) <= 0) cycle
355 vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
356 swtpdt = this%fmi%gwfsat(n)
357 swt = this%fmi%gwfsatold(n,
delt)
358 idiag = this%dis%con%ia(n)
359 const1 = this%distcoef(n)
364 volfracm = this%get_volfracm(n)
365 rhobm = this%bulk_density(n)
366 call mst_srb_term(this%isrb, volfracm, rhobm, vcell, tled, cnew(n), &
367 cold(n), swtpdt, swt, const1, const2, &
368 hcofval=hhcof, rhsval=rrhs)
371 call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
372 rhs(n) = rhs(n) + rrhs
381 subroutine mst_srb_term(isrb, volfracm, rhobm, vcell, tled, cnew, cold, &
382 swnew, swold, const1, const2, rate, hcofval, rhsval)
384 integer(I4B),
intent(in) :: isrb
385 real(DP),
intent(in) :: volfracm
386 real(DP),
intent(in) :: rhobm
387 real(DP),
intent(in) :: vcell
388 real(DP),
intent(in) :: tled
389 real(DP),
intent(in) :: cnew
390 real(DP),
intent(in) :: cold
391 real(DP),
intent(in) :: swnew
392 real(DP),
intent(in) :: swold
393 real(DP),
intent(in) :: const1
394 real(DP),
intent(in) :: const2
395 real(DP),
intent(out),
optional :: rate
396 real(DP),
intent(out),
optional :: hcofval
397 real(DP),
intent(out),
optional :: rhsval
410 term = -volfracm * rhobm * vcell * tled * const1
411 if (
present(hcofval)) hcofval = term * swnew
412 if (
present(rhsval)) rhsval = term * swold * cold
413 if (
present(rate)) rate = term * swnew * cnew - term * swold * cold
417 cavg =
dhalf * (cold + cnew)
434 term = -volfracm * rhobm * vcell * tled
435 cbaravg = (cbarold + cbarnew) *
dhalf
436 swavg = (swnew + swold) *
dhalf
437 if (
present(hcofval))
then
438 hcofval = term * derv * swavg
440 if (
present(rhsval))
then
441 rhsval = term * derv * swavg * cold - term * cbaravg * (swnew - swold)
443 if (
present(rate))
then
444 rate = term * derv * swavg * (cnew - cold) &
445 + term * cbaravg * (swnew - swold)
460 integer,
intent(in) :: nodes
461 real(DP),
intent(in),
dimension(nodes) :: cold
462 integer(I4B),
intent(in) :: nja
464 integer(I4B),
intent(in),
dimension(nja) :: idxglo
465 real(DP),
intent(inout),
dimension(nodes) :: rhs
466 real(DP),
intent(in),
dimension(nodes) :: cnew
467 integer(I4B),
intent(in) :: kiter
469 integer(I4B) :: n, idiag
470 real(DP) :: hhcof, rrhs
478 real(DP) :: decay_rate
483 do n = 1, this%dis%nodes
486 if (this%ibound(n) <= 0) cycle
491 vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
492 swnew = this%fmi%gwfsat(n)
493 distcoef = this%distcoef(n)
494 idiag = this%dis%con%ia(n)
495 volfracm = this%get_volfracm(n)
496 rhobm = this%bulk_density(n)
497 term = this%decay_sorbed(n) * volfracm * rhobm * swnew * vcell
500 select case (this%idcy)
502 select case (this%isrb)
507 if (cnew(n) >
dzero)
then
508 hhcof = -term * distcoef
525 if (distcoef >
dzero)
then
526 select case (this%isrb)
528 csrbold = cold(n) * distcoef
529 csrbnew = cnew(n) * distcoef
539 this%decayslast(n), &
540 kiter, csrbold, csrbnew,
delt)
541 this%decayslast(n) = decay_rate
542 rrhs = decay_rate * volfracm * rhobm * swnew * vcell
547 call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
548 rhs(n) = rhs(n) + rrhs
557 subroutine mst_cq(this, nodes, cnew, cold, flowja)
560 integer(I4B),
intent(in) :: nodes
561 real(DP),
intent(in),
dimension(nodes) :: cnew
562 real(DP),
intent(in),
dimension(nodes) :: cold
563 real(DP),
dimension(:),
contiguous,
intent(inout) :: flowja
566 call this%mst_cq_sto(nodes, cnew, cold, flowja)
570 call this%mst_cq_dcy(nodes, cnew, cold, flowja)
575 call this%mst_cq_srb(nodes, cnew, cold, flowja)
580 call this%mst_cq_dcy_srb(nodes, cnew, cold, flowja)
585 call this%mst_calc_csrb(cnew)
598 integer(I4B),
intent(in) :: nodes
599 real(DP),
intent(in),
dimension(nodes) :: cnew
600 real(DP),
intent(in),
dimension(nodes) :: cold
601 real(DP),
dimension(:),
contiguous,
intent(inout) :: flowja
604 integer(I4B) :: idiag
607 real(DP) :: vnew, vold
608 real(DP) :: hhcof, rrhs
615 this%ratesto(n) =
dzero
618 if (this%ibound(n) <= 0) cycle
621 vnew = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n)) * &
622 this%fmi%gwfsat(n) * this%thetam(n)
624 if (this%fmi%igwfstrgss /= 0) vold = vold + this%fmi%gwfstrgss(n) *
delt
625 if (this%fmi%igwfstrgsy /= 0) vold = vold + this%fmi%gwfstrgsy(n) *
delt
629 rrhs = -vold * tled * cold(n)
630 rate = hhcof * cnew(n) - rrhs
631 this%ratesto(n) = rate
632 idiag = this%dis%con%ia(n)
633 flowja(idiag) = flowja(idiag) + rate
646 integer(I4B),
intent(in) :: nodes
647 real(DP),
intent(in),
dimension(nodes) :: cnew
648 real(DP),
intent(in),
dimension(nodes) :: cold
649 real(DP),
dimension(:),
contiguous,
intent(inout) :: flowja
652 integer(I4B) :: idiag
655 real(DP) :: hhcof, rrhs
657 real(DP) :: decay_rate
665 this%ratedcy(n) =
dzero
666 if (this%ibound(n) <= 0) cycle
669 vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
670 swtpdt = this%fmi%gwfsat(n)
677 if (cnew(n) >
dzero)
then
678 hhcof = -this%decay(n) * vcell * swtpdt * this%thetam(n)
682 0, cold(n), cnew(n),
delt)
683 rrhs = decay_rate * vcell * swtpdt * this%thetam(n)
685 rate = hhcof * cnew(n) - rrhs
686 this%ratedcy(n) = rate
687 idiag = this%dis%con%ia(n)
688 flowja(idiag) = flowja(idiag) + rate
702 integer(I4B),
intent(in) :: nodes
703 real(DP),
intent(in),
dimension(nodes) :: cnew
704 real(DP),
intent(in),
dimension(nodes) :: cold
705 real(DP),
dimension(:),
contiguous,
intent(inout) :: flowja
708 integer(I4B) :: idiag
711 real(DP) :: swt, swtpdt
725 this%ratesrb(n) =
dzero
728 if (this%ibound(n) <= 0) cycle
731 vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
732 swtpdt = this%fmi%gwfsat(n)
733 swt = this%fmi%gwfsatold(n,
delt)
734 volfracm = this%get_volfracm(n)
735 rhobm = this%bulk_density(n)
736 const1 = this%distcoef(n)
741 call mst_srb_term(this%isrb, volfracm, rhobm, vcell, tled, cnew(n), &
742 cold(n), swtpdt, swt, const1, const2, &
744 this%ratesrb(n) = rate
745 idiag = this%dis%con%ia(n)
746 flowja(idiag) = flowja(idiag) + rate
760 integer(I4B),
intent(in) :: nodes
761 real(DP),
intent(in),
dimension(nodes) :: cnew
762 real(DP),
intent(in),
dimension(nodes) :: cold
763 real(DP),
dimension(:),
contiguous,
intent(inout) :: flowja
766 integer(I4B) :: idiag
768 real(DP) :: hhcof, rrhs
778 real(DP) :: decay_rate
785 this%ratedcys(n) =
dzero
788 if (this%ibound(n) <= 0) cycle
793 vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
794 swnew = this%fmi%gwfsat(n)
795 distcoef = this%distcoef(n)
796 volfracm = this%get_volfracm(n)
797 rhobm = this%bulk_density(n)
798 term = this%decay_sorbed(n) * volfracm * rhobm * swnew * vcell
801 select case (this%idcy)
805 select case (this%isrb)
810 if (cnew(n) >
dzero)
then
811 hhcof = -term * distcoef
828 if (distcoef >
dzero)
then
829 select case (this%isrb)
831 csrbold = cold(n) * distcoef
832 csrbnew = cnew(n) * distcoef
841 this%decayslast(n), &
842 0, csrbold, csrbnew,
delt)
843 rrhs = decay_rate * volfracm * rhobm * swnew * vcell
848 rate = hhcof * cnew(n) - rrhs
849 this%ratedcys(n) = rate
850 idiag = this%dis%con%ia(n)
851 flowja(idiag) = flowja(idiag) + rate
861 real(DP),
intent(in),
dimension(:) :: cnew
870 if (this%ibound(n) > 0)
then
871 distcoef = this%distcoef(n)
872 select case (this%isrb)
874 csrb = cnew(n) * distcoef
888 subroutine mst_bd(this, isuppress_output, model_budget)
894 integer(I4B),
intent(in) :: isuppress_output
895 type(
budgettype),
intent(inout) :: model_budget
902 call model_budget%addentry(rin, rout,
delt,
budtxt(1), &
903 isuppress_output, rowlabel=this%packName)
908 call model_budget%addentry(rin, rout,
delt,
budtxt(2), &
909 isuppress_output, rowlabel=this%packName)
915 call model_budget%addentry(rin, rout,
delt,
budtxt(3), &
916 isuppress_output, rowlabel=this%packName)
922 call model_budget%addentry(rin, rout,
delt,
budtxt(4), &
923 isuppress_output, rowlabel=this%packName)
934 integer(I4B),
intent(in) :: icbcfl
935 integer(I4B),
intent(in) :: icbcun
937 integer(I4B) :: ibinun
938 integer(I4B) :: iprint, nvaluesp, nwidthp
939 character(len=1) :: cdatafmp =
' ', editdesc =
' '
943 if (this%ipakcb < 0)
then
945 elseif (this%ipakcb == 0)
then
950 if (icbcfl == 0) ibinun = 0
953 if (ibinun /= 0)
then
958 call this%dis%record_array(this%ratesto, this%iout, iprint, -ibinun, &
959 budtxt(1), cdatafmp, nvaluesp, &
960 nwidthp, editdesc, dinact)
964 call this%dis%record_array(this%ratedcy, this%iout, iprint, -ibinun, &
965 budtxt(2), cdatafmp, nvaluesp, &
966 nwidthp, editdesc, dinact)
970 call this%dis%record_array(this%ratesrb, this%iout, iprint, -ibinun, &
971 budtxt(3), cdatafmp, nvaluesp, &
972 nwidthp, editdesc, dinact)
976 call this%dis%record_array(this%ratedcys, this%iout, iprint, -ibinun, &
977 budtxt(4), cdatafmp, nvaluesp, &
978 nwidthp, editdesc, dinact)
987 integer(I4B),
intent(in) :: idvsave
989 character(len=1) :: cdatafmp =
' ', editdesc =
' '
990 integer(I4B) :: ibinun
991 integer(I4B) :: iprint
992 integer(I4B) :: nvaluesp
993 integer(I4B) :: nwidthp
997 if (this%ioutsorbate /= 0)
then
1002 if (idvsave == 0) ibinun = 0
1005 if (ibinun /= 0)
then
1008 if (this%ioutsorbate /= 0)
then
1009 ibinun = this%ioutsorbate
1010 call this%dis%record_array(this%csrb, this%iout, iprint, ibinun, &
1011 ' SORBATE', cdatafmp, nvaluesp, &
1012 nwidthp, editdesc, dinact)
1029 if (this%inunit > 0)
then
1048 this%ibound => null()
1055 call this%NumericalPackageType%da()
1069 call this%NumericalPackageType%allocate_scalars()
1073 call mem_allocate(this%ioutsorbate,
'IOUTSORBATE', this%memoryPath)
1078 this%ioutsorbate = 0
1092 integer(I4B),
intent(in) :: nodes
1098 call mem_allocate(this%porosity, nodes,
'POROSITY', this%memoryPath)
1099 call mem_allocate(this%thetam, nodes,
'THETAM', this%memoryPath)
1100 call mem_allocate(this%volfracim, nodes,
'VOLFRACIM', this%memoryPath)
1101 call mem_allocate(this%ratesto, nodes,
'RATESTO', this%memoryPath)
1105 call mem_allocate(this%ratedcy, 1,
'RATEDCY', this%memoryPath)
1106 call mem_allocate(this%decay, 1,
'DECAY', this%memoryPath)
1107 call mem_allocate(this%decaylast, 1,
'DECAYLAST', this%memoryPath)
1109 call mem_allocate(this%ratedcy, this%dis%nodes,
'RATEDCY', this%memoryPath)
1110 call mem_allocate(this%decay, nodes,
'DECAY', this%memoryPath)
1111 call mem_allocate(this%decaylast, nodes,
'DECAYLAST', this%memoryPath)
1114 call mem_allocate(this%ratedcys, this%dis%nodes,
'RATEDCYS', &
1116 call mem_allocate(this%decayslast, this%dis%nodes,
'DECAYSLAST', &
1119 call mem_allocate(this%ratedcys, 1,
'RATEDCYS', this%memoryPath)
1120 call mem_allocate(this%decayslast, 1,
'DECAYSLAST', this%memoryPath)
1122 call mem_allocate(this%decay_sorbed, 1,
'DECAY_SORBED', &
1127 call mem_allocate(this%bulk_density, 1,
'BULK_DENSITY', this%memoryPath)
1129 call mem_allocate(this%distcoef, 1,
'DISTCOEF', this%memoryPath)
1130 call mem_allocate(this%ratesrb, 1,
'RATESRB', this%memoryPath)
1131 call mem_allocate(this%csrb, 1,
'CSRB', this%memoryPath)
1133 call mem_allocate(this%bulk_density, nodes,
'BULK_DENSITY', this%memoryPath)
1134 call mem_allocate(this%distcoef, nodes,
'DISTCOEF', this%memoryPath)
1135 call mem_allocate(this%ratesrb, nodes,
'RATESRB', this%memoryPath)
1136 call mem_allocate(this%csrb, nodes,
'CSRB', this%memoryPath)
1140 call mem_allocate(this%sp2, nodes,
'SP2', this%memoryPath)
1146 this%porosity(n) =
dzero
1147 this%thetam(n) =
dzero
1148 this%volfracim(n) =
dzero
1149 this%ratesto(n) =
dzero
1151 do n = 1,
size(this%decay)
1152 this%decay(n) =
dzero
1153 this%ratedcy(n) =
dzero
1154 this%decaylast(n) =
dzero
1156 do n = 1,
size(this%bulk_density)
1157 this%bulk_density(n) =
dzero
1158 this%distcoef(n) =
dzero
1159 this%ratesrb(n) =
dzero
1160 this%csrb(n) =
dzero
1162 do n = 1,
size(this%sp2)
1165 do n = 1,
size(this%ratedcys)
1166 this%ratedcys(n) =
dzero
1167 this%decayslast(n) =
dzero
1186 character(len=LENVARNAME),
dimension(3) :: sorption_method = &
1187 &[character(len=LENVARNAME) ::
'LINEAR',
'FREUNDLICH',
'LANGMUIR']
1188 character(len=linelength) :: fname
1191 call mem_set_value(this%ipakcb,
'SAVE_FLOWS', this%input_mempath, &
1193 call mem_set_value(this%idcy,
'ORDER1_DECAY', this%input_mempath, &
1195 call mem_set_value(this%idcy,
'ORDER0_DECAY', this%input_mempath, &
1197 call mem_set_value(this%isrb,
'SORPTION', this%input_mempath, &
1198 sorption_method, found%sorption)
1199 call mem_set_value(fname,
'SORBATEFILE', this%input_mempath, &
1203 if (found%save_flows) this%ipakcb = -1
1206 if (found%sorption)
then
1207 if (this%isrb == izero)
then
1208 call store_error(
'Unknown sorption type was specified. &
1209 &Sorption must be specified as LINEAR, &
1210 &FREUNDLICH, or LANGMUIR.')
1214 if (found%sorbatefile)
then
1216 call openfile(this%ioutsorbate, this%iout, fname,
'DATA(BINARY)', &
1221 if (this%iout > 0)
then
1222 call this%log_options(found, fname)
1232 character(len=*),
intent(in) :: sorbate_fname
1234 character(len=*),
parameter :: fmtisvflow = &
1235 "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
1236 &WHENEVER ICBCFL IS NOT ZERO.')"
1237 character(len=*),
parameter :: fmtlinear = &
1238 &
"(4x,'LINEAR SORPTION IS ACTIVE. ')"
1239 character(len=*),
parameter :: fmtfreundlich = &
1240 &
"(4x,'FREUNDLICH SORPTION IS ACTIVE. ')"
1241 character(len=*),
parameter :: fmtlangmuir = &
1242 &
"(4x,'LANGMUIR SORPTION IS ACTIVE. ')"
1243 character(len=*),
parameter :: fmtidcy1 = &
1244 &
"(4x,'FIRST-ORDER DECAY IS ACTIVE. ')"
1245 character(len=*),
parameter :: fmtidcy2 = &
1246 &
"(4x,'ZERO-ORDER DECAY IS ACTIVE. ')"
1247 character(len=*),
parameter :: fmtfileout = &
1248 "(4x,'MST ',1x,a,1x,' WILL BE SAVED TO FILE: ',a,/4x,&
1249 &'OPENED ON UNIT: ',I7)"
1251 write (this%iout,
'(1x,a)')
'PROCESSING MOBILE STORAGE AND TRANSFER OPTIONS'
1252 if (found%save_flows)
then
1253 write (this%iout, fmtisvflow)
1255 if (found%order1_decay)
then
1256 write (this%iout, fmtidcy1)
1258 if (found%order0_decay)
then
1259 write (this%iout, fmtidcy2)
1261 if (found%sorption)
then
1262 select case (this%isrb)
1264 write (this%iout, fmtlinear)
1266 write (this%iout, fmtfreundlich)
1268 write (this%iout, fmtlangmuir)
1271 if (found%sorbatefile)
then
1272 write (this%iout, fmtfileout) &
1273 'SORBATE', sorbate_fname, this%ioutsorbate
1275 write (this%iout,
'(1x,a)')
'END OF MOBILE STORAGE AND TRANSFER OPTIONS'
1290 character(len=LINELENGTH) :: errmsg
1292 integer(I4B) :: n, asize
1293 integer(I4B),
dimension(:),
pointer,
contiguous :: map
1297 if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1301 call get_isize(
'BULK_DENSITY', this%input_mempath, asize)
1304 'BULK_DENSITY', this%memoryPath)
1305 call get_isize(
'DISTCOEF', this%input_mempath, asize)
1311 call get_isize(
'DECAY', this%input_mempath, asize)
1313 call mem_reallocate(this%decay, this%dis%nodes,
'DECAY', this%memoryPath)
1315 call get_isize(
'DECAY_SORBED', this%input_mempath, asize)
1318 'DECAY_SORBED', this%memoryPath)
1321 call get_isize(
'SP2', this%input_mempath, asize)
1323 call mem_reallocate(this%sp2, this%dis%nodes,
'SP2', this%memoryPath)
1327 call mem_set_value(this%porosity,
'POROSITY', this%input_mempath, map, &
1329 call mem_set_value(this%decay,
'DECAY', this%input_mempath, map, &
1331 call mem_set_value(this%decay_sorbed,
'DECAY_SORBED', this%input_mempath, &
1332 map, found%decay_sorbed)
1333 call mem_set_value(this%bulk_density,
'BULK_DENSITY', this%input_mempath, &
1334 map, found%bulk_density)
1335 call mem_set_value(this%distcoef,
'DISTCOEF', this%input_mempath, map, &
1337 call mem_set_value(this%sp2,
'SP2', this%input_mempath, map, &
1341 if (this%iout > 0)
then
1342 call this%log_data(found)
1346 if (.not. found%porosity)
then
1347 write (errmsg,
'(a)')
'POROSITY not specified in GRIDDATA block.'
1354 if (.not. found%bulk_density)
then
1355 write (errmsg,
'(a)')
'Sorption is active but BULK_DENSITY &
1356 ¬ specified. BULK_DENSITY must be specified in GRIDDATA block.'
1359 if (.not. found%distcoef)
then
1360 write (errmsg,
'(a)')
'Sorption is active but distribution &
1361 &coefficient not specified. DISTCOEF must be specified in &
1366 if (.not. found%sp2)
then
1367 write (errmsg,
'(a)')
'Freundlich or langmuir sorption is active &
1368 &but SP2 not specified. SP2 must be specified in &
1374 if (found%bulk_density)
then
1375 write (
warnmsg,
'(a)')
'Sorption is not active but &
1376 &BULK_DENSITY was specified. BULK_DENSITY will have no affect on &
1377 &simulation results.'
1379 write (this%iout,
'(1x,a)')
'WARNING. '//
warnmsg
1381 if (found%distcoef)
then
1382 write (
warnmsg,
'(a)')
'Sorption is not active but &
1383 &distribution coefficient was specified. DISTCOEF will have &
1384 &no affect on simulation results.'
1386 write (this%iout,
'(1x,a)')
'WARNING. '//
warnmsg
1389 write (
warnmsg,
'(a)')
'Sorption is not active but &
1390 &SP2 was specified. SP2 will have &
1391 &no affect on simulation results.'
1393 write (this%iout,
'(1x,a)')
'WARNING. '//
warnmsg
1399 if (.not. found%decay)
then
1400 write (errmsg,
'(a)')
'First or zero order decay is &
1401 &active but the first rate coefficient is not specified. DECAY &
1402 &must be specified in GRIDDATA block.'
1405 if (.not. found%decay_sorbed)
then
1411 write (errmsg,
'(a)')
'DECAY_SORBED not provided in GRIDDATA &
1412 &block but decay and sorption are active. Specify DECAY_SORBED &
1413 &in GRIDDATA block.'
1418 if (found%decay)
then
1419 write (
warnmsg,
'(a)')
'First- or zero-order decay &
1420 &is not active but decay was specified. DECAY will &
1421 &have no affect on simulation results.'
1423 write (this%iout,
'(1x,a)')
'WARNING. '//
warnmsg
1425 if (found%decay_sorbed)
then
1426 write (
warnmsg,
'(a)')
'First- or zero-order decay &
1427 &is not active but DECAY_SORBED was specified. &
1428 &DECAY_SORBED will have no affect on simulation results.'
1430 write (this%iout,
'(1x,a)')
'WARNING. '//
warnmsg
1440 do n = 1,
size(this%porosity)
1441 this%thetam(n) = this%porosity(n)
1452 write (this%iout,
'(1x,a)')
'PROCESSING GRIDDATA'
1453 if (found%porosity)
then
1454 write (this%iout,
'(4x,a)')
'MOBILE DOMAIN POROSITY set from input file'
1456 if (found%decay)
then
1457 write (this%iout,
'(4x,a)')
'DECAY RATE set from input file'
1459 if (found%decay_sorbed)
then
1460 write (this%iout,
'(4x,a)')
'DECAY SORBED RATE set from input file'
1462 if (found%bulk_density)
then
1463 write (this%iout,
'(4x,a)')
'BULK DENSITY set from input file'
1465 if (found%distcoef)
then
1466 write (this%iout,
'(4x,a)')
'DISTRIBUTION COEFFICIENT set from input file'
1469 write (this%iout,
'(4x,a)')
'SECOND SORPTION PARAM set from input file'
1471 write (this%iout,
'(1x,a)')
'END PROCESSING GRIDDATA'
1482 real(DP),
dimension(:),
intent(in) :: volfracim
1487 do n = 1, this%dis%nodes
1488 this%volfracim(n) = this%volfracim(n) + volfracim(n)
1493 do n = 1, this%dis%nodes
1494 this%thetam(n) = this%get_volfracm(n) * this%porosity(n)
1505 integer(I4B),
intent(in) :: node
1507 real(dp) :: volfracm
1509 volfracm =
done - this%volfracim(node)
1518 real(dp),
intent(in) :: conc
1519 real(dp),
intent(in) :: kf
1520 real(dp),
intent(in) :: a
1524 if (conc >
dzero)
then
1537 real(dp),
intent(in) :: conc
1538 real(dp),
intent(in) :: kl
1539 real(dp),
intent(in) :: sbar
1543 if (conc >
dzero)
then
1544 cbar = (kl * sbar * conc) / (
done + kl * conc)
1556 real(dp),
intent(in) :: conc
1557 real(dp),
intent(in) :: kf
1558 real(dp),
intent(in) :: a
1562 if (conc >
dzero)
then
1563 derv = kf * a * conc**(a -
done)
1575 real(dp),
intent(in) :: conc
1576 real(dp),
intent(in) :: kl
1577 real(dp),
intent(in) :: sbar
1581 if (conc >
dzero)
then
1582 derv = (kl * sbar) / (
done + kl * conc)**
dtwo
1592 real(dp),
intent(in) :: conc
1593 real(dp),
intent(in) :: kf
1594 real(dp),
intent(in) :: a
1598 if (conc >
dzero)
then
1599 kd = kf * conc**(a -
done)
1609 real(dp),
intent(in) :: conc
1610 real(dp),
intent(in) :: kl
1611 real(dp),
intent(in) :: sbar
1615 if (conc >
dzero)
then
1616 kd = (kl * sbar) / (
done + kl * conc)
1631 cold, cnew, delt)
result(decay_rate)
1633 real(dp),
intent(in) :: decay_rate_usr
1634 real(dp),
intent(in) :: decay_rate_last
1635 integer(I4B),
intent(in) :: kiter
1636 real(dp),
intent(in) :: cold
1637 real(dp),
intent(in) :: cnew
1638 real(dp),
intent(in) :: delt
1640 real(dp) :: decay_rate
1643 if (decay_rate_usr <
dzero)
then
1646 decay_rate = decay_rate_usr
1652 if (kiter == 1)
then
1653 decay_rate = min(decay_rate_usr, cold / delt)
1655 decay_rate = decay_rate_last
1656 if (cnew <
dzero)
then
1657 decay_rate = decay_rate_last + cnew / delt
1658 else if (cnew > cold)
then
1659 decay_rate = decay_rate_last + cnew / delt
1661 decay_rate = min(decay_rate_usr, decay_rate)
1663 decay_rate = max(decay_rate,
dzero)
This module contains the BudgetModule.
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
This module contains simulation constants.
integer(i4b), parameter linelength
maximum length of a standard line
@ mnormal
normal output mode
real(dp), parameter dhnoflo
real no flow constant
integer(i4b), parameter lenvarname
maximum length of a variable name
real(dp), parameter dhalf
real constant 1/2
integer(i4b), parameter izero
integer constant zero
real(dp), parameter dzero
real constant zero
integer(i4b), parameter maxcharlen
maximum length of char string
real(dp), parameter dtwo
real constant 2
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
real(dp), parameter done
real constant 1
– @ brief Mobile Storage and Transfer (MST) Module
subroutine mst_srb_term(isrb, volfracm, rhobm, vcell, tled, cnew, cold, swnew, swold, const1, const2, rate, hcofval, rhsval)
@ brief Calculate sorption terms
subroutine mst_cq_srb(this, nodes, cnew, cold, flowja)
@ brief Calculate sorption terms for package
subroutine mst_calc_csrb(this, cnew)
@ brief Calculate sorbed concentration
subroutine addto_volfracim(this, volfracim)
@ brief Add volfrac values to volfracim
real(dp) function get_zero_order_decay(decay_rate_usr, decay_rate_last, kiter, cold, cnew, delt)
@ brief Calculate zero-order decay rate and constrain if necessary
character(len=lenbudtxt), dimension(nbditems) budtxt
subroutine mst_cq(this, nodes, cnew, cold, flowja)
@ brief Calculate flows for package
@ decay_zero_order
Zeroth-order decay.
@ decay_first_order
First-order decay.
@ sorption_freund
Freundlich sorption between aqueous and solid phases.
@ sorption_lang
Langmuir sorption between aqueous and solid phases.
@ sorption_linear
Linear sorption between aqueous and solid phases.
@ decay_off
Decay (or production) of mass inactive (default)
@ sorption_off
Sorption is inactive (default)
subroutine mst_bd(this, isuppress_output, model_budget)
@ brief Calculate budget terms for package
subroutine allocate_arrays(this, nodes)
@ brief Allocate arrays for package
real(dp) function get_freundlich_derivative(conc, kf, a)
@ brief Calculate sorption derivative using Freundlich
real(dp) function get_langmuir_kd(conc, kl, sbar)
@ brief Get effective Langmuir distribution coefficient
subroutine mst_cq_dcy_srb(this, nodes, cnew, cold, flowja)
@ brief Calculate decay-sorption terms for package
real(dp) function get_freundlich_conc(conc, kf, a)
@ brief Calculate sorption concentration using Freundlich
subroutine, public mst_cr(mstobj, name_model, input_mempath, inunit, iout, fmi)
@ brief Create a new package object
real(dp) function get_volfracm(this, node)
@ brief Return mobile domain volume fraction
subroutine mst_ot_dv(this, idvsave)
Save sorbate concentration array to binary file.
subroutine mst_ot_flow(this, icbcfl, icbcun)
@ brief Output flow terms for package
subroutine mst_ar(this, dis, ibound)
@ brief Allocate and read method for package
subroutine mst_fc_srb(this, nodes, cold, nja, matrix_sln, idxglo, rhs, cnew)
@ brief Fill sorption coefficient method for package
subroutine mst_cq_dcy(this, nodes, cnew, cold, flowja)
@ brief Calculate decay terms for package
subroutine allocate_scalars(this)
@ brief Allocate scalar variables for package
subroutine log_options(this, found, sorbate_fname)
Log user options to list file.
subroutine source_options(this)
@ brief Source options for package
subroutine mst_fc_dcy_srb(this, nodes, cold, nja, matrix_sln, idxglo, rhs, cnew, kiter)
@ brief Fill sorption-decay coefficient method for package
subroutine source_data(this)
@ brief Source data for package
subroutine log_data(this, found)
Log user data to list file.
subroutine mst_da(this)
@ brief Deallocate
subroutine mst_fc_dcy(this, nodes, cold, cnew, nja, matrix_sln, idxglo, rhs, kiter)
@ brief Fill decay coefficient method for package
integer(i4b), parameter nbditems
subroutine mst_fc(this, nodes, cold, nja, matrix_sln, idxglo, cnew, rhs, kiter)
@ brief Fill coefficient method for package
real(dp) function get_freundlich_kd(conc, kf, a)
@ brief Get effective Freundlich distribution coefficient
subroutine mst_fc_sto(this, nodes, cold, nja, matrix_sln, idxglo, rhs)
@ brief Fill storage coefficient method for package
real(dp) function get_langmuir_derivative(conc, kl, sbar)
@ brief Calculate sorption derivative using Langmuir
subroutine mst_cq_sto(this, nodes, cnew, cold, flowja)
@ brief Calculate storage terms for package
real(dp) function get_langmuir_conc(conc, kl, sbar)
@ brief Calculate sorption concentration using Langmuir
This module defines variable data types.
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the base numerical package type.
This module contains simulation methods.
subroutine, public store_warning(msg, substring)
Store warning message.
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
character(len=maxcharlen) warnmsg
warning message string
real(dp), pointer, public delt
length of the current time step
Derived type for the Budget object.
@ brief Mobile storage and transfer