MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
gwt-ist.f90
Go to the documentation of this file.
1 !> -- @ brief Immobile Storage and Transfer (IST) Module
2 !!
3 !! The GwtIstModule is contains the GwtIstType, which is the
4 !! derived type responsible for adding the effects of an
5 !! immobile domain. In addition to representing transfer
6 !! of mass between the mobile and immobile domain, the IST
7 !! Package also represents the following processes within
8 !! the immobile domain
9 !! 1. Changes in dissolved solute mass
10 !! 2. Decay of dissolved solute mass
11 !! 3. Sorption
12 !! 4. Decay of sorbed solute mass
13 !<
15 
16  use kindmodule, only: dp, i4b
17  use constantsmodule, only: done, dzero, dhalf, lenftype, &
20  use bndmodule, only: bndtype
21  use budgetmodule, only: budgettype
22  use tspfmimodule, only: tspfmitype
29  !
30  implicit none
31  !
32  private
33  public :: ist_create
34  !
35  character(len=LENFTYPE) :: ftype = 'IST'
36  character(len=LENPACKAGENAME) :: text = ' IMMOBILE DOMAIN'
37  integer(I4B), parameter :: nbditems = 5
38  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt
39  data budtxt/' STORAGE-AQUEOUS', ' STORAGE-SORBED', &
40  ' DECAY-AQUEOUS', ' DECAY-SORBED', &
41  ' MOBILE-DOMAIN'/
42 
43  !> @ brief Immobile storage and transfer
44  !!
45  !! Data and methods for handling the effects of an
46  !! immobile domain. Note that there can be as many of these
47  !! domains as necessary. Each immobile domain represents
48  !! changes in immobile solute storage, decay of dissolved
49  !! immobile solute mass, sorption within the immobile domain,
50  !! and decay of immobile domain sorbed mass. The immobile
51  !! domain also includes exchange with the mobile domain.
52  !<
53  type, extends(bndtype) :: gwtisttype
54 
55  type(tspfmitype), pointer :: fmi => null() !< flow model interface
56  type(gwtmsttype), pointer :: mst => null() !< mobile storage and transfer
57  type(budgettype), pointer :: budget => null() !< budget
58  type(outputcontroldatatype), pointer :: ocd => null() !< output control data
59  character(len=LINELENGTH) :: lstfmt !< lst file CIM print format
60  integer(I4B), pointer :: icimout => null() !< unit number for binary cim output
61  integer(I4B), pointer :: ibudgetout => null() !< binary budget output file
62  integer(I4B), pointer :: ibudcsv => null() !< unit number for csv budget output file
63  integer(I4B), pointer :: ioutsorbate => null() !< unit number for sorbate concentration output
64  integer(I4B), pointer :: idcy => null() !< order of decay rate (0:none, 1:first, 2:zero)
65  integer(I4B), pointer :: isrb => null() !< sorption active flag (0:off, 1:on); only linear is supported in ist
66  integer(I4B), pointer :: kiter => null() !< picard iteration counter
67  real(dp), pointer, contiguous :: cim(:) => null() !< concentration for immobile domain
68  real(dp), pointer, contiguous :: cimnew(:) => null() !< immobile concentration at end of current time step
69  real(dp), pointer, contiguous :: cimold(:) => null() !< immobile concentration at end of last time step
70  real(dp), pointer, contiguous :: cimsrb(:) => null() !< sorbate concentration in immobile domain
71  real(dp), pointer, contiguous :: zetaim(:) => null() !< mass transfer rate to immobile domain
72  real(dp), pointer, contiguous :: porosity(:) => null() !< immobile domain porosity defined as volume of immobile voids per volume of immobile domain
73  real(dp), pointer, contiguous :: volfrac(:) => null() !< volume fraction of the immobile domain defined as volume of immobile domain per aquifer volume
74  real(dp), pointer, contiguous :: bulk_density(:) => null() !< bulk density of immobile domain defined as mass of solids in immobile domain per volume of immobile domain
75  real(dp), pointer, contiguous :: distcoef(:) => null() !< distribution coefficient
76  real(dp), pointer, contiguous :: sp2(:) => null() !< second sorption parameter
77  real(dp), pointer, contiguous :: decay(:) => null() !< first or zero order rate constant for liquid
78  real(dp), pointer, contiguous :: decaylast(:) => null() !< decay rate used for last iteration (needed for zero order decay)
79  real(dp), pointer, contiguous :: decayslast(:) => null() !< sorbed decay rate used for last iteration (needed for zero order decay)
80  real(dp), pointer, contiguous :: decay_sorbed(:) => null() !< first or zero order rate constant for sorbed mass
81  real(dp), pointer, contiguous :: strg(:) => null() !< mass transfer rate
82  real(dp) :: budterm(2, nbditems) !< immmobile domain mass summaries
83  class(isothermtype), pointer :: isotherm => null() !< pointer to isotherm object
84 
85  contains
86 
87  procedure :: bnd_ar => ist_ar
88  procedure :: bnd_rp => ist_rp
89  procedure :: bnd_ad => ist_ad
90  procedure :: bnd_fc => ist_fc
91  procedure :: bnd_cq => ist_cq
92  procedure :: bnd_bd => ist_bd
93  procedure :: bnd_ot_model_flows => ist_ot_model_flows
94  procedure :: bnd_ot_dv => ist_ot_dv
97  procedure :: bnd_ot_bdsummary => ist_ot_bdsummary
98  procedure :: bnd_da => ist_da
99  procedure :: allocate_scalars
100  procedure :: read_options
101  procedure :: read_dimensions => ist_read_dimensions
102  procedure :: source_options
103  procedure :: source_data
104  procedure :: log_options
105  procedure :: log_data
106  procedure :: get_thetaim
107  procedure :: ist_calc_csrb
108  procedure, private :: ist_allocate_arrays
109 
110  end type gwtisttype
111 
112 contains
113 
114  !> @ brief Create a new package object
115  !!
116  !! Create a new IST object
117  !!
118  !<
119  subroutine ist_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, &
120  mempath, fmi, mst)
121  ! -- dummy
122  class(bndtype), pointer :: packobj !< BndType pointer that will point to new IST Package
123  integer(I4B), intent(in) :: id !< name of the model
124  integer(I4B), intent(in) :: ibcnum !< consecutive package number
125  integer(I4B), intent(in) :: inunit !< unit number of package input file
126  integer(I4B), intent(in) :: iout !< unit number of model listing file
127  character(len=*), intent(in) :: namemodel !< name of the model
128  character(len=*), intent(in) :: pakname !< name of the package
129  character(len=*), intent(in) :: mempath
130  type(tspfmitype), pointer, intent(in) :: fmi
131  type(gwtmsttype), pointer, intent(in) :: mst
132  ! -- local
133  type(gwtisttype), pointer :: istobj
134  !
135  ! -- allocate the object and assign values to object variables
136  allocate (istobj)
137  packobj => istobj
138  !
139  ! -- create name and memory path
140  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
141  packobj%text = text
142  !
143  ! -- allocate scalars
144  call packobj%allocate_scalars()
145  !
146  ! -- initialize package
147  call packobj%pack_initialize()
148  !
149  ! -- store values
150  packobj%inunit = inunit
151  packobj%iout = iout
152  packobj%id = id
153  packobj%ibcnum = ibcnum
154  packobj%ncolbnd = 1
155  packobj%iscloc = 1
156  !
157  ! -- Point IST specific variables
158  istobj%fmi => fmi
159  istobj%mst => mst
160  end subroutine ist_create
161 
162  !> @ brief Allocate and read method for package
163  !!
164  !! Method to allocate and read static data for the package.
165  !!
166  !<
167  subroutine ist_ar(this)
168  ! -- modules
170  use budgetmodule, only: budget_cr
171  ! -- dummy
172  class(gwtisttype), intent(inout) :: this !< GwtIstType object
173  ! -- local
174  integer(I4B) :: n
175  !
176  ! -- Allocate arrays
177  call this%ist_allocate_arrays()
178  !
179  ! -- Now that arrays are allocated, check in the cimnew array to
180  ! the output control manager for subsequent printing/saving
181  call this%ocd%init_dbl('CIM', this%cimnew, this%dis, 'PRINT LAST ', &
182  'COLUMNS 10 WIDTH 11 DIGITS 4 GENERAL ', &
183  this%iout, dhnoflo)
184 
185  ! -- apply user override if provided
186  if (this%lstfmt /= '') then
187  call this%ocd%set_prnfmt(trim(this%lstfmt)//" ", 0)
188  end if
189  !
190  ! -- source the data block
191  call this%source_data()
192  !
193  ! -- set cimnew to the cim start values read from input
194  do n = 1, this%dis%nodes
195  this%cimnew(n) = this%cim(n)
196  end do
197  !
198  ! -- add volfrac to the volfracim accumulator in mst package
199  call this%mst%addto_volfracim(this%volfrac)
200  !
201  ! -- setup the immobile domain budget
202  call budget_cr(this%budget, this%memoryPath)
203  call this%budget%budget_df(nbditems, 'MASS', 'M', bdzone=this%packName)
204  call this%budget%set_ibudcsv(this%ibudcsv)
205  !
206  ! -- Create isotherm object if sorption is active
207  this%isotherm => create_isotherm(this%isrb, this%distcoef, this%sp2)
208  !
209  ! -- Perform a check to ensure that sorption and decay are set
210  ! consistently between the MST and IST packages.
211  if (this%idcy /= this%mst%idcy) then
212  call store_error('DECAY must be activated consistently between the &
213  &MST and IST Packages. Activate or deactivate DECAY for both &
214  &Packages.')
215  end if
216  if (this%isrb /= this%mst%isrb) then
217  call store_error('SORPTION must be activated consistently between the &
218  &MST and IST Packages. Activate or deactivate SORPTION for both &
219  &Packages. If activated, the same type of sorption (LINEAR, &
220  &FREUNDLICH, or LANGMUIR) must be specified in the options block of &
221  &both the MST and IST Packages.')
222  end if
223  if (count_errors() > 0) then
224  call store_error_filename(this%input_fname)
225  end if
226  end subroutine ist_ar
227 
228  !> @ brief Read and prepare method for package
229  !!
230  !! Method to read and prepare package data
231  !!
232  !<
233  subroutine ist_rp(this)
234  ! -- dummy
235  class(gwtisttype), intent(inout) :: this !< GwtIstType object
236  ! -- local
237  ! -- format
238  end subroutine ist_rp
239 
240  !> @ brief Advance the ist package
241  !!
242  !! Advance the IST Package and handle the adaptive time stepping
243  !! feature by copying from new to old or old to new accordingly
244  !!
245  !<
246  subroutine ist_ad(this)
247  ! -- modules
249  ! -- dummy variables
250  class(gwtisttype) :: this !< GwtIstType object
251  ! -- local variables
252  integer(I4B) :: n
253  !
254  ! -- Call parent advance
255  call this%BndType%bnd_ad()
256  !
257  ! -- set independent kiter counter to zero
258  this%kiter = 0
259  !
260  ! -- copy cimnew into cimold or vice versa if this is a repeat of
261  ! a failed time step
262  if (ifailedstepretry == 0) then
263  do n = 1, this%dis%nodes
264  this%cimold(n) = this%cimnew(n)
265  end do
266  else
267  do n = 1, this%dis%nodes
268  this%cimnew(n) = this%cimold(n)
269  end do
270  end if
271  end subroutine ist_ad
272 
273  !> @ brief Fill coefficient method for package
274  !<
275  subroutine ist_fc(this, rhs, ia, idxglo, matrix_sln)
276  ! modules
277  use tdismodule, only: delt
278  ! dummy
279  class(gwtisttype) :: this !< GwtIstType object
280  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
281  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
282  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
283  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
284  ! local
285  integer(I4B) :: n, idiag
286  real(DP) :: tled
287  real(DP) :: hhcof, rrhs
288  real(DP) :: swt, swtpdt
289  real(DP) :: vcell
290  real(DP) :: thetaim
291  real(DP) :: zetaim
292  real(DP) :: kdnew
293  real(DP) :: kdold
294  real(DP) :: volfracim
295  real(DP) :: rhobim
296  real(DP) :: lambda1im
297  real(DP) :: lambda2im
298  real(DP) :: gamma1im
299  real(DP) :: gamma2im
300  real(DP) :: cimold
301  real(DP) :: f
302  real(DP) :: cimsrbold
303  real(DP) :: cimsrbnew
304  real(DP), dimension(9) :: ddterm
305 
306  ! set variables
307  tled = done / delt
308  this%kiter = this%kiter + 1
309 
310  ! loop through each node and calculate immobile domain contribution
311  ! to hcof and rhs
312  do n = 1, this%dis%nodes
313 
314  ! skip if transport inactive
315  if (this%ibound(n) <= 0) cycle
316 
317  ! calculate new and old water volumes
318  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
319  swtpdt = this%fmi%gwfsat(n)
320  swt = this%fmi%gwfsatold(n, delt)
321  thetaim = this%get_thetaim(n)
322  idiag = ia(n)
323 
324  ! set exchange coefficient
325  zetaim = this%zetaim(n)
326 
327  ! Add dual domain mass transfer contributions to rhs and hcof
328  ! dcimsrbdc = DZERO
329  kdnew = dzero
330  kdold = dzero
331  volfracim = dzero
332  rhobim = dzero
333  lambda1im = dzero
334  lambda2im = dzero
335  gamma1im = dzero
336  gamma2im = dzero
337 
338  ! set variables for decay of aqueous solute
339  if (this%idcy == 1) lambda1im = this%decay(n)
340  if (this%idcy == 2) then
341  gamma1im = get_zero_order_decay(this%decay(n), this%decaylast(n), &
342  this%kiter, this%cimold(n), &
343  this%cimnew(n), delt)
344  this%decaylast(n) = gamma1im
345  end if
346 
347  ! setup sorption variables
348  if (this%isrb > 0) then
349 
350  ! initialize sorption variables
351  volfracim = this%volfrac(n)
352  rhobim = this%bulk_density(n)
353 
354  ! set isotherm dependent sorption variables
355  cimsrbnew = this%isotherm%value(this%cimnew, n)
356  cimsrbold = this%isotherm%value(this%cimold, n)
357  select case (this%isrb)
358  case (sorption_linear)
359  kdnew = this%distcoef(n)
360  kdold = this%distcoef(n)
361  case (sorption_freund)
362  kdnew = get_freundlich_kd(this%cimnew(n), this%distcoef(n), &
363  this%sp2(n))
364  kdold = get_freundlich_kd(this%cimold(n), this%distcoef(n), &
365  this%sp2(n))
366  case (sorption_lang)
367  kdnew = get_langmuir_kd(this%cimnew(n), this%distcoef(n), &
368  this%sp2(n))
369  kdold = get_langmuir_kd(this%cimold(n), this%distcoef(n), &
370  this%sp2(n))
371  end select
372 
373  ! set decay of sorbed solute parameters
374  if (this%idcy == 1) then
375  lambda2im = this%decay_sorbed(n)
376  else if (this%idcy == 2) then
377  gamma2im = get_zero_order_decay(this%decay_sorbed(n), &
378  this%decayslast(n), &
379  this%kiter, cimsrbold, &
380  cimsrbnew, delt)
381  this%decayslast(n) = gamma2im
382  end if
383  end if
384 
385  ! calculate dual domain terms and get the hcof and rhs contributions
386  call get_ddterm(thetaim, vcell, delt, swtpdt, &
387  volfracim, rhobim, kdnew, kdold, lambda1im, lambda2im, &
388  gamma1im, gamma2im, zetaim, ddterm, f)
389  cimold = this%cimold(n)
390  call get_hcofrhs(ddterm, f, cimold, hhcof, rrhs)
391 
392  ! update solution accumulators
393  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
394  rhs(n) = rhs(n) + rrhs
395 
396  end do
397 
398  end subroutine ist_fc
399 
400  !> @ brief Calculate package flows.
401  !<
402  subroutine ist_cq(this, x, flowja, iadv)
403  ! modules
404  use tdismodule, only: delt
405  use constantsmodule, only: dzero
406  ! dummy
407  class(gwtisttype), intent(inout) :: this !< GwtIstType object
408  real(DP), dimension(:), intent(in) :: x !< current dependent-variable value
409  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
410  integer(I4B), optional, intent(in) :: iadv !< flag that indicates if this is an advance package
411  ! local
412  integer(I4B) :: idiag
413  integer(I4B) :: n
414  real(DP) :: rate
415  real(DP) :: swt, swtpdt
416  real(DP) :: hhcof, rrhs
417  real(DP) :: vcell
418  real(DP) :: thetaim
419  real(DP) :: zetaim
420  real(DP) :: kdnew
421  real(DP) :: kdold
422  real(DP) :: volfracim
423  real(DP) :: rhobim
424  real(DP) :: lambda1im
425  real(DP) :: lambda2im
426  real(DP) :: gamma1im
427  real(DP) :: gamma2im
428  real(DP) :: cimnew
429  real(DP) :: cimold
430  real(DP) :: f
431  real(DP) :: cimsrbold
432  real(DP) :: cimsrbnew
433  real(DP), dimension(9) :: ddterm
434  ! formats
435 
436  ! initialize
437  this%budterm(:, :) = dzero
438 
439  ! Calculate immobile domain transfer rate
440  do n = 1, this%dis%nodes
441 
442  ! skip if transport inactive
443  rate = dzero
444  cimnew = dzero
445  if (this%ibound(n) > 0) then
446 
447  ! calculate new and old water volumes
448  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
449  swtpdt = this%fmi%gwfsat(n)
450  swt = this%fmi%gwfsatold(n, delt)
451  thetaim = this%get_thetaim(n)
452 
453  ! set exchange coefficient
454  zetaim = this%zetaim(n)
455 
456  ! Calculate exchange with immobile domain
457  rate = dzero
458  hhcof = dzero
459  rrhs = dzero
460  kdnew = dzero
461  kdold = dzero
462  volfracim = dzero
463  rhobim = dzero
464  lambda1im = dzero
465  lambda2im = dzero
466  gamma1im = dzero
467  gamma2im = dzero
468  if (this%idcy == 1) lambda1im = this%decay(n)
469  if (this%idcy == 2) then
470  gamma1im = get_zero_order_decay(this%decay(n), this%decaylast(n), 0, &
471  this%cimold(n), this%cimnew(n), delt)
472  end if
473 
474  ! setup sorption variables
475  if (this%isrb > 0) then
476 
477  ! initialize sorption variables
478  volfracim = this%volfrac(n)
479  rhobim = this%bulk_density(n)
480 
481  ! set isotherm dependent sorption variables
482  cimsrbnew = this%isotherm%value(this%cimnew, n)
483  cimsrbold = this%isotherm%value(this%cimold, n)
484  select case (this%isrb)
485  case (sorption_linear)
486  kdnew = this%distcoef(n)
487  kdold = this%distcoef(n)
488  case (sorption_freund)
489  kdnew = get_freundlich_kd(this%cimnew(n), this%distcoef(n), &
490  this%sp2(n))
491  kdold = get_freundlich_kd(this%cimold(n), this%distcoef(n), &
492  this%sp2(n))
493  case (sorption_lang)
494  kdnew = get_langmuir_kd(this%cimnew(n), this%distcoef(n), &
495  this%sp2(n))
496  kdold = get_langmuir_kd(this%cimold(n), this%distcoef(n), &
497  this%sp2(n))
498  end select
499 
500  ! set decay of sorbed solute parameters
501  if (this%idcy == 1) then
502  lambda2im = this%decay_sorbed(n)
503  else if (this%idcy == 2) then
504  gamma2im = get_zero_order_decay(this%decay_sorbed(n), &
505  this%decayslast(n), &
506  0, cimsrbold, &
507  cimsrbnew, delt)
508  end if
509  end if
510 
511  ! calculate the terms and then get the hcof and rhs contributions
512  call get_ddterm(thetaim, vcell, delt, swtpdt, &
513  volfracim, rhobim, kdnew, kdold, lambda1im, lambda2im, &
514  gamma1im, gamma2im, zetaim, ddterm, f)
515  cimold = this%cimold(n)
516  call get_hcofrhs(ddterm, f, cimold, hhcof, rrhs)
517 
518  ! calculate rate from hcof and rhs
519  rate = hhcof * x(n) - rrhs
520 
521  ! calculate immobile domain concentration
522  cimnew = get_ddconc(ddterm, f, cimold, x(n))
523 
524  ! accumulate the budget terms
525  call accumulate_budterm(this%budterm, ddterm, cimnew, cimold, x(n), &
526  this%idcy)
527  end if
528 
529  ! store rate and add to flowja
530  this%strg(n) = rate
531  idiag = this%dis%con%ia(n)
532  flowja(idiag) = flowja(idiag) + rate
533 
534  ! store immobile domain concentration
535  this%cimnew(n) = cimnew
536 
537  end do
538 
539  ! calculate csrb
540  if (this%isrb /= 0) then
541  call this%ist_calc_csrb(this%cimnew)
542  end if
543 
544  end subroutine ist_cq
545 
546  !> @ brief Calculate immobile sorbed concentration
547  !<
548  subroutine ist_calc_csrb(this, cim)
549  ! -- dummy
550  class(gwtisttype) :: this !< GwtMstType object
551  real(DP), intent(in), dimension(:) :: cim !< immobile domain aqueous concentration at end of this time step
552  ! -- local
553  integer(I4B) :: n
554  real(DP) :: csrb
555 
556  ! Calculate sorbed concentration
557  do n = 1, size(cim)
558  csrb = dzero
559  if (this%ibound(n) > 0 .and. this%isrb /= sorption_off) then
560  csrb = this%isotherm%value(cim, n)
561  end if
562  this%cimsrb(n) = csrb
563  end do
564 
565  end subroutine ist_calc_csrb
566 
567  !> @ brief Add package flows to model budget.
568  !!
569  !! Add the flow between IST package and the model (ratin and ratout) to the
570  !! model budget.
571  !!
572  !<
573  subroutine ist_bd(this, model_budget)
574  ! -- modules
575  use tdismodule, only: delt
577  ! -- dummy
578  class(gwtisttype) :: this !< GwtIstType object
579  type(budgettype), intent(inout) :: model_budget !< model budget object
580  ! -- local
581  real(DP) :: ratin
582  real(DP) :: ratout
583  integer(I4B) :: isuppress_output
584  isuppress_output = 0
585  call rate_accumulator(this%strg(:), ratin, ratout)
586  call model_budget%addentry(ratin, ratout, delt, this%text, &
587  isuppress_output, this%packName)
588  end subroutine ist_bd
589 
590  !> @ brief Output model flow terms.
591  !!
592  !! Output flow terms between the IST package and model to a binary file and/or
593  !! print flows to the model listing file.
594  !!
595  !<
596  subroutine ist_ot_model_flows(this, icbcfl, ibudfl, icbcun, imap)
597  ! -- modules
598  use constantsmodule, only: dzero
599  ! -- dummy
600  class(gwtisttype) :: this !< GwtIstType object
601  integer(I4B), intent(in) :: icbcfl !< flag for cell-by-cell output
602  integer(I4B), intent(in) :: ibudfl !< flag indication if cell-by-cell data should be saved
603  integer(I4B), intent(in) :: icbcun !< unit number for cell-by-cell output
604  integer(I4B), dimension(:), optional, intent(in) :: imap !< mapping vector
605  ! -- local
606  integer(I4B) :: n
607  integer(I4B) :: ibinun
608  integer(I4B) :: nbound
609  integer(I4B) :: naux
610  real(DP) :: rate
611  real(DP), dimension(0) :: auxrow
612  !
613  ! -- Set unit number for binary output
614  if (this%ipakcb < 0) then
615  ibinun = icbcun
616  elseif (this%ipakcb == 0) then
617  ibinun = 0
618  else
619  ibinun = this%ipakcb
620  end if
621  if (icbcfl == 0) ibinun = 0
622  !
623  ! -- Record the storage rate if requested
624  !
625  ! -- If cell-by-cell flows will be saved as a list, write header.
626  if (ibinun /= 0) then
627  nbound = this%dis%nodes
628  naux = 0
629  call this%dis%record_srcdst_list_header(this%text, this%name_model, &
630  this%name_model, this%name_model, &
631  this%packName, naux, this%auxname, &
632  ibinun, nbound, this%iout)
633  end if
634  !
635  ! -- Calculate immobile domain rhs and hcof
636  do n = 1, this%dis%nodes
637  !
638  ! -- skip if transport inactive
639  rate = dzero
640  if (this%ibound(n) > 0) then
641  !
642  ! -- set rate from this%strg
643  rate = this%strg(n)
644  end if
645  !
646  ! -- If saving cell-by-cell flows in list, write flow
647  if (ibinun /= 0) then
648  call this%dis%record_mf6_list_entry(ibinun, n, n, rate, &
649  naux, auxrow, &
650  olconv=.true., &
651  olconv2=.true.)
652  end if
653  !
654  end do
655  end subroutine ist_ot_model_flows
656 
657  !> @ brief Output dependent variables.
658  !<
659  subroutine ist_ot_dv(this, idvsave, idvprint)
660  ! dummy variables
661  class(gwtisttype) :: this !< BndType object
662  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
663  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
664 
665  call this%output_immobile_concentration(idvsave, idvprint)
666  call this%output_immobile_sorbate_concentration(idvsave, idvprint)
667 
668  end subroutine ist_ot_dv
669 
670  !> @ brief Output immobile domain aqueous concentration.
671  !<
672  subroutine output_immobile_concentration(this, idvsave, idvprint)
673  ! modules
674  use tdismodule, only: kstp, endofperiod
675  ! dummy variables
676  class(gwtisttype) :: this !< BndType object
677  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
678  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
679  ! local
680  integer(I4B) :: ipflg
681  integer(I4B) :: ibinun
682  !
683  ! Save cim to a binary file. ibinun is a flag where 1 indicates that
684  ! cim should be written to a binary file if a binary file is open for it.
685  ipflg = 0
686  ibinun = 1
687  if (idvsave == 0) ibinun = 0
688  if (ibinun /= 0) then
689  call this%ocd%ocd_ot(ipflg, kstp, endofperiod, this%iout, &
690  iprint_opt=0, isav_opt=ibinun)
691  end if
692  !
693  ! Print immobile domain concentrations to listing file
694  if (idvprint /= 0) then
695  call this%ocd%ocd_ot(ipflg, kstp, endofperiod, this%iout, &
696  iprint_opt=idvprint, isav_opt=0)
697  end if
698 
699  end subroutine output_immobile_concentration
700 
701  !> @ brief Output immobile domain sorbate concentration.
702  !<
703  subroutine output_immobile_sorbate_concentration(this, idvsave, idvprint)
704  ! modules
705  ! dummy
706  class(gwtisttype) :: this !< BndType object
707  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
708  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
709  ! local
710  character(len=1) :: cdatafmp = ' ', editdesc = ' '
711  integer(I4B) :: ibinun
712  integer(I4B) :: iprint, nvaluesp, nwidthp
713  real(DP) :: dinact
714 
715  ! Save cimsrb to a binary file. ibinun is a flag where 1 indicates that
716  ! cim should be written to a binary file if a binary file is open for it.
717  ! Set unit number for sorbate output
718  if (this%ioutsorbate /= 0) then
719  ibinun = 1
720  else
721  ibinun = 0
722  end if
723  if (idvsave == 0) ibinun = 0
724 
725  ! save sorbate concentration array
726  if (ibinun /= 0) then
727  iprint = 0
728  dinact = dhnoflo
729  if (this%ioutsorbate /= 0) then
730  ibinun = this%ioutsorbate
731  call this%dis%record_array(this%cimsrb, this%iout, iprint, ibinun, &
732  ' SORBATE', cdatafmp, nvaluesp, &
733  nwidthp, editdesc, dinact)
734  end if
735  end if
736 
738 
739  !> @ brief Output IST package budget summary.
740  !!
741  !! Output advanced boundary package budget summary. This method only needs
742  !! to be overridden for advanced packages that save budget summaries
743  !! to the model listing file.
744  !!
745  !<
746  subroutine ist_ot_bdsummary(this, kstp, kper, iout, ibudfl)
747  ! -- modules
748  use tdismodule, only: delt, totim
749  ! -- dummy variables
750  class(gwtisttype) :: this !< GwtIstType object
751  integer(I4B), intent(in) :: kstp !< time step number
752  integer(I4B), intent(in) :: kper !< period number
753  integer(I4B), intent(in) :: iout !< flag and unit number for the model listing file
754  integer(I4B), intent(in) :: ibudfl !< flag indicating budget should be written
755  ! -- local
756  integer(I4B) :: isuppress_output = 0
757  !
758  ! -- Fill budget terms
759  call this%budget%reset()
760  call this%budget%addentry(this%budterm, delt, budtxt, isuppress_output)
761  !
762  ! -- Write budget to list file
763  call this%budget%finalize_step(delt)
764  if (ibudfl /= 0) then
765  call this%budget%budget_ot(kstp, kper, iout)
766  end if
767  !
768  ! -- Write budget csv
769  call this%budget%writecsv(totim)
770  end subroutine ist_ot_bdsummary
771 
772  !> @ brief Deallocate package memory
773  !!
774  !! Deallocate package scalars and arrays.
775  !!
776  !<
777  subroutine ist_da(this)
778  ! -- modules
780  ! -- dummy
781  class(gwtisttype) :: this !< GwtIstType object
782  !
783  ! -- Deallocate arrays if package was active
784  if (this%inunit > 0) then
785  call mem_deallocate(this%icimout)
786  call mem_deallocate(this%ibudgetout)
787  call mem_deallocate(this%ibudcsv)
788  call mem_deallocate(this%ioutsorbate)
789  call mem_deallocate(this%idcy)
790  call mem_deallocate(this%isrb)
791  call mem_deallocate(this%kiter)
792  call mem_deallocate(this%cim)
793  call mem_deallocate(this%cimnew)
794  call mem_deallocate(this%cimold)
795  call mem_deallocate(this%cimsrb)
796  call mem_deallocate(this%zetaim)
797  call mem_deallocate(this%porosity)
798  call mem_deallocate(this%volfrac)
799  call mem_deallocate(this%bulk_density)
800  call mem_deallocate(this%distcoef)
801  call mem_deallocate(this%sp2)
802  call mem_deallocate(this%decay)
803  call mem_deallocate(this%decaylast)
804  call mem_deallocate(this%decayslast)
805  call mem_deallocate(this%decay_sorbed)
806  call mem_deallocate(this%strg)
807  this%fmi => null()
808  this%mst => null()
809  end if
810  !
811  ! -- Scalars
812  !
813  ! -- Objects
814  call this%budget%budget_da()
815  deallocate (this%budget)
816  call this%ocd%ocd_da()
817  deallocate (this%ocd)
818  if (associated(this%isotherm)) then
819  deallocate (this%isotherm)
820  nullify (this%isotherm)
821  end if
822  !
823  ! -- deallocate parent
824  call this%BndType%bnd_da()
825  end subroutine ist_da
826 
827  !> @ brief Allocate package scalars
828  !!
829  !! Allocate and initialize package scalars.
830  !!
831  !<
832  subroutine allocate_scalars(this)
833  ! -- modules
836  ! -- dummy
837  class(gwtisttype) :: this !< GwtIstType object
838  ! -- local
839  !
840  ! -- call standard BndType allocate scalars
841  call this%BndType%allocate_scalars()
842  !
843  ! -- Allocate
844  call mem_allocate(this%icimout, 'ICIMOUT', this%memoryPath)
845  call mem_allocate(this%ibudgetout, 'IBUDGETOUT', this%memoryPath)
846  call mem_allocate(this%ibudcsv, 'IBUDCSV', this%memoryPath)
847  call mem_allocate(this%ioutsorbate, 'IOUTSORBATE', this%memoryPath)
848  call mem_allocate(this%isrb, 'ISRB', this%memoryPath)
849  call mem_allocate(this%idcy, 'IDCY', this%memoryPath)
850  call mem_allocate(this%kiter, 'KITER', this%memoryPath)
851  !
852  ! -- Initialize
853  this%lstfmt = ''
854  this%icimout = 0
855  this%ibudgetout = 0
856  this%ibudcsv = 0
857  this%ioutsorbate = 0
858  this%isrb = 0
859  this%idcy = 0
860  this%kiter = 0
861  !
862  ! -- Create the ocd object, which is used to manage printing and saving
863  ! of the immobile domain concentrations
864  call ocd_cr(this%ocd)
865  end subroutine allocate_scalars
866 
867  !> @ brief Allocate package arrays
868  !!
869  !! Allocate and initialize package arrays.
870  !!
871  !<
872  subroutine ist_allocate_arrays(this)
873  ! -- modules
875  ! -- dummy
876  class(gwtisttype), intent(inout) :: this !< GwtIstType object
877  ! -- local
878  integer(I4B) :: n
879  !
880  ! -- call standard BndType allocate scalars
881  ! nbound and maxbound are 0 in order to keep memory footprint low
882  call this%BndType%allocate_arrays()
883  !
884  ! -- allocate ist arrays of size nodes
885  call mem_allocate(this%strg, this%dis%nodes, 'STRG', this%memoryPath)
886  call mem_allocate(this%cim, this%dis%nodes, 'CIM', this%memoryPath)
887  call mem_allocate(this%cimnew, this%dis%nodes, 'CIMNEW', this%memoryPath)
888  call mem_allocate(this%cimold, this%dis%nodes, 'CIMOLD', this%memoryPath)
889  call mem_allocate(this%porosity, this%dis%nodes, 'POROSITY', this%memoryPath)
890  call mem_allocate(this%zetaim, this%dis%nodes, 'ZETAIM', this%memoryPath)
891  call mem_allocate(this%volfrac, this%dis%nodes, 'VOLFRAC', this%memoryPath)
892  if (this%isrb == 0) then
893  call mem_allocate(this%bulk_density, 1, 'BULK_DENSITY', this%memoryPath)
894  call mem_allocate(this%distcoef, 1, 'DISTCOEF', this%memoryPath)
895  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
896  call mem_allocate(this%cimsrb, 1, 'SORBATE', this%memoryPath)
897  else
898  call mem_allocate(this%bulk_density, this%dis%nodes, 'BULK_DENSITY', &
899  this%memoryPath)
900  call mem_allocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
901  this%memoryPath)
902  call mem_allocate(this%cimsrb, this%dis%nodes, 'SORBATE', &
903  this%memoryPath)
904  if (this%isrb == 1) then
905  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
906  else
907  call mem_allocate(this%sp2, this%dis%nodes, 'SP2', this%memoryPath)
908  end if
909  end if
910  if (this%idcy == 0) then
911  call mem_allocate(this%decay, 1, 'DECAY', this%memoryPath)
912  call mem_allocate(this%decaylast, 1, 'DECAYLAST', this%memoryPath)
913  else
914  call mem_allocate(this%decay, this%dis%nodes, 'DECAY', this%memoryPath)
915  call mem_allocate(this%decaylast, this%dis%nodes, 'DECAYLAST', &
916  this%memoryPath)
917  end if
918  if (this%isrb == 0 .and. this%idcy == 0) then
919  call mem_allocate(this%decayslast, 1, 'DECAYSLAST', this%memoryPath)
920  else
921  call mem_allocate(this%decayslast, this%dis%nodes, 'DECAYSLAST', &
922  this%memoryPath)
923  end if
924  call mem_allocate(this%decay_sorbed, 1, 'DECAY_SORBED', this%memoryPath)
925  !
926  ! -- initialize
927  do n = 1, this%dis%nodes
928  this%strg(n) = dzero
929  this%cim(n) = dzero
930  this%cimnew(n) = dzero
931  this%cimold(n) = dzero
932  this%porosity(n) = dzero
933  this%zetaim(n) = dzero
934  this%volfrac(n) = dzero
935  end do
936  do n = 1, size(this%bulk_density)
937  this%bulk_density(n) = dzero
938  this%distcoef(n) = dzero
939  this%cimsrb(n) = dzero
940  end do
941  do n = 1, size(this%sp2)
942  this%sp2(n) = dzero
943  end do
944  do n = 1, size(this%decay)
945  this%decay(n) = dzero
946  this%decaylast(n) = dzero
947  end do
948  do n = 1, size(this%decayslast)
949  this%decayslast(n) = dzero
950  end do
951  !
952  ! -- Set pointers
953  this%ocd%dis => this%dis
954  end subroutine ist_allocate_arrays
955 
956  !> @ brief Source options for package
957  !!
958  !! Method to source options for the package.
959  !<
960  subroutine source_options(this)
961  ! -- modules
964  use openspecmodule, only: access, form
968  ! -- dummy
969  class(gwtisttype), intent(inout) :: this
970  ! -- locals
971  character(len=LINELENGTH) :: prnfmt
972  integer(I4B), pointer :: columns, width, digits
973  type(gwtistparamfoundtype) :: found
974  character(len=LENVARNAME), dimension(3) :: sorption_method = &
975  &[character(len=LENVARNAME) :: 'LINEAR', 'FREUNDLICH', 'LANGMUIR']
976  character(len=linelength) :: sorbate_fname, cim6_fname, budget_fname, &
977  budgetcsv_fname
978  allocate (columns)
979  allocate (width)
980  allocate (digits)
981  !
982  ! -- update defaults with memory sourced values
983  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', this%input_mempath, &
984  found%save_flows)
985  call mem_set_value(budget_fname, 'BUDGETFILE', this%input_mempath, &
986  found%budgetfile)
987  call mem_set_value(budgetcsv_fname, 'BUDGETCSVFILE', this%input_mempath, &
988  found%budgetcsvfile)
989  call mem_set_value(this%isrb, 'SORPTION', this%input_mempath, &
990  sorption_method, found%sorption)
991  call mem_set_value(this%idcy, 'ORDER1_DECAY', this%input_mempath, &
992  found%order1_decay)
993  call mem_set_value(this%idcy, 'ORDER0_DECAY', this%input_mempath, &
994  found%order0_decay)
995  call mem_set_value(cim6_fname, 'CIMFILE', this%input_mempath, &
996  found%cimfile)
997  call mem_set_value(sorbate_fname, 'SORBATEFILE', this%input_mempath, &
998  found%sorbatefile)
999  call mem_set_value(columns, 'COLUMNS', this%input_mempath, &
1000  found%columns)
1001  call mem_set_value(width, 'WIDTH', this%input_mempath, &
1002  found%width)
1003  call mem_set_value(digits, 'DIGITS', this%input_mempath, &
1004  found%digits)
1005  call mem_set_value(prnfmt, 'FORMAT', this%input_mempath, &
1006  found%format)
1007 
1008  ! -- found side effects
1009  if (found%save_flows) this%ipakcb = -1
1010  if (found%budgetfile) then
1011  call assign_iounit(this%ibudgetout, this%inunit, "BUDGET fileout")
1012  call openfile(this%ibudgetout, this%iout, budget_fname, 'DATA(BINARY)', &
1013  form, access, 'REPLACE', mode_opt=mnormal)
1014  end if
1015  if (found%budgetcsvfile) then
1016  call assign_iounit(this%ibudcsv, this%inunit, "BUDGETCSV fileout")
1017  call openfile(this%ibudcsv, this%iout, budgetcsv_fname, 'CSV', &
1018  filstat_opt='REPLACE')
1019  end if
1020  if (found%sorption) then
1021  if (this%isrb == 0) then
1022  call store_error('Unknown sorption type was specified. &
1023  &Sorption must be specified as LINEAR, &
1024  &FREUNDLICH, or LANGMUIR.')
1025  call store_error_filename(this%input_fname)
1026  end if
1027  end if
1028  if (found%order1_decay) this%idcy = 1
1029  if (found%order0_decay) this%idcy = 2
1030  if (found%cimfile) then
1031  call this%ocd%set_ocfile(cim6_fname, this%iout)
1032  end if
1033  if (found%columns .and. found%width .and. &
1034  found%digits .and. found%format) then
1035  write (this%lstfmt, '(a,i0,a,i0,a,i0,a)') 'COLUMNS ', columns, &
1036  ' WIDTH ', width, ' DIGITS ', digits, ' '//trim(prnfmt)
1037  end if
1038  if (found%sorbatefile) then
1039  this%ioutsorbate = getunit()
1040  call openfile(this%ioutsorbate, this%iout, sorbate_fname, &
1041  'DATA(BINARY)', form, access, 'REPLACE', mode_opt=mnormal)
1042  end if
1043  !
1044  ! -- log options
1045  if (this%iout > 0) then
1046  call this%log_options(found, cim6_fname, budget_fname, &
1047  budgetcsv_fname, sorbate_fname)
1048  end if
1049 
1050  deallocate (columns)
1051  deallocate (width)
1052  deallocate (digits)
1053  end subroutine source_options
1054 
1055  !> @brief Log user options to list file
1056  !<
1057  subroutine log_options(this, found, cim6_fname, budget_fname, &
1058  budgetcsv_fname, sorbate_fname)
1060  class(gwtisttype), intent(inout) :: this
1061  type(gwtistparamfoundtype), intent(in) :: found
1062  character(len=*), intent(in) :: cim6_fname
1063  character(len=*), intent(in) :: budget_fname
1064  character(len=*), intent(in) :: budgetcsv_fname
1065  character(len=*), intent(in) :: sorbate_fname
1066  ! -- formats
1067  character(len=*), parameter :: fmtisvflow = &
1068  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
1069  &WHENEVER ICBCFL IS NOT ZERO.')"
1070  character(len=*), parameter :: fmtlinear = &
1071  &"(4x,'LINEAR SORPTION IS SELECTED. ')"
1072  character(len=*), parameter :: fmtfreundlich = &
1073  &"(4x,'FREUNDLICH SORPTION IS ACTIVE. ')"
1074  character(len=*), parameter :: fmtlangmuir = &
1075  &"(4x,'LANGMUIR SORPTION IS ACTIVE. ')"
1076  character(len=*), parameter :: fmtidcy1 = &
1077  &"(4x,'FIRST-ORDER DECAY IS ACTIVE. ')"
1078  character(len=*), parameter :: fmtidcy2 = &
1079  &"(4x,'ZERO-ORDER DECAY IS ACTIVE. ')"
1080  character(len=*), parameter :: fmtistbin = &
1081  "(4x, 'IST ', 1x, a, 1x, ' WILL BE SAVED TO FILE: ', a, &
1082  &/4x, 'OPENED ON UNIT: ', I0)"
1083 
1084  write (this%iout, '(1x,a)') 'PROCESSING IMMOBILE STORAGE AND TRANSFER &
1085  &OPTIONS'
1086  if (found%save_flows) then
1087  write (this%iout, fmtisvflow)
1088  end if
1089  if (found%budgetfile) then
1090  write (this%iout, fmtistbin) 'BUDGET', trim(adjustl(budget_fname)), &
1091  this%ibudgetout
1092  end if
1093  if (found%budgetcsvfile) then
1094  write (this%iout, fmtistbin) 'BUDGET CSV', trim(adjustl(budgetcsv_fname)), &
1095  this%ibudcsv
1096  end if
1097  if (found%sorption) then
1098  select case (this%isrb)
1099  case (sorption_linear)
1100  write (this%iout, fmtlinear)
1101  case (sorption_freund)
1102  write (this%iout, fmtfreundlich)
1103  case (sorption_lang)
1104  write (this%iout, fmtlangmuir)
1105  end select
1106  end if
1107  if (found%order1_decay) then
1108  write (this%iout, fmtidcy1)
1109  end if
1110  if (found%order0_decay) then
1111  write (this%iout, fmtidcy2)
1112  end if
1113  if (found%sorbatefile) then
1114  write (this%iout, fmtistbin) &
1115  'SORBATE', sorbate_fname, this%ioutsorbate
1116  end if
1117  write (this%iout, '(1x,a)') 'END OF IMMOBILE STORAGE AND TRANSFER &
1118  &OPTIONS'
1119  end subroutine log_options
1120 
1121  !> @ brief Read options for package
1122  !!
1123  !! Read options for boundary packages.
1124  !!
1125  !<
1126  subroutine read_options(this)
1127  ! -- modules
1128  ! -- dummy
1129  class(gwtisttype), intent(inout) :: this !< GwtIstType object
1130 
1131  ! -- source options
1132  call this%source_options()
1133  end subroutine read_options
1134 
1135  !> @ brief Source data for package
1136  !!
1137  !! Method to source data for the package.
1138  !<
1139  subroutine source_data(this)
1140  ! -- modules
1141  use constantsmodule, only: linelength
1142  use simvariablesmodule, only: errmsg, warnmsg
1148  ! -- dummy
1149  class(gwtisttype) :: this
1150  ! -- locals
1151  !character(len=LINELENGTH) :: errmsg
1152  type(gwtistparamfoundtype) :: found
1153  integer(I4B) :: asize
1154  integer(I4B), dimension(:), pointer, contiguous :: map
1155  !
1156  ! -- set map to convert user input data into reduced data
1157  map => null()
1158  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1159  !
1160  ! -- reallocate
1161  if (this%isrb == 0) then
1162  call get_isize('BULK_DENSITY', this%input_mempath, asize)
1163  if (asize > 0) &
1164  call mem_reallocate(this%bulk_density, this%dis%nodes, &
1165  'BULK_DENSITY', this%memoryPath)
1166  call get_isize('DISTCOEF', this%input_mempath, asize)
1167  if (asize > 0) &
1168  call mem_reallocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
1169  this%memoryPath)
1170  end if
1171  if (this%idcy == 0) then
1172  call get_isize('DECAY', this%input_mempath, asize)
1173  if (asize > 0) &
1174  call mem_reallocate(this%decay, this%dis%nodes, 'DECAY', this%memoryPath)
1175  end if
1176  call get_isize('DECAY_SORBED', this%input_mempath, asize)
1177  if (asize > 0) then
1178  call mem_reallocate(this%decay_sorbed, this%dis%nodes, &
1179  'DECAY_SORBED', this%memoryPath)
1180  end if
1181  if (this%isrb < 2) then
1182  call get_isize('SP2', this%input_mempath, asize)
1183  if (asize > 0) &
1184  call mem_reallocate(this%sp2, this%dis%nodes, 'SP2', this%memoryPath)
1185  end if
1186  !
1187  ! -- update defaults with memory sourced values
1188  call mem_set_value(this%porosity, 'POROSITY', this%input_mempath, map, &
1189  found%porosity)
1190  call mem_set_value(this%volfrac, 'VOLFRAC', this%input_mempath, map, &
1191  found%volfrac)
1192  call mem_set_value(this%zetaim, 'ZETAIM', this%input_mempath, map, &
1193  found%zetaim)
1194  call mem_set_value(this%cim, 'CIM', this%input_mempath, map, &
1195  found%cim)
1196  call mem_set_value(this%decay, 'DECAY', this%input_mempath, map, &
1197  found%decay)
1198  call mem_set_value(this%decay_sorbed, 'DECAY_SORBED', this%input_mempath, &
1199  map, found%decay_sorbed)
1200  call mem_set_value(this%bulk_density, 'BULK_DENSITY', this%input_mempath, &
1201  map, found%bulk_density)
1202  call mem_set_value(this%distcoef, 'DISTCOEF', this%input_mempath, map, &
1203  found%distcoef)
1204  call mem_set_value(this%sp2, 'SP2', this%input_mempath, map, &
1205  found%sp2)
1206 
1207  ! -- log options
1208  if (this%iout > 0) then
1209  call this%log_data(found)
1210  end if
1211 
1212  ! -- Check for required sorption variables
1213  if (this%isrb > 0) then
1214  if (.not. found%bulk_density) then
1215  write (errmsg, '(a)') 'Sorption is active but BULK_DENSITY &
1216  &not specified. BULK_DENSITY must be specified in GRIDDATA block.'
1217  call store_error(errmsg)
1218  end if
1219  if (.not. found%distcoef) then
1220  write (errmsg, '(a)') 'Sorption is active but distribution &
1221  &coefficient not specified. DISTCOEF must be specified in &
1222  &GRIDDATA block.'
1223  call store_error(errmsg)
1224  end if
1225  if (this%isrb > 1) then
1226  if (.not. found%sp2) then
1227  write (errmsg, '(a)') 'Freundlich or langmuir sorption is active &
1228  &but SP2 not specified. SP2 must be specified in &
1229  &GRIDDATA block.'
1230  call store_error(errmsg)
1231  end if
1232  end if
1233  else
1234  if (found%bulk_density) then
1235  write (warnmsg, '(a)') 'Sorption is not active but &
1236  &BULK_DENSITY was specified. BULK_DENSITY will have no affect on &
1237  &simulation results.'
1238  call store_warning(warnmsg)
1239  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1240  end if
1241  if (found%distcoef) then
1242  write (warnmsg, '(a)') 'Sorption is not active but &
1243  &distribution coefficient was specified. DISTCOEF will have &
1244  &no affect on simulation results.'
1245  call store_warning(warnmsg)
1246  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1247  end if
1248  if (found%sp2) then
1249  write (warnmsg, '(a)') 'Sorption is not active but &
1250  &SP2 was specified. SP2 will have &
1251  &no affect on simulation results.'
1252  call store_warning(warnmsg)
1253  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1254  end if
1255  end if
1256 
1257  ! -- Check for required decay/production rate coefficients
1258  if (this%idcy > 0) then
1259  if (.not. found%decay) then
1260  write (errmsg, '(a)') 'First or zero order decay is &
1261  &active but the first rate coefficient is not specified. DECAY &
1262  &must be specified in GRIDDATA block.'
1263  call store_error(errmsg)
1264  end if
1265  if (.not. found%decay_sorbed) then
1266  !
1267  ! -- If DECAY_SORBED not specified and sorption is active, then
1268  ! terminate with an error
1269  write (errmsg, '(a)') 'DECAY_SORBED not provided in GRIDDATA &
1270  &block but decay and sorption are active. Specify DECAY_SORBED &
1271  &in GRIDDATA block.'
1272  call store_error(errmsg)
1273  end if
1274  else
1275  if (found%decay) then
1276  write (warnmsg, '(a)') 'First- or zero-order decay &
1277  &is not active but decay was specified. DECAY will &
1278  &have no affect on simulation results.'
1279  call store_warning(warnmsg)
1280  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1281  end if
1282  if (found%decay_sorbed) then
1283  write (warnmsg, '(a)') 'First- or zero-order decay &
1284  &is not active but DECAY_SORBED was specified. &
1285  &DECAY_SORBED will have no affect on simulation results.'
1286  call store_warning(warnmsg)
1287  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1288  end if
1289  end if
1290 
1291  ! -- Check for required dual domain arrays or warn if they are specified
1292  ! but won't be used.
1293  if (.not. found%cim) then
1294  write (this%iout, '(1x,a)') 'Warning. Dual domain is active but &
1295  &initial immobile domain concentration was not specified. &
1296  &Setting CIM to zero.'
1297  end if
1298  if (.not. found%zetaim) then
1299  write (errmsg, '(a)') 'Dual domain is active but dual &
1300  &domain mass transfer rate (ZETAIM) was not specified. ZETAIM &
1301  &must be specified in GRIDDATA block.'
1302  call store_error(errmsg)
1303  end if
1304  if (.not. found%porosity) then
1305  write (errmsg, '(a)') 'Dual domain is active but &
1306  &immobile domain POROSITY was not specified. POROSITY &
1307  &must be specified in GRIDDATA block.'
1308  call store_error(errmsg)
1309  end if
1310  if (.not. found%volfrac) then
1311  write (errmsg, '(a)') 'Dual domain is active but &
1312  &immobile domain VOLFRAC was not specified. VOLFRAC &
1313  &must be specified in GRIDDATA block. This is a new &
1314  &requirement for MODFLOW versions later than version &
1315  &6.4.1.'
1316  call store_error(errmsg)
1317  end if
1318 
1319  ! -- terminate if errors
1320  if (count_errors() > 0) then
1321  call store_error_filename(this%input_fname)
1322  end if
1323  end subroutine source_data
1324 
1325  !> @brief Log user data to list file
1326  !<
1327  subroutine log_data(this, found)
1329  class(gwtisttype), intent(inout) :: this
1330  type(gwtistparamfoundtype), intent(in) :: found
1331 
1332  write (this%iout, '(1x,a)') 'PROCESSING GRIDDATA'
1333  if (found%porosity) then
1334  write (this%iout, '(4x,a)') 'MOBILE DOMAIN POROSITY set from input file'
1335  end if
1336  if (found%volfrac) then
1337  write (this%iout, '(4x,a)') 'VOLFRAC set from input file'
1338  end if
1339  if (found%zetaim) then
1340  write (this%iout, '(4x,a)') 'ZETAIM set from input file'
1341  end if
1342  if (found%cim) then
1343  write (this%iout, '(4x,a)') 'CIM set from input file'
1344  end if
1345  if (found%decay) then
1346  write (this%iout, '(4x,a)') 'DECAY RATE set from input file'
1347  end if
1348  if (found%decay_sorbed) then
1349  write (this%iout, '(4x,a)') 'DECAY SORBED RATE set from input file'
1350  end if
1351  if (found%bulk_density) then
1352  write (this%iout, '(4x,a)') 'BULK DENSITY set from input file'
1353  end if
1354  if (found%distcoef) then
1355  write (this%iout, '(4x,a)') 'DISTRIBUTION COEFFICIENT set from input file'
1356  end if
1357  if (found%sp2) then
1358  write (this%iout, '(4x,a)') 'SECOND SORPTION PARAM set from input file'
1359  end if
1360  write (this%iout, '(1x,a)') 'END PROCESSING GRIDDATA'
1361  end subroutine log_data
1362 
1363  !> @ brief Read dimensions for package
1364  !!
1365  !! Read dimensions for package.
1366  !!
1367  !<
1368  subroutine ist_read_dimensions(this)
1369  ! -- dummy
1370  class(gwtisttype), intent(inout) :: this !< GwtIstType object
1371  ! -- local
1372  ! -- format
1373  end subroutine ist_read_dimensions
1374 
1375  !> @ brief Return thetaim
1376  !!
1377  !! Calculate and return thetaim, volume of immobile voids per aquifer volume
1378  !!
1379  !<
1380  function get_thetaim(this, node) result(thetaim)
1381  ! -- modules
1382  ! -- dummy
1383  class(gwtisttype) :: this !< GwtIstType object
1384  integer(I4B), intent(in) :: node !< node number
1385  ! -- return
1386  real(dp) :: thetaim
1387  !
1388  thetaim = this%volfrac(node) * this%porosity(node)
1389  end function get_thetaim
1390 
1391  !> @ brief Calculate immobile domain equation terms
1392  !!
1393  !! This subroutine calculates the immobile domain (or dual domain) terms.
1394  !! The resulting ddterm and f terms are described in the GWT model report.
1395  !! The terms are concentration coefficients used in the balance equation
1396  !! for the immobile domain.
1397  !!
1398  !<
1399  subroutine get_ddterm(thetaim, vcell, delt, swtpdt, &
1400  volfracim, rhobim, kdnew, kdold, &
1401  lambda1im, lambda2im, &
1402  gamma1im, gamma2im, zetaim, ddterm, f)
1403  ! -- dummy
1404  real(DP), intent(in) :: thetaim !< immobile domain porosity
1405  real(DP), intent(in) :: vcell !< volume of cell
1406  real(DP), intent(in) :: delt !< length of time step
1407  real(DP), intent(in) :: swtpdt !< cell saturation at end of time step
1408  real(DP), intent(in) :: volfracim !< volume fraction of immobile domain
1409  real(DP), intent(in) :: rhobim !< bulk density for the immobile domain (fim * rhob)
1410  real(DP), intent(in) :: kdnew !< effective distribution coefficient for new time
1411  real(DP), intent(in) :: kdold !< effective distribution coefficient for old time
1412  real(DP), intent(in) :: lambda1im !< first-order decay rate in aqueous phase
1413  real(DP), intent(in) :: lambda2im !< first-order decay rate in sorbed phase
1414  real(DP), intent(in) :: gamma1im !< zero-order decay rate in aqueous phase
1415  real(DP), intent(in) :: gamma2im !< zero-order decay rate in sorbed phase
1416  real(DP), intent(in) :: zetaim !< transfer coefficient between mobile and immobile domains
1417  real(DP), dimension(:), intent(inout) :: ddterm !< nine terms comprising the balance equation of the immobile domain
1418  real(DP), intent(inout) :: f !< the f term used to calculate the immobile domain concentration
1419  ! -- local
1420  real(DP) :: tled
1421  !
1422  ! -- initialize
1423  tled = done / delt
1424  !
1425  ! -- Calculate terms. These terms correspond to the concentration
1426  ! coefficients in equation 7-4 of the GWT model report. However,
1427  ! an updated equation is presented as 9-9 in the supplemental technical
1428  ! information guide (mf6suptechinfo.pdf)
1429  ddterm(1) = thetaim * vcell * tled
1430  ddterm(2) = thetaim * vcell * tled
1431  ddterm(3) = volfracim * rhobim * vcell * kdnew * tled
1432  ddterm(4) = volfracim * rhobim * vcell * kdold * tled
1433  ddterm(5) = thetaim * lambda1im * vcell
1434  ddterm(6) = lambda2im * volfracim * rhobim * kdnew * vcell
1435  ddterm(7) = thetaim * gamma1im * vcell
1436  ddterm(8) = gamma2im * volfracim * rhobim * vcell
1437  ddterm(9) = vcell * swtpdt * zetaim
1438  !
1439  ! -- calculate denominator term, f
1440  f = ddterm(1) + ddterm(3) + ddterm(5) + ddterm(6) + ddterm(9)
1441  end subroutine get_ddterm
1442 
1443  !> @ brief Calculate the hcof and rhs terms for immobile domain
1444  !!
1445  !! This subroutine calculates the hcof and rhs terms that must
1446  !! be added to the solution system of equations
1447  !!
1448  !<
1449  subroutine get_hcofrhs(ddterm, f, cimold, hcof, rhs)
1450  ! -- dummy
1451  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1452  real(DP), intent(in) :: f !< the f term used to calculate the immobile domain concentration
1453  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1454  real(DP), intent(inout) :: hcof !< calculated contribution for the a-matrix diagonal position
1455  real(DP), intent(inout) :: rhs !< calculated contribution for the solution right-hand side
1456  !
1457  ! -- calculate hcof
1458  hcof = ddterm(9)**2 / f - ddterm(9)
1459  !
1460  ! -- calculate rhs, and switch the sign because this term needs to
1461  ! be moved to the left hand side
1462  rhs = (ddterm(2) + ddterm(4)) * cimold - ddterm(7) - ddterm(8)
1463  rhs = rhs * ddterm(9) / f
1464  rhs = -rhs
1465  end subroutine get_hcofrhs
1466 
1467  !> @ brief Calculate the immobile domain concentration
1468  !!
1469  !! This function calculates the concentration of the immobile domain.
1470  !!
1471  !<
1472  function get_ddconc(ddterm, f, cimold, cnew) result(cimnew)
1473  ! -- dummy
1474  real(dp), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1475  real(dp), intent(in) :: f !< the f term used to calculate the immobile domain concentration
1476  real(dp), intent(in) :: cimold !< immobile domain concentration at end of last time step
1477  real(dp), intent(in) :: cnew !< concentration of the mobile domain at the end of the time step
1478  ! -- result
1479  real(dp) :: cimnew !< calculated concentration of the immobile domain
1480  ! -- local
1481  !
1482  ! -- calculate ddconc
1483  cimnew = (ddterm(2) + ddterm(4)) * cimold + ddterm(9) * cnew - ddterm(7) &
1484  - ddterm(8)
1485  cimnew = cimnew / f
1486  end function get_ddconc
1487 
1488  !> @ brief Calculate the immobile domain budget terms
1489  !!
1490  !! This subroutine calculates and accumulates the immobile domain
1491  !! budget terms into the budterm accumulator
1492  !!
1493  !<
1494  subroutine accumulate_budterm(budterm, ddterm, cimnew, cimold, cnew, idcy)
1495  ! -- modules
1496  ! -- dummy
1497  real(DP), dimension(:, :), intent(inout) :: budterm !<
1498  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1499  real(DP), intent(in) :: cimnew !< immobile domain concenration at the end of this time step
1500  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1501  real(DP), intent(in) :: cnew !< mobile domain concentration at the end of this time step
1502  integer(I4B), intent(in) :: idcy !< order of decay rate (0:none, 1:first, 2:zero)
1503  ! -- local
1504  real(DP) :: rate
1505  integer(I4B) :: i
1506  !
1507  ! -- calculate STORAGE-AQUEOUS
1508  i = 1
1509  rate = -ddterm(1) * cimnew + ddterm(2) * cimold
1510  if (rate > dzero) then
1511  budterm(1, i) = budterm(1, i) + rate
1512  else
1513  budterm(2, i) = budterm(2, i) - rate
1514  end if
1515  !
1516  ! -- calculate STORAGE-SORBED
1517  i = 2
1518  rate = -ddterm(3) * cimnew + ddterm(4) * cimold
1519  if (rate > dzero) then
1520  budterm(1, i) = budterm(1, i) + rate
1521  else
1522  budterm(2, i) = budterm(2, i) - rate
1523  end if
1524  !
1525  ! -- calculate DECAY-AQUEOUS
1526  i = 3
1527  rate = dzero
1528  if (idcy == 1) then
1529  rate = -ddterm(5) * cimnew
1530  else if (idcy == 2) then
1531  rate = -ddterm(7)
1532  else
1533  rate = dzero
1534  end if
1535  if (rate > dzero) then
1536  budterm(1, i) = budterm(1, i) + rate
1537  else
1538  budterm(2, i) = budterm(2, i) - rate
1539  end if
1540  !
1541  ! -- calculate DECAY-SORBED
1542  i = 4
1543  if (idcy == 1) then
1544  rate = -ddterm(6) * cimnew
1545  else if (idcy == 2) then
1546  rate = -ddterm(8)
1547  else
1548  rate = dzero
1549  end if
1550  if (rate > dzero) then
1551  budterm(1, i) = budterm(1, i) + rate
1552  else
1553  budterm(2, i) = budterm(2, i) - rate
1554  end if
1555  !
1556  ! -- calculate MOBILE-DOMAIN
1557  i = 5
1558  rate = ddterm(9) * cnew - ddterm(9) * cimnew
1559  if (rate > dzero) then
1560  budterm(1, i) = budterm(1, i) + rate
1561  else
1562  budterm(2, i) = budterm(2, i) - rate
1563  end if
1564  !
1565  end subroutine accumulate_budterm
1566 
1567  !> @ brief Get effective Freundlich distribution coefficient
1568  !<
1569  function get_freundlich_kd(conc, kf, a) result(kd)
1570  ! -- dummy
1571  real(dp), intent(in) :: conc !< solute concentration
1572  real(dp), intent(in) :: kf !< freundlich constant
1573  real(dp), intent(in) :: a !< freundlich exponent
1574  ! -- return
1575  real(dp) :: kd !< effective distribution coefficient
1576  !
1577  if (conc > dzero) then
1578  kd = kf * conc**(a - done)
1579  else
1580  kd = dzero
1581  end if
1582  end function get_freundlich_kd
1583 
1584  !> @ brief Get effective Langmuir distribution coefficient
1585  !<
1586  function get_langmuir_kd(conc, kl, sbar) result(kd)
1587  ! -- dummy
1588  real(dp), intent(in) :: conc !< solute concentration
1589  real(dp), intent(in) :: kl !< langmuir constant
1590  real(dp), intent(in) :: sbar !< langmuir sorption sites
1591  ! -- return
1592  real(dp) :: kd !< effective distribution coefficient
1593  !
1594  if (conc > dzero) then
1595  kd = (kl * sbar) / (done + kl * conc)
1596  else
1597  kd = dzero
1598  end if
1599  end function get_langmuir_kd
1600 
1601 end module gwtistmodule
This module contains the base boundary package.
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public budget_cr(this, name_model)
@ brief Create a new budget object
Definition: Budget.f90:84
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
@ mnormal
normal output mode
Definition: Constants.f90:206
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
– @ brief Immobile Storage and Transfer (IST) Module
Definition: gwt-ist.f90:14
subroutine ist_da(this)
@ brief Deallocate package memory
Definition: gwt-ist.f90:778
subroutine log_data(this, found)
Log user data to list file.
Definition: gwt-ist.f90:1328
subroutine ist_calc_csrb(this, cim)
@ brief Calculate immobile sorbed concentration
Definition: gwt-ist.f90:549
subroutine ist_rp(this)
@ brief Read and prepare method for package
Definition: gwt-ist.f90:234
real(dp) function get_ddconc(ddterm, f, cimold, cnew)
@ brief Calculate the immobile domain concentration
Definition: gwt-ist.f90:1473
real(dp) function get_thetaim(this, node)
@ brief Return thetaim
Definition: gwt-ist.f90:1381
subroutine source_data(this)
@ brief Source data for package
Definition: gwt-ist.f90:1140
subroutine get_hcofrhs(ddterm, f, cimold, hcof, rhs)
@ brief Calculate the hcof and rhs terms for immobile domain
Definition: gwt-ist.f90:1450
subroutine get_ddterm(thetaim, vcell, delt, swtpdt, volfracim, rhobim, kdnew, kdold, lambda1im, lambda2im, gamma1im, gamma2im, zetaim, ddterm, f)
@ brief Calculate immobile domain equation terms
Definition: gwt-ist.f90:1403
integer(i4b), parameter nbditems
Definition: gwt-ist.f90:37
subroutine output_immobile_sorbate_concentration(this, idvsave, idvprint)
@ brief Output immobile domain sorbate concentration.
Definition: gwt-ist.f90:704
real(dp) function get_freundlich_kd(conc, kf, a)
@ brief Get effective Freundlich distribution coefficient
Definition: gwt-ist.f90:1570
subroutine source_options(this)
@ brief Source options for package
Definition: gwt-ist.f90:961
real(dp) function get_langmuir_kd(conc, kl, sbar)
@ brief Get effective Langmuir distribution coefficient
Definition: gwt-ist.f90:1587
subroutine ist_ot_model_flows(this, icbcfl, ibudfl, icbcun, imap)
@ brief Output model flow terms.
Definition: gwt-ist.f90:597
subroutine ist_read_dimensions(this)
@ brief Read dimensions for package
Definition: gwt-ist.f90:1369
subroutine log_options(this, found, cim6_fname, budget_fname, budgetcsv_fname, sorbate_fname)
Log user options to list file.
Definition: gwt-ist.f90:1059
subroutine allocate_scalars(this)
@ brief Allocate package scalars
Definition: gwt-ist.f90:833
character(len=lenftype) ftype
Definition: gwt-ist.f90:35
subroutine ist_ar(this)
@ brief Allocate and read method for package
Definition: gwt-ist.f90:168
subroutine ist_ot_bdsummary(this, kstp, kper, iout, ibudfl)
@ brief Output IST package budget summary.
Definition: gwt-ist.f90:747
subroutine ist_ot_dv(this, idvsave, idvprint)
@ brief Output dependent variables.
Definition: gwt-ist.f90:660
subroutine ist_allocate_arrays(this)
@ brief Allocate package arrays
Definition: gwt-ist.f90:873
subroutine ist_cq(this, x, flowja, iadv)
@ brief Calculate package flows.
Definition: gwt-ist.f90:403
subroutine ist_ad(this)
@ brief Advance the ist package
Definition: gwt-ist.f90:247
subroutine accumulate_budterm(budterm, ddterm, cimnew, cimold, cnew, idcy)
@ brief Calculate the immobile domain budget terms
Definition: gwt-ist.f90:1495
subroutine ist_bd(this, model_budget)
@ brief Add package flows to model budget.
Definition: gwt-ist.f90:574
subroutine ist_fc(this, rhs, ia, idxglo, matrix_sln)
@ brief Fill coefficient method for package
Definition: gwt-ist.f90:276
character(len=lenbudtxt), dimension(nbditems) budtxt
Definition: gwt-ist.f90:38
subroutine read_options(this)
@ brief Read options for package
Definition: gwt-ist.f90:1127
subroutine, public ist_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath, fmi, mst)
@ brief Create a new package object
Definition: gwt-ist.f90:121
character(len=lenpackagename) text
Definition: gwt-ist.f90:36
subroutine output_immobile_concentration(this, idvsave, idvprint)
@ brief Output immobile domain aqueous concentration.
Definition: gwt-ist.f90:673
– @ brief Mobile Storage and Transfer (MST) Module
Definition: gwt-mst.f90:10
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
Definition: gwt-mst.f90:1451
subroutine, public assign_iounit(iounit, errunit, description)
@ brief assign io unit number
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
integer(i4b), parameter sorption_off
Sorption is inactive (default)
Definition: IsothermEnum.f90:6
integer(i4b), parameter sorption_freund
Freundlich sorption between aqueous and solid phases.
Definition: IsothermEnum.f90:8
integer(i4b), parameter sorption_lang
Langmuir sorption between aqueous and solid phases.
Definition: IsothermEnum.f90:9
integer(i4b), parameter sorption_linear
Linear sorption between aqueous and solid phases.
Definition: IsothermEnum.f90:7
class(isothermtype) function, pointer, public create_isotherm(isotherm_type, distcoef, sp2)
Create an isotherm object based on type and parameters.
This module defines variable data types.
Definition: kind.f90:8
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Output control data module.
subroutine, public ocd_cr(ocdobj)
@ brief Create a new output control data type.
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
integer(i4b) ifailedstepretry
current retry for this time step
character(len=maxcharlen) warnmsg
warning message string
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
@ brief BndType
Derived type for the Budget object.
Definition: Budget.f90:39
@ brief Immobile storage and transfer
Definition: gwt-ist.f90:53
@ brief Mobile storage and transfer
Definition: gwt-mst.f90:49