MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwt-mst.f90
Go to the documentation of this file.
1 !> -- @ brief Mobile Storage and Transfer (MST) Module
2 !!
3 !! The GwtMstModule contains the GwtMstType, which is the
4 !! derived type responsible for adding the effects of
5 !! 1. Changes in dissolved solute mass
6 !! 2. Decay of dissolved solute mass
7 !! 3. Sorption
8 !! 4. Decay of sorbed solute mass
9 !<
11 
12  use kindmodule, only: dp, i4b
15  use simmodule, only: store_error, count_errors, &
19  use basedismodule, only: disbasetype
20  use tspfmimodule, only: tspfmitype
21 
22  implicit none
23  public :: gwtmsttype
24  public :: mst_cr
25  !
26  integer(I4B), parameter :: nbditems = 4
27  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt
28  data budtxt/' STORAGE-AQUEOUS', ' DECAY-AQUEOUS', &
29  ' STORAGE-SORBED', ' DECAY-SORBED'/
30 
31  !> @ brief Mobile storage and transfer
32  !!
33  !! Data and methods for handling changes in solute storage,
34  !! decay of dissolved solute mass, sorption, and decay of
35  !! sorbed mass.
36  !<
37  type, extends(numericalpackagetype) :: gwtmsttype
38  !
39  ! -- storage
40  real(dp), dimension(:), pointer, contiguous :: porosity => null() !< mobile porosity defined as volume mobile voids per volume of mobile domain
41  real(dp), dimension(:), pointer, contiguous :: thetam => null() !< mobile porosity defined as volume mobile voids per volume of aquifer
42  real(dp), dimension(:), pointer, contiguous :: volfracim => null() !< sum of all immobile domain volume fractions
43  real(dp), dimension(:), pointer, contiguous :: ratesto => null() !< rate of mobile storage
44  !
45  ! -- decay
46  integer(I4B), pointer :: idcy => null() !< order of decay rate (0:none, 1:first, 2:zero)
47  real(dp), dimension(:), pointer, contiguous :: decay => null() !< first or zero order decay rate (aqueous)
48  real(dp), dimension(:), pointer, contiguous :: decay_sorbed => null() !< first or zero order decay rate (sorbed)
49  real(dp), dimension(:), pointer, contiguous :: ratedcy => null() !< rate of decay
50  real(dp), dimension(:), pointer, contiguous :: decaylast => null() !< decay rate used for last iteration (needed for zero order decay)
51  real(dp), dimension(:), pointer, contiguous :: decayslast => null() !< sorbed decay rate used for last iteration (needed for zero order decay)
52  !
53  ! -- sorption
54  integer(I4B), pointer :: isrb => null() !< sorption active flag (0:off, 1:linear, 2:freundlich, 3:langmuir)
55  real(dp), dimension(:), pointer, contiguous :: bulk_density => null() !< bulk density of mobile domain; mass of mobile domain solid per aquifer volume
56  real(dp), dimension(:), pointer, contiguous :: distcoef => null() !< kd distribution coefficient
57  real(dp), dimension(:), pointer, contiguous :: sp2 => null() !< second sorption parameter
58  real(dp), dimension(:), pointer, contiguous :: ratesrb => null() !< rate of sorption
59  real(dp), dimension(:), pointer, contiguous :: ratedcys => null() !< rate of sorbed mass decay
60  !
61  ! -- misc
62  integer(I4B), dimension(:), pointer, contiguous :: ibound => null() !< pointer to model ibound
63  type(tspfmitype), pointer :: fmi => null() !< pointer to fmi object
64 
65  contains
66 
67  procedure :: mst_ar
68  procedure :: mst_fc
69  procedure :: mst_fc_sto
70  procedure :: mst_fc_dcy
71  procedure :: mst_fc_srb
72  procedure :: mst_fc_dcy_srb
73  procedure :: mst_cq
74  procedure :: mst_cq_sto
75  procedure :: mst_cq_dcy
76  procedure :: mst_cq_srb
77  procedure :: mst_cq_dcy_srb
78  procedure :: mst_bd
79  procedure :: mst_ot_flow
80  procedure :: mst_da
81  procedure :: allocate_scalars
82  procedure :: addto_volfracim
83  procedure :: get_volfracm
84  procedure, private :: allocate_arrays
85  procedure, private :: read_options
86  procedure, private :: read_data
87 
88  end type gwtmsttype
89 
90 contains
91 
92  !> @ brief Create a new package object
93  !!
94  !! Create a new MST object
95  !!
96  !<
97  subroutine mst_cr(mstobj, name_model, inunit, iout, fmi)
98  ! -- dummy
99  type(gwtmsttype), pointer :: mstobj !< unallocated new mst object to create
100  character(len=*), intent(in) :: name_model !< name of the model
101  integer(I4B), intent(in) :: inunit !< unit number of WEL package input file
102  integer(I4B), intent(in) :: iout !< unit number of model listing file
103  type(tspfmitype), intent(in), target :: fmi !< fmi package for this GWT model
104  !
105  ! -- Create the object
106  allocate (mstobj)
107  !
108  ! -- create name and memory path
109  call mstobj%set_names(1, name_model, 'MST', 'MST')
110  !
111  ! -- Allocate scalars
112  call mstobj%allocate_scalars()
113  !
114  ! -- Set variables
115  mstobj%inunit = inunit
116  mstobj%iout = iout
117  mstobj%fmi => fmi
118  !
119  ! -- Initialize block parser
120  call mstobj%parser%Initialize(mstobj%inunit, mstobj%iout)
121  !
122  ! -- Return
123  return
124  end subroutine mst_cr
125 
126  !> @ brief Allocate and read method for package
127  !!
128  !! Method to allocate and read static data for the package.
129  !!
130  !<
131  subroutine mst_ar(this, dis, ibound)
132  ! -- modules
133  ! -- dummy
134  class(gwtmsttype), intent(inout) :: this !< GwtMstType object
135  class(disbasetype), pointer, intent(in) :: dis !< pointer to dis package
136  integer(I4B), dimension(:), pointer, contiguous :: ibound !< pointer to GWT ibound array
137  ! -- local
138  ! -- formats
139  character(len=*), parameter :: fmtmst = &
140  "(1x,/1x,'MST -- MOBILE STORAGE AND TRANSFER PACKAGE, VERSION 1, &
141  &7/29/2020 INPUT READ FROM UNIT ', i0, //)"
142  !
143  ! --print a message identifying the mobile storage and transfer package.
144  write (this%iout, fmtmst) this%inunit
145  !
146  ! -- Read options
147  call this%read_options()
148  !
149  ! -- store pointers to arguments that were passed in
150  this%dis => dis
151  this%ibound => ibound
152  !
153  ! -- Allocate arrays
154  call this%allocate_arrays(dis%nodes)
155  !
156  ! -- read the data block
157  call this%read_data()
158  !
159  ! -- Return
160  return
161  end subroutine mst_ar
162 
163  !> @ brief Fill coefficient method for package
164  !!
165  !! Method to calculate and fill coefficients for the package.
166  !!
167  !<
168  subroutine mst_fc(this, nodes, cold, nja, matrix_sln, idxglo, cnew, &
169  rhs, kiter)
170  ! -- modules
171  ! -- dummy
172  class(gwtmsttype) :: this !< GwtMstType object
173  integer, intent(in) :: nodes !< number of nodes
174  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
175  integer(I4B), intent(in) :: nja !< number of GWT connections
176  class(matrixbasetype), pointer :: matrix_sln !< solution matrix
177  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
178  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
179  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
180  integer(I4B), intent(in) :: kiter !< solution outer iteration number
181  ! -- local
182  !
183  ! -- storage contribution
184  call this%mst_fc_sto(nodes, cold, nja, matrix_sln, idxglo, rhs)
185  !
186  ! -- decay contribution
187  if (this%idcy /= 0) then
188  call this%mst_fc_dcy(nodes, cold, cnew, nja, matrix_sln, idxglo, &
189  rhs, kiter)
190  end if
191  !
192  ! -- sorption contribution
193  if (this%isrb /= 0) then
194  call this%mst_fc_srb(nodes, cold, nja, matrix_sln, idxglo, rhs, cnew)
195  end if
196  !
197  ! -- decay sorbed contribution
198  if (this%isrb /= 0 .and. this%idcy /= 0) then
199  call this%mst_fc_dcy_srb(nodes, cold, nja, matrix_sln, idxglo, rhs, &
200  cnew, kiter)
201  end if
202  !
203  ! -- Return
204  return
205  end subroutine mst_fc
206 
207  !> @ brief Fill storage coefficient method for package
208  !!
209  !! Method to calculate and fill storage coefficients for the package.
210  !!
211  !<
212  subroutine mst_fc_sto(this, nodes, cold, nja, matrix_sln, idxglo, rhs)
213  ! -- modules
214  use tdismodule, only: delt
215  ! -- dummy
216  class(gwtmsttype) :: this !< GwtMstType object
217  integer, intent(in) :: nodes !< number of nodes
218  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
219  integer(I4B), intent(in) :: nja !< number of GWT connections
220  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
221  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
222  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
223  ! -- local
224  integer(I4B) :: n, idiag
225  real(DP) :: tled
226  real(DP) :: hhcof, rrhs
227  real(DP) :: vnew, vold
228  !
229  ! -- set variables
230  tled = done / delt
231  !
232  ! -- loop through and calculate storage contribution to hcof and rhs
233  do n = 1, this%dis%nodes
234  !
235  ! -- skip if transport inactive
236  if (this%ibound(n) <= 0) cycle
237  !
238  ! -- calculate new and old water volumes
239  vnew = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n)) * &
240  this%fmi%gwfsat(n) * this%thetam(n)
241  vold = vnew
242  if (this%fmi%igwfstrgss /= 0) vold = vold + this%fmi%gwfstrgss(n) * delt
243  if (this%fmi%igwfstrgsy /= 0) vold = vold + this%fmi%gwfstrgsy(n) * delt
244  !
245  ! -- add terms to diagonal and rhs accumulators
246  hhcof = -vnew * tled
247  rrhs = -vold * tled * cold(n)
248  idiag = this%dis%con%ia(n)
249  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
250  rhs(n) = rhs(n) + rrhs
251  end do
252  !
253  ! -- Return
254  return
255  end subroutine mst_fc_sto
256 
257  !> @ brief Fill decay coefficient method for package
258  !!
259  !! Method to calculate and fill decay coefficients for the package.
260  !!
261  !<
262  subroutine mst_fc_dcy(this, nodes, cold, cnew, nja, matrix_sln, &
263  idxglo, rhs, kiter)
264  ! -- modules
265  use tdismodule, only: delt
266  ! -- dummy
267  class(gwtmsttype) :: this !< GwtMstType object
268  integer, intent(in) :: nodes !< number of nodes
269  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
270  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
271  integer(I4B), intent(in) :: nja !< number of GWT connections
272  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
273  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
274  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
275  integer(I4B), intent(in) :: kiter !< solution outer iteration number
276  ! -- local
277  integer(I4B) :: n, idiag
278  real(DP) :: hhcof, rrhs
279  real(DP) :: swtpdt
280  real(DP) :: vcell
281  real(DP) :: decay_rate
282  !
283  ! -- loop through and calculate decay contribution to hcof and rhs
284  do n = 1, this%dis%nodes
285  !
286  ! -- skip if transport inactive
287  if (this%ibound(n) <= 0) cycle
288  !
289  ! -- calculate new and old water volumes
290  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
291  swtpdt = this%fmi%gwfsat(n)
292  !
293  ! -- add decay rate terms to accumulators
294  idiag = this%dis%con%ia(n)
295  if (this%idcy == 1) then
296  !
297  ! -- first order decay rate is a function of concentration, so add
298  ! to left hand side
299  hhcof = -this%decay(n) * vcell * swtpdt * this%thetam(n)
300  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
301  elseif (this%idcy == 2) then
302  !
303  ! -- Call function to get zero-order decay rate, which may be changed
304  ! from the user-specified rate to prevent negative concentrations
305  decay_rate = get_zero_order_decay(this%decay(n), this%decaylast(n), &
306  kiter, cold(n), cnew(n), delt)
307  this%decaylast(n) = decay_rate
308  rrhs = decay_rate * vcell * swtpdt * this%thetam(n)
309  rhs(n) = rhs(n) + rrhs
310  end if
311  !
312  end do
313  !
314  ! -- Return
315  return
316  end subroutine mst_fc_dcy
317 
318  !> @ brief Fill sorption coefficient method for package
319  !!
320  !! Method to calculate and fill sorption coefficients for the package.
321  !!
322  !<
323  subroutine mst_fc_srb(this, nodes, cold, nja, matrix_sln, idxglo, rhs, &
324  cnew)
325  ! -- modules
326  use tdismodule, only: delt
327  ! -- dummy
328  class(gwtmsttype) :: this !< GwtMstType object
329  integer, intent(in) :: nodes !< number of nodes
330  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
331  integer(I4B), intent(in) :: nja !< number of GWT connections
332  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
333  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
334  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
335  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
336  ! -- local
337  integer(I4B) :: n, idiag
338  real(DP) :: tled
339  real(DP) :: hhcof, rrhs
340  real(DP) :: swt, swtpdt
341  real(DP) :: vcell
342  real(DP) :: const1
343  real(DP) :: const2
344  real(DP) :: volfracm
345  real(DP) :: rhobm
346  !
347  ! -- set variables
348  tled = done / delt
349  !
350  ! -- loop through and calculate sorption contribution to hcof and rhs
351  do n = 1, this%dis%nodes
352  !
353  ! -- skip if transport inactive
354  if (this%ibound(n) <= 0) cycle
355  !
356  ! -- assign variables
357  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
358  swtpdt = this%fmi%gwfsat(n)
359  swt = this%fmi%gwfsatold(n, delt)
360  idiag = this%dis%con%ia(n)
361  const1 = this%distcoef(n)
362  const2 = 0.
363  if (this%isrb > 1) const2 = this%sp2(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)
369  !
370  ! -- Add hhcof to diagonal and rrhs to right-hand side
371  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
372  rhs(n) = rhs(n) + rrhs
373  !
374  end do
375  !
376  ! -- Return
377  return
378  end subroutine mst_fc_srb
379 
380  !> @ brief Calculate sorption terms
381  !!
382  !! Subroutine to calculate sorption terms
383  !!
384  !<
385  subroutine mst_srb_term(isrb, volfracm, rhobm, vcell, tled, cnew, cold, &
386  swnew, swold, const1, const2, rate, hcofval, rhsval)
387  ! -- modules
388  ! -- dummy
389  integer(I4B), intent(in) :: isrb !< sorption flag 1, 2, 3 are linear, freundlich, and langmuir
390  real(DP), intent(in) :: volfracm !< volume fraction of mobile domain (fhat_m)
391  real(DP), intent(in) :: rhobm !< bulk density of mobile domain (rhob_m)
392  real(DP), intent(in) :: vcell !< volume of cell
393  real(DP), intent(in) :: tled !< one over time step length
394  real(DP), intent(in) :: cnew !< concentration at end of this time step
395  real(DP), intent(in) :: cold !< concentration at end of last time step
396  real(DP), intent(in) :: swnew !< cell saturation at end of this time step
397  real(DP), intent(in) :: swold !< cell saturation at end of last time step
398  real(DP), intent(in) :: const1 !< distribution coefficient or freundlich or langmuir constant
399  real(DP), intent(in) :: const2 !< zero, freundlich exponent, or langmuir sorption sites
400  real(DP), intent(out), optional :: rate !< calculated sorption rate
401  real(DP), intent(out), optional :: hcofval !< diagonal contribution to solution coefficient matrix
402  real(DP), intent(out), optional :: rhsval !< contribution to solution right-hand-side
403  ! -- local
404  real(DP) :: term
405  real(DP) :: derv
406  real(DP) :: cbarnew
407  real(DP) :: cbarold
408  real(DP) :: cavg
409  real(DP) :: cbaravg
410  real(DP) :: swavg
411  !
412  ! -- Calculate based on type of sorption
413  if (isrb == 1) then
414  ! -- linear
415  term = -volfracm * rhobm * vcell * tled * const1
416  if (present(hcofval)) hcofval = term * swnew
417  if (present(rhsval)) rhsval = term * swold * cold
418  if (present(rate)) rate = term * swnew * cnew - term * swold * cold
419  else
420  !
421  ! -- calculate average aqueous concentration
422  cavg = dhalf * (cold + cnew)
423  !
424  ! -- set values based on isotherm
425  if (isrb == 2) then
426  ! -- freundlich
427  cbarnew = get_freundlich_conc(cnew, const1, const2)
428  cbarold = get_freundlich_conc(cold, const1, const2)
429  derv = get_freundlich_derivative(cavg, const1, const2)
430  else if (isrb == 3) then
431  ! -- langmuir
432  cbarnew = get_langmuir_conc(cnew, const1, const2)
433  cbarold = get_langmuir_conc(cold, const1, const2)
434  derv = get_langmuir_derivative(cavg, const1, const2)
435  end if
436  !
437  ! -- calculate hcof, rhs, and rate for freundlich and langmuir
438  term = -volfracm * rhobm * vcell * tled
439  cbaravg = (cbarold + cbarnew) * dhalf
440  swavg = (swnew + swold) * dhalf
441  if (present(hcofval)) then
442  hcofval = term * derv * swavg
443  end if
444  if (present(rhsval)) then
445  rhsval = term * derv * swavg * cold - term * cbaravg * (swnew - swold)
446  end if
447  if (present(rate)) then
448  rate = term * derv * swavg * (cnew - cold) &
449  + term * cbaravg * (swnew - swold)
450  end if
451  end if
452  return
453  end subroutine mst_srb_term
454 
455  !> @ brief Fill sorption-decay coefficient method for package
456  !!
457  !! Method to calculate and fill sorption-decay coefficients for the package.
458  !!
459  !<
460  subroutine mst_fc_dcy_srb(this, nodes, cold, nja, matrix_sln, idxglo, &
461  rhs, cnew, kiter)
462  ! -- modules
463  use tdismodule, only: delt
464  ! -- dummy
465  class(gwtmsttype) :: this !< GwtMstType object
466  integer, intent(in) :: nodes !< number of nodes
467  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
468  integer(I4B), intent(in) :: nja !< number of GWT connections
469  class(matrixbasetype), pointer :: matrix_sln !< solution coefficient matrix
470  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
471  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
472  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
473  integer(I4B), intent(in) :: kiter !< solution outer iteration number
474  ! -- local
475  integer(I4B) :: n, idiag
476  real(DP) :: hhcof, rrhs
477  real(DP) :: vcell
478  real(DP) :: swnew
479  real(DP) :: distcoef
480  real(DP) :: volfracm
481  real(DP) :: rhobm
482  real(DP) :: term
483  real(DP) :: csrb
484  real(DP) :: decay_rate
485  real(DP) :: csrbold
486  real(DP) :: csrbnew
487  !
488  ! -- loop through and calculate sorption contribution to hcof and rhs
489  do n = 1, this%dis%nodes
490  !
491  ! -- skip if transport inactive
492  if (this%ibound(n) <= 0) cycle
493  !
494  ! -- set variables
495  hhcof = dzero
496  rrhs = dzero
497  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
498  swnew = this%fmi%gwfsat(n)
499  distcoef = this%distcoef(n)
500  idiag = this%dis%con%ia(n)
501  volfracm = this%get_volfracm(n)
502  rhobm = this%bulk_density(n)
503  term = this%decay_sorbed(n) * volfracm * rhobm * swnew * vcell
504  !
505  ! -- add sorbed mass decay rate terms to accumulators
506  if (this%idcy == 1) then
507  !
508  if (this%isrb == 1) then
509  !
510  ! -- first order decay rate is a function of concentration, so add
511  ! to left hand side
512  hhcof = -term * distcoef
513  else if (this%isrb == 2) then
514  !
515  ! -- nonlinear Freundlich sorption, so add to RHS
516  csrb = get_freundlich_conc(cnew(n), distcoef, this%sp2(n))
517  rrhs = term * csrb
518  else if (this%isrb == 3) then
519  !
520  ! -- nonlinear Lanmuir sorption, so add to RHS
521  csrb = get_langmuir_conc(cnew(n), distcoef, this%sp2(n))
522  rrhs = term * csrb
523  end if
524  elseif (this%idcy == 2) then
525  !
526  ! -- Call function to get zero-order decay rate, which may be changed
527  ! from the user-specified rate to prevent negative concentrations
528  if (distcoef > dzero) then
529 
530  if (this%isrb == 1) then
531  csrbold = cold(n) * distcoef
532  csrbnew = cnew(n) * distcoef
533  else if (this%isrb == 2) then
534  csrbold = get_freundlich_conc(cold(n), distcoef, this%sp2(n))
535  csrbnew = get_freundlich_conc(cnew(n), distcoef, this%sp2(n))
536  else if (this%isrb == 3) then
537  csrbold = get_langmuir_conc(cold(n), distcoef, this%sp2(n))
538  csrbnew = get_langmuir_conc(cnew(n), distcoef, this%sp2(n))
539  end if
540 
541  decay_rate = get_zero_order_decay(this%decay_sorbed(n), &
542  this%decayslast(n), &
543  kiter, csrbold, csrbnew, delt)
544  this%decayslast(n) = decay_rate
545  rrhs = decay_rate * volfracm * rhobm * swnew * vcell
546  end if
547 
548  end if
549  !
550  ! -- Add hhcof to diagonal and rrhs to right-hand side
551  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
552  rhs(n) = rhs(n) + rrhs
553  !
554  end do
555  !
556  ! -- Return
557  return
558  end subroutine mst_fc_dcy_srb
559 
560  !> @ brief Calculate flows for package
561  !!
562  !! Method to calculate flows for the package.
563  !!
564  !<
565  subroutine mst_cq(this, nodes, cnew, cold, flowja)
566  ! -- modules
567  ! -- dummy
568  class(gwtmsttype) :: this !< GwtMstType object
569  integer(I4B), intent(in) :: nodes !< number of nodes
570  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
571  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
572  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
573  ! -- local
574  !
575  ! - storage
576  call this%mst_cq_sto(nodes, cnew, cold, flowja)
577  !
578  ! -- decay
579  if (this%idcy /= 0) then
580  call this%mst_cq_dcy(nodes, cnew, cold, flowja)
581  end if
582  !
583  ! -- sorption
584  if (this%isrb /= 0) then
585  call this%mst_cq_srb(nodes, cnew, cold, flowja)
586  end if
587  !
588  ! -- decay sorbed
589  if (this%isrb /= 0 .and. this%idcy /= 0) then
590  call this%mst_cq_dcy_srb(nodes, cnew, cold, flowja)
591  end if
592  !
593  ! -- Return
594  return
595  end subroutine mst_cq
596 
597  !> @ brief Calculate storage terms for package
598  !!
599  !! Method to calculate storage terms for the package.
600  !!
601  !<
602  subroutine mst_cq_sto(this, nodes, cnew, cold, flowja)
603  ! -- modules
604  use tdismodule, only: delt
605  ! -- dummy
606  class(gwtmsttype) :: this !< GwtMstType object
607  integer(I4B), intent(in) :: nodes !< number of nodes
608  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
609  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
610  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
611  ! -- local
612  integer(I4B) :: n
613  integer(I4B) :: idiag
614  real(DP) :: rate
615  real(DP) :: tled
616  real(DP) :: vnew, vold
617  real(DP) :: hhcof, rrhs
618  !
619  ! -- initialize
620  tled = done / delt
621  !
622  ! -- Calculate storage change
623  do n = 1, nodes
624  this%ratesto(n) = dzero
625  !
626  ! -- skip if transport inactive
627  if (this%ibound(n) <= 0) cycle
628  !
629  ! -- calculate new and old water volumes
630  vnew = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n)) * &
631  this%fmi%gwfsat(n) * this%thetam(n)
632  vold = vnew
633  if (this%fmi%igwfstrgss /= 0) vold = vold + this%fmi%gwfstrgss(n) * delt
634  if (this%fmi%igwfstrgsy /= 0) vold = vold + this%fmi%gwfstrgsy(n) * delt
635  !
636  ! -- calculate rate
637  hhcof = -vnew * tled
638  rrhs = -vold * tled * cold(n)
639  rate = hhcof * cnew(n) - rrhs
640  this%ratesto(n) = rate
641  idiag = this%dis%con%ia(n)
642  flowja(idiag) = flowja(idiag) + rate
643  end do
644  !
645  ! -- Return
646  return
647  end subroutine mst_cq_sto
648 
649  !> @ brief Calculate decay terms for package
650  !!
651  !! Method to calculate decay terms for the package.
652  !!
653  !<
654  subroutine mst_cq_dcy(this, nodes, cnew, cold, flowja)
655  ! -- modules
656  use tdismodule, only: delt
657  ! -- dummy
658  class(gwtmsttype) :: this !< GwtMstType object
659  integer(I4B), intent(in) :: nodes !< number of nodes
660  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
661  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
662  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
663  ! -- local
664  integer(I4B) :: n
665  integer(I4B) :: idiag
666  real(DP) :: rate
667  real(DP) :: swtpdt
668  real(DP) :: hhcof, rrhs
669  real(DP) :: vcell
670  real(DP) :: decay_rate
671  !
672  ! -- initialize
673  !
674  ! -- Calculate decay change
675  do n = 1, nodes
676  !
677  ! -- skip if transport inactive
678  this%ratedcy(n) = dzero
679  if (this%ibound(n) <= 0) cycle
680  !
681  ! -- calculate new and old water volumes
682  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
683  swtpdt = this%fmi%gwfsat(n)
684  !
685  ! -- calculate decay gains and losses
686  rate = dzero
687  hhcof = dzero
688  rrhs = dzero
689  if (this%idcy == 1) then
690  hhcof = -this%decay(n) * vcell * swtpdt * this%thetam(n)
691  elseif (this%idcy == 2) then
692  decay_rate = get_zero_order_decay(this%decay(n), this%decaylast(n), &
693  0, cold(n), cnew(n), delt)
694  rrhs = decay_rate * vcell * swtpdt * this%thetam(n)
695  end if
696  rate = hhcof * cnew(n) - rrhs
697  this%ratedcy(n) = rate
698  idiag = this%dis%con%ia(n)
699  flowja(idiag) = flowja(idiag) + rate
700  !
701  end do
702  !
703  ! -- Return
704  return
705  end subroutine mst_cq_dcy
706 
707  !> @ brief Calculate sorption terms for package
708  !!
709  !! Method to calculate sorption terms for the package.
710  !!
711  !<
712  subroutine mst_cq_srb(this, nodes, cnew, cold, flowja)
713  ! -- modules
714  use tdismodule, only: delt
715  ! -- dummy
716  class(gwtmsttype) :: this !< GwtMstType object
717  integer(I4B), intent(in) :: nodes !< number of nodes
718  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
719  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
720  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
721  ! -- local
722  integer(I4B) :: n
723  integer(I4B) :: idiag
724  real(DP) :: rate
725  real(DP) :: tled
726  real(DP) :: swt, swtpdt
727  real(DP) :: vcell
728  real(DP) :: volfracm
729  real(DP) :: rhobm
730  real(DP) :: const1
731  real(DP) :: const2
732  !
733  ! -- initialize
734  tled = done / delt
735  !
736  ! -- Calculate sorption change
737  do n = 1, nodes
738  !
739  ! -- initialize rates
740  this%ratesrb(n) = dzero
741  !
742  ! -- skip if transport inactive
743  if (this%ibound(n) <= 0) cycle
744  !
745  ! -- assign variables
746  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
747  swtpdt = this%fmi%gwfsat(n)
748  swt = this%fmi%gwfsatold(n, delt)
749  volfracm = this%get_volfracm(n)
750  rhobm = this%bulk_density(n)
751  const1 = this%distcoef(n)
752  const2 = 0.
753  if (this%isrb > 1) const2 = this%sp2(n)
754  call mst_srb_term(this%isrb, volfracm, rhobm, vcell, tled, cnew(n), &
755  cold(n), swtpdt, swt, const1, const2, &
756  rate=rate)
757  this%ratesrb(n) = rate
758  idiag = this%dis%con%ia(n)
759  flowja(idiag) = flowja(idiag) + rate
760  !
761  end do
762  !
763  ! -- Return
764  return
765  end subroutine mst_cq_srb
766 
767  !> @ brief Calculate decay-sorption terms for package
768  !!
769  !! Method to calculate decay-sorption terms for the package.
770  !!
771  !<
772  subroutine mst_cq_dcy_srb(this, nodes, cnew, cold, flowja)
773  ! -- modules
774  use tdismodule, only: delt
775  ! -- dummy
776  class(gwtmsttype) :: this !< GwtMstType object
777  integer(I4B), intent(in) :: nodes !< number of nodes
778  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
779  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
780  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
781  ! -- local
782  integer(I4B) :: n
783  integer(I4B) :: idiag
784  real(DP) :: rate
785  real(DP) :: hhcof, rrhs
786  real(DP) :: vcell
787  real(DP) :: swnew
788  real(DP) :: distcoef
789  real(DP) :: volfracm
790  real(DP) :: rhobm
791  real(DP) :: term
792  real(DP) :: csrb
793  real(DP) :: csrbnew
794  real(DP) :: csrbold
795  real(DP) :: decay_rate
796  !
797  ! -- Calculate sorbed decay change
798  ! This routine will only be called if sorption and decay are active
799  do n = 1, nodes
800  !
801  ! -- initialize rates
802  this%ratedcys(n) = dzero
803  !
804  ! -- skip if transport inactive
805  if (this%ibound(n) <= 0) cycle
806  !
807  ! -- set variables
808  hhcof = dzero
809  rrhs = dzero
810  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
811  swnew = this%fmi%gwfsat(n)
812  distcoef = this%distcoef(n)
813  volfracm = this%get_volfracm(n)
814  rhobm = this%bulk_density(n)
815  term = this%decay_sorbed(n) * volfracm * rhobm * swnew * vcell
816  !
817  ! -- add sorbed mass decay rate terms to accumulators
818  if (this%idcy == 1) then
819  !
820  if (this%isrb == 1) then
821  !
822  ! -- first order decay rate is a function of concentration, so add
823  ! to left hand side
824  hhcof = -term * distcoef
825  else if (this%isrb == 2) then
826  !
827  ! -- nonlinear Freundlich sorption, so add to RHS
828  csrb = get_freundlich_conc(cnew(n), distcoef, this%sp2(n))
829  rrhs = term * csrb
830  else if (this%isrb == 3) then
831  !
832  ! -- nonlinear Lanmuir sorption, so add to RHS
833  csrb = get_langmuir_conc(cnew(n), distcoef, this%sp2(n))
834  rrhs = term * csrb
835  end if
836  elseif (this%idcy == 2) then
837  !
838  ! -- Call function to get zero-order decay rate, which may be changed
839  ! from the user-specified rate to prevent negative concentrations
840  if (distcoef > dzero) then
841  if (this%isrb == 1) then
842  csrbold = cold(n) * distcoef
843  csrbnew = cnew(n) * distcoef
844  else if (this%isrb == 2) then
845  csrbold = get_freundlich_conc(cold(n), distcoef, this%sp2(n))
846  csrbnew = get_freundlich_conc(cnew(n), distcoef, this%sp2(n))
847  else if (this%isrb == 3) then
848  csrbold = get_langmuir_conc(cold(n), distcoef, this%sp2(n))
849  csrbnew = get_langmuir_conc(cnew(n), distcoef, this%sp2(n))
850  end if
851  decay_rate = get_zero_order_decay(this%decay_sorbed(n), &
852  this%decayslast(n), &
853  0, csrbold, csrbnew, delt)
854  rrhs = decay_rate * volfracm * rhobm * swnew * vcell
855  end if
856  end if
857  !
858  ! -- calculate rate
859  rate = hhcof * cnew(n) - rrhs
860  this%ratedcys(n) = rate
861  idiag = this%dis%con%ia(n)
862  flowja(idiag) = flowja(idiag) + rate
863  !
864  end do
865  !
866  ! -- Return
867  return
868  end subroutine mst_cq_dcy_srb
869 
870  !> @ brief Calculate budget terms for package
871  !!
872  !! Method to calculate budget terms for the package.
873  !!
874  !<
875  subroutine mst_bd(this, isuppress_output, model_budget)
876  ! -- modules
877  use tdismodule, only: delt
879  ! -- dummy
880  class(gwtmsttype) :: this !< GwtMstType object
881  integer(I4B), intent(in) :: isuppress_output !< flag to suppress output
882  type(budgettype), intent(inout) :: model_budget !< model budget object
883  ! -- local
884  real(DP) :: rin
885  real(DP) :: rout
886  !
887  ! -- sto
888  call rate_accumulator(this%ratesto, rin, rout)
889  call model_budget%addentry(rin, rout, delt, budtxt(1), &
890  isuppress_output, rowlabel=this%packName)
891  !
892  ! -- dcy
893  if (this%idcy /= 0) then
894  call rate_accumulator(this%ratedcy, rin, rout)
895  call model_budget%addentry(rin, rout, delt, budtxt(2), &
896  isuppress_output, rowlabel=this%packName)
897  end if
898  !
899  ! -- srb
900  if (this%isrb /= 0) then
901  call rate_accumulator(this%ratesrb, rin, rout)
902  call model_budget%addentry(rin, rout, delt, budtxt(3), &
903  isuppress_output, rowlabel=this%packName)
904  end if
905  !
906  ! -- srb dcy
907  if (this%isrb /= 0 .and. this%idcy /= 0) then
908  call rate_accumulator(this%ratedcys, rin, rout)
909  call model_budget%addentry(rin, rout, delt, budtxt(4), &
910  isuppress_output, rowlabel=this%packName)
911  end if
912  !
913  ! -- Return
914  return
915  end subroutine mst_bd
916 
917  !> @ brief Output flow terms for package
918  !!
919  !! Method to output terms for the package.
920  !!
921  !<
922  subroutine mst_ot_flow(this, icbcfl, icbcun)
923  ! -- dummy
924  class(gwtmsttype) :: this !< GwtMstType object
925  integer(I4B), intent(in) :: icbcfl !< flag and unit number for cell-by-cell output
926  integer(I4B), intent(in) :: icbcun !< flag indication if cell-by-cell data should be saved
927  ! -- local
928  integer(I4B) :: ibinun
929  !character(len=16), dimension(2) :: aname
930  integer(I4B) :: iprint, nvaluesp, nwidthp
931  character(len=1) :: cdatafmp = ' ', editdesc = ' '
932  real(DP) :: dinact
933  !
934  ! -- Set unit number for binary output
935  if (this%ipakcb < 0) then
936  ibinun = icbcun
937  elseif (this%ipakcb == 0) then
938  ibinun = 0
939  else
940  ibinun = this%ipakcb
941  end if
942  if (icbcfl == 0) ibinun = 0
943  !
944  ! -- Record the storage rate if requested
945  if (ibinun /= 0) then
946  iprint = 0
947  dinact = dzero
948  !
949  ! -- sto
950  call this%dis%record_array(this%ratesto, this%iout, iprint, -ibinun, &
951  budtxt(1), cdatafmp, nvaluesp, &
952  nwidthp, editdesc, dinact)
953  !
954  ! -- dcy
955  if (this%idcy /= 0) &
956  call this%dis%record_array(this%ratedcy, this%iout, iprint, -ibinun, &
957  budtxt(2), cdatafmp, nvaluesp, &
958  nwidthp, editdesc, dinact)
959  !
960  ! -- srb
961  if (this%isrb /= 0) &
962  call this%dis%record_array(this%ratesrb, this%iout, iprint, -ibinun, &
963  budtxt(3), cdatafmp, nvaluesp, &
964  nwidthp, editdesc, dinact)
965  !
966  ! -- dcy srb
967  if (this%isrb /= 0 .and. this%idcy /= 0) &
968  call this%dis%record_array(this%ratedcys, this%iout, iprint, -ibinun, &
969  budtxt(4), cdatafmp, nvaluesp, &
970  nwidthp, editdesc, dinact)
971  end if
972  !
973  ! -- Return
974  return
975  end subroutine mst_ot_flow
976 
977  !> @ brief Deallocate
978  !!
979  !! Method to deallocate memory for the package.
980  !!
981  !<
982  subroutine mst_da(this)
983  ! -- modules
985  ! -- dummy
986  class(gwtmsttype) :: this !< GwtMstType object
987  !
988  ! -- Deallocate arrays if package was active
989  if (this%inunit > 0) then
990  call mem_deallocate(this%porosity)
991  call mem_deallocate(this%thetam)
992  call mem_deallocate(this%volfracim)
993  call mem_deallocate(this%ratesto)
994  call mem_deallocate(this%idcy)
995  call mem_deallocate(this%decay)
996  call mem_deallocate(this%decay_sorbed)
997  call mem_deallocate(this%ratedcy)
998  call mem_deallocate(this%decaylast)
999  call mem_deallocate(this%decayslast)
1000  call mem_deallocate(this%isrb)
1001  call mem_deallocate(this%bulk_density)
1002  call mem_deallocate(this%distcoef)
1003  call mem_deallocate(this%sp2)
1004  call mem_deallocate(this%ratesrb)
1005  call mem_deallocate(this%ratedcys)
1006  this%ibound => null()
1007  this%fmi => null()
1008  end if
1009  !
1010  ! -- Scalars
1011  !
1012  ! -- deallocate parent
1013  call this%NumericalPackageType%da()
1014  !
1015  ! -- Return
1016  return
1017  end subroutine mst_da
1018 
1019  !> @ brief Allocate scalar variables for package
1020  !!
1021  !! Method to allocate scalar variables for the package.
1022  !!
1023  !<
1024  subroutine allocate_scalars(this)
1025  ! -- modules
1027  ! -- dummy
1028  class(gwtmsttype) :: this !< GwtMstType object
1029  ! -- local
1030  !
1031  ! -- Allocate scalars in NumericalPackageType
1032  call this%NumericalPackageType%allocate_scalars()
1033  !
1034  ! -- Allocate
1035  call mem_allocate(this%isrb, 'ISRB', this%memoryPath)
1036  call mem_allocate(this%idcy, 'IDCY', this%memoryPath)
1037  !
1038  ! -- Initialize
1039  this%isrb = 0
1040  this%idcy = 0
1041  !
1042  ! -- Return
1043  return
1044  end subroutine allocate_scalars
1045 
1046  !> @ brief Allocate arrays for package
1047  !!
1048  !! Method to allocate arrays for the package.
1049  !!
1050  !<
1051  subroutine allocate_arrays(this, nodes)
1052  ! -- modules
1054  use constantsmodule, only: dzero
1055  ! -- dummy
1056  class(gwtmsttype) :: this !< GwtMstType object
1057  integer(I4B), intent(in) :: nodes !< number of nodes
1058  ! -- local
1059  integer(I4B) :: n
1060  !
1061  ! -- Allocate
1062  ! -- sto
1063  call mem_allocate(this%porosity, nodes, 'POROSITY', this%memoryPath)
1064  call mem_allocate(this%thetam, nodes, 'THETAM', this%memoryPath)
1065  call mem_allocate(this%volfracim, nodes, 'VOLFRACIM', this%memoryPath)
1066  call mem_allocate(this%ratesto, nodes, 'RATESTO', this%memoryPath)
1067  !
1068  ! -- dcy
1069  if (this%idcy == 0) then
1070  call mem_allocate(this%ratedcy, 1, 'RATEDCY', this%memoryPath)
1071  call mem_allocate(this%decay, 1, 'DECAY', this%memoryPath)
1072  call mem_allocate(this%decaylast, 1, 'DECAYLAST', this%memoryPath)
1073  else
1074  call mem_allocate(this%ratedcy, this%dis%nodes, 'RATEDCY', this%memoryPath)
1075  call mem_allocate(this%decay, nodes, 'DECAY', this%memoryPath)
1076  call mem_allocate(this%decaylast, nodes, 'DECAYLAST', this%memoryPath)
1077  end if
1078  if (this%idcy /= 0 .and. this%isrb /= 0) then
1079  call mem_allocate(this%ratedcys, this%dis%nodes, 'RATEDCYS', &
1080  this%memoryPath)
1081  call mem_allocate(this%decayslast, this%dis%nodes, 'DECAYSLAST', &
1082  this%memoryPath)
1083  else
1084  call mem_allocate(this%ratedcys, 1, 'RATEDCYS', this%memoryPath)
1085  call mem_allocate(this%decayslast, 1, 'DECAYSLAST', this%memoryPath)
1086  end if
1087  call mem_allocate(this%decay_sorbed, 1, 'DECAY_SORBED', &
1088  this%memoryPath)
1089  !
1090  ! -- srb
1091  if (this%isrb == 0) then
1092  call mem_allocate(this%bulk_density, 1, 'BULK_DENSITY', this%memoryPath)
1093  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
1094  call mem_allocate(this%distcoef, 1, 'DISTCOEF', this%memoryPath)
1095  call mem_allocate(this%ratesrb, 1, 'RATESRB', this%memoryPath)
1096  else
1097  call mem_allocate(this%bulk_density, nodes, 'BULK_DENSITY', this%memoryPath)
1098  call mem_allocate(this%distcoef, nodes, 'DISTCOEF', this%memoryPath)
1099  call mem_allocate(this%ratesrb, nodes, 'RATESRB', this%memoryPath)
1100  if (this%isrb == 1) then
1101  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
1102  else
1103  call mem_allocate(this%sp2, nodes, 'SP2', this%memoryPath)
1104  end if
1105  end if
1106  !
1107  ! -- Initialize
1108  do n = 1, nodes
1109  this%porosity(n) = dzero
1110  this%thetam(n) = dzero
1111  this%volfracim(n) = dzero
1112  this%ratesto(n) = dzero
1113  end do
1114  do n = 1, size(this%decay)
1115  this%decay(n) = dzero
1116  this%ratedcy(n) = dzero
1117  this%decaylast(n) = dzero
1118  end do
1119  do n = 1, size(this%bulk_density)
1120  this%bulk_density(n) = dzero
1121  this%distcoef(n) = dzero
1122  this%ratesrb(n) = dzero
1123  end do
1124  do n = 1, size(this%sp2)
1125  this%sp2(n) = dzero
1126  end do
1127  do n = 1, size(this%ratedcys)
1128  this%ratedcys(n) = dzero
1129  this%decayslast(n) = dzero
1130  end do
1131  !
1132  ! -- Return
1133  return
1134  end subroutine allocate_arrays
1135 
1136  !> @ brief Read options for package
1137  !!
1138  !! Method to read options for the package.
1139  !!
1140  !<
1141  subroutine read_options(this)
1142  ! -- modules
1143  use constantsmodule, only: linelength
1144  ! -- dummy
1145  class(gwtmsttype) :: this !< GwtMstType object
1146  ! -- local
1147  character(len=LINELENGTH) :: keyword, keyword2
1148  integer(I4B) :: ierr
1149  logical :: isfound, endOfBlock
1150  ! -- formats
1151  character(len=*), parameter :: fmtisvflow = &
1152  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
1153  &WHENEVER ICBCFL IS NOT ZERO.')"
1154  character(len=*), parameter :: fmtisrb = &
1155  &"(4x,'LINEAR SORPTION IS ACTIVE. ')"
1156  character(len=*), parameter :: fmtfreundlich = &
1157  &"(4x,'FREUNDLICH SORPTION IS ACTIVE. ')"
1158  character(len=*), parameter :: fmtlangmuir = &
1159  &"(4x,'LANGMUIR SORPTION IS ACTIVE. ')"
1160  character(len=*), parameter :: fmtidcy1 = &
1161  &"(4x,'FIRST-ORDER DECAY IS ACTIVE. ')"
1162  character(len=*), parameter :: fmtidcy2 = &
1163  &"(4x,'ZERO-ORDER DECAY IS ACTIVE. ')"
1164  !
1165  ! -- get options block
1166  call this%parser%GetBlock('OPTIONS', isfound, ierr, blockrequired=.false., &
1167  supportopenclose=.true.)
1168  !
1169  ! -- parse options block if detected
1170  if (isfound) then
1171  write (this%iout, '(1x,a)') 'PROCESSING MOBILE STORAGE AND TRANSFER OPTIONS'
1172  do
1173  call this%parser%GetNextLine(endofblock)
1174  if (endofblock) exit
1175  call this%parser%GetStringCaps(keyword)
1176  select case (keyword)
1177  case ('SAVE_FLOWS')
1178  this%ipakcb = -1
1179  write (this%iout, fmtisvflow)
1180  case ('SORBTION', 'SORPTION')
1181  this%isrb = 1
1182  call this%parser%GetStringCaps(keyword2)
1183  if (trim(adjustl(keyword2)) == 'LINEAR') this%isrb = 1
1184  if (trim(adjustl(keyword2)) == 'FREUNDLICH') this%isrb = 2
1185  if (trim(adjustl(keyword2)) == 'LANGMUIR') this%isrb = 3
1186  select case (this%isrb)
1187  case (1)
1188  write (this%iout, fmtisrb)
1189  case (2)
1190  write (this%iout, fmtfreundlich)
1191  case (3)
1192  write (this%iout, fmtlangmuir)
1193  end select
1194  case ('FIRST_ORDER_DECAY')
1195  this%idcy = 1
1196  write (this%iout, fmtidcy1)
1197  case ('ZERO_ORDER_DECAY')
1198  this%idcy = 2
1199  write (this%iout, fmtidcy2)
1200  case default
1201  write (errmsg, '(a,a)') 'Unknown MST option: ', trim(keyword)
1202  call store_error(errmsg)
1203  call this%parser%StoreErrorUnit()
1204  end select
1205  end do
1206  write (this%iout, '(1x,a)') 'END OF MOBILE STORAGE AND TRANSFER OPTIONS'
1207  end if
1208  !
1209  ! -- Return
1210  return
1211  end subroutine read_options
1212 
1213  !> @ brief Read data for package
1214  !!
1215  !! Method to read data for the package.
1216  !!
1217  !<
1218  subroutine read_data(this)
1219  ! -- modules
1220  use constantsmodule, only: linelength
1222  ! -- dummy
1223  class(gwtmsttype) :: this !< GwtMstType object
1224  ! -- local
1225  character(len=LINELENGTH) :: keyword
1226  character(len=:), allocatable :: line
1227  integer(I4B) :: istart, istop, lloc, ierr, n
1228  logical :: isfound, endOfBlock
1229  logical, dimension(6) :: lname
1230  character(len=24), dimension(6) :: aname
1231  ! -- formats
1232  ! -- data
1233  data aname(1)/' MOBILE DOMAIN POROSITY'/
1234  data aname(2)/' BULK DENSITY'/
1235  data aname(3)/'DISTRIBUTION COEFFICIENT'/
1236  data aname(4)/' DECAY RATE'/
1237  data aname(5)/' DECAY SORBED RATE'/
1238  data aname(6)/' SECOND SORPTION PARAM'/
1239  !
1240  ! -- initialize
1241  isfound = .false.
1242  lname(:) = .false.
1243  !
1244  ! -- get griddata block
1245  call this%parser%GetBlock('GRIDDATA', isfound, ierr)
1246  if (isfound) then
1247  write (this%iout, '(1x,a)') 'PROCESSING GRIDDATA'
1248  do
1249  call this%parser%GetNextLine(endofblock)
1250  if (endofblock) exit
1251  call this%parser%GetStringCaps(keyword)
1252  call this%parser%GetRemainingLine(line)
1253  lloc = 1
1254  select case (keyword)
1255  case ('POROSITY')
1256  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1257  this%parser%iuactive, this%porosity, &
1258  aname(1))
1259  lname(1) = .true.
1260  case ('BULK_DENSITY')
1261  if (this%isrb == 0) &
1262  call mem_reallocate(this%bulk_density, this%dis%nodes, &
1263  'BULK_DENSITY', trim(this%memoryPath))
1264  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1265  this%parser%iuactive, &
1266  this%bulk_density, aname(2))
1267  lname(2) = .true.
1268  case ('DISTCOEF')
1269  if (this%isrb == 0) &
1270  call mem_reallocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
1271  trim(this%memoryPath))
1272  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1273  this%parser%iuactive, this%distcoef, &
1274  aname(3))
1275  lname(3) = .true.
1276  case ('DECAY')
1277  if (this%idcy == 0) &
1278  call mem_reallocate(this%decay, this%dis%nodes, 'DECAY', &
1279  trim(this%memoryPath))
1280  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1281  this%parser%iuactive, this%decay, &
1282  aname(4))
1283  lname(4) = .true.
1284  case ('DECAY_SORBED')
1285  call mem_reallocate(this%decay_sorbed, this%dis%nodes, &
1286  'DECAY_SORBED', trim(this%memoryPath))
1287  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1288  this%parser%iuactive, &
1289  this%decay_sorbed, aname(5))
1290  lname(5) = .true.
1291  case ('SP2')
1292  if (this%isrb < 2) &
1293  call mem_reallocate(this%sp2, this%dis%nodes, 'SP2', &
1294  trim(this%memoryPath))
1295  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1296  this%parser%iuactive, this%sp2, &
1297  aname(6))
1298  lname(6) = .true.
1299  case default
1300  write (errmsg, '(a,a)') 'Unknown GRIDDATA tag: ', trim(keyword)
1301  call store_error(errmsg)
1302  call this%parser%StoreErrorUnit()
1303  end select
1304  end do
1305  write (this%iout, '(1x,a)') 'END PROCESSING GRIDDATA'
1306  else
1307  write (errmsg, '(a)') 'Required GRIDDATA block not found.'
1308  call store_error(errmsg)
1309  call this%parser%StoreErrorUnit()
1310  end if
1311  !
1312  ! -- Check for required porosity
1313  if (.not. lname(1)) then
1314  write (errmsg, '(a)') 'POROSITY not specified in GRIDDATA block.'
1315  call store_error(errmsg)
1316  end if
1317  !
1318  ! -- Check for required sorption variables
1319  if (this%isrb > 0) then
1320  if (.not. lname(2)) then
1321  write (errmsg, '(a)') 'Sorption is active but BULK_DENSITY &
1322  &not specified. BULK_DENSITY must be specified in GRIDDATA block.'
1323  call store_error(errmsg)
1324  end if
1325  if (.not. lname(3)) then
1326  write (errmsg, '(a)') 'Sorption is active but distribution &
1327  &coefficient not specified. DISTCOEF must be specified in &
1328  &GRIDDATA block.'
1329  call store_error(errmsg)
1330  end if
1331  if (this%isrb > 1) then
1332  if (.not. lname(6)) then
1333  write (errmsg, '(a)') 'Freundlich or langmuir sorption is active &
1334  &but SP2 not specified. SP2 must be specified in &
1335  &GRIDDATA block.'
1336  call store_error(errmsg)
1337  end if
1338  end if
1339  else
1340  if (lname(2)) then
1341  write (warnmsg, '(a)') 'Sorption is not active but &
1342  &BULK_DENSITY was specified. BULK_DENSITY will have no affect on &
1343  &simulation results.'
1344  call store_warning(warnmsg)
1345  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1346  end if
1347  if (lname(3)) then
1348  write (warnmsg, '(a)') 'Sorption is not active but &
1349  &distribution coefficient was specified. DISTCOEF will have &
1350  &no affect on simulation results.'
1351  call store_warning(warnmsg)
1352  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1353  end if
1354  if (lname(6)) then
1355  write (warnmsg, '(a)') 'Sorption is not active but &
1356  &SP2 was specified. SP2 will have &
1357  &no affect on simulation results.'
1358  call store_warning(warnmsg)
1359  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1360  end if
1361  end if
1362  !
1363  ! -- Check for required decay/production rate coefficients
1364  if (this%idcy > 0) then
1365  if (.not. lname(4)) then
1366  write (errmsg, '(a)') 'First or zero order decay is &
1367  &active but the first rate coefficient is not specified. DECAY &
1368  &must be specified in GRIDDATA block.'
1369  call store_error(errmsg)
1370  end if
1371  if (.not. lname(5)) then
1372  !
1373  ! -- If DECAY_SORBED not specified and sorption is active, then
1374  ! terminate with an error
1375  if (this%isrb > 0) then
1376  write (errmsg, '(a)') 'DECAY_SORBED not provided in GRIDDATA &
1377  &block but decay and sorption are active. Specify DECAY_SORBED &
1378  &in GRIDDATA block.'
1379  call store_error(errmsg)
1380  end if
1381  end if
1382  else
1383  if (lname(4)) then
1384  write (warnmsg, '(a)') 'First- or zero-order decay &
1385  &is not active but decay was specified. DECAY will &
1386  &have no affect on simulation results.'
1387  call store_warning(warnmsg)
1388  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1389  end if
1390  if (lname(5)) then
1391  write (warnmsg, '(a)') 'First- or zero-order decay &
1392  &is not active but DECAY_SORBED was specified. &
1393  &DECAY_SORBED will have no affect on simulation results.'
1394  call store_warning(warnmsg)
1395  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1396  end if
1397  end if
1398  !
1399  ! -- terminate if errors
1400  if (count_errors() > 0) then
1401  call this%parser%StoreErrorUnit()
1402  end if
1403  !
1404  ! -- initialize thetam from porosity
1405  do n = 1, size(this%porosity)
1406  this%thetam(n) = this%porosity(n)
1407  end do
1408  !
1409  ! -- Return
1410  return
1411  end subroutine read_data
1412 
1413  !> @ brief Add volfrac values to volfracim
1414  !!
1415  !! Method to add immobile domain volume fracions, which are stored as a
1416  !! cumulative value in volfracim.
1417  !!
1418  !<
1419  subroutine addto_volfracim(this, volfracim)
1420  ! -- modules
1421  ! -- dummy
1422  class(gwtmsttype) :: this !< GwtMstType object
1423  real(DP), dimension(:), intent(in) :: volfracim !< immobile domain volume fraction that contributes to total immobile volume fraction
1424  ! -- local
1425  integer(I4B) :: n
1426  !
1427  ! -- Add to volfracim
1428  do n = 1, this%dis%nodes
1429  this%volfracim(n) = this%volfracim(n) + volfracim(n)
1430  end do
1431  !
1432  ! -- An immobile domain is adding a volume fraction, so update thetam
1433  ! accordingly.
1434  do n = 1, this%dis%nodes
1435  this%thetam(n) = this%get_volfracm(n) * this%porosity(n)
1436  end do
1437  !
1438  ! -- Return
1439  return
1440  end subroutine addto_volfracim
1441 
1442  !> @ brief Return mobile domain volume fraction
1443  !!
1444  !! Calculate and return the volume fraction of the aquifer that is mobile
1445  !!
1446  !<
1447  function get_volfracm(this, node) result(volfracm)
1448  ! -- modules
1449  ! -- dummy
1450  class(gwtmsttype) :: this !< GwtMstType object
1451  integer(I4B), intent(in) :: node !< node number
1452  ! -- return
1453  real(dp) :: volfracm
1454  !
1455  volfracm = done - this%volfracim(node)
1456  !
1457  ! -- Return
1458  return
1459  end function get_volfracm
1460 
1461  !> @ brief Calculate sorption concentration using Freundlich
1462  !!
1463  !! Function to calculate sorption concentration using Freundlich
1464  !!
1465  !<
1466  function get_freundlich_conc(conc, kf, a) result(cbar)
1467  ! -- dummy
1468  real(dp), intent(in) :: conc !< solute concentration
1469  real(dp), intent(in) :: kf !< freundlich constant
1470  real(dp), intent(in) :: a !< freundlich exponent
1471  ! -- return
1472  real(dp) :: cbar
1473  !
1474  if (conc > dzero) then
1475  cbar = kf * conc**a
1476  else
1477  cbar = dzero
1478  end if
1479  return
1480  end function
1481 
1482  !> @ brief Calculate sorption concentration using Langmuir
1483  !!
1484  !! Function to calculate sorption concentration using Langmuir
1485  !!
1486  !<
1487  function get_langmuir_conc(conc, kl, sbar) result(cbar)
1488  ! -- dummy
1489  real(dp), intent(in) :: conc !< solute concentration
1490  real(dp), intent(in) :: kl !< langmuir constant
1491  real(dp), intent(in) :: sbar !< langmuir sorption sites
1492  ! -- return
1493  real(dp) :: cbar
1494  !
1495  if (conc > dzero) then
1496  cbar = (kl * sbar * conc) / (done + kl * conc)
1497  else
1498  cbar = dzero
1499  end if
1500  return
1501  end function
1502 
1503  !> @ brief Calculate sorption derivative using Freundlich
1504  !!
1505  !! Function to calculate sorption derivative using Freundlich
1506  !!
1507  !<
1508  function get_freundlich_derivative(conc, kf, a) result(derv)
1509  ! -- dummy
1510  real(dp), intent(in) :: conc !< solute concentration
1511  real(dp), intent(in) :: kf !< freundlich constant
1512  real(dp), intent(in) :: a !< freundlich exponent
1513  ! -- return
1514  real(dp) :: derv
1515  !
1516  if (conc > dzero) then
1517  derv = kf * a * conc**(a - done)
1518  else
1519  derv = dzero
1520  end if
1521  return
1522  end function
1523 
1524  !> @ brief Calculate sorption derivative using Langmuir
1525  !!
1526  !! Function to calculate sorption derivative using Langmuir
1527  !!
1528  !<
1529  function get_langmuir_derivative(conc, kl, sbar) result(derv)
1530  ! -- dummy
1531  real(dp), intent(in) :: conc !< solute concentration
1532  real(dp), intent(in) :: kl !< langmuir constant
1533  real(dp), intent(in) :: sbar !< langmuir sorption sites
1534  ! -- return
1535  real(dp) :: derv
1536  !
1537  if (conc > dzero) then
1538  derv = (kl * sbar) / (done + kl * conc)**dtwo
1539  else
1540  derv = dzero
1541  end if
1542  return
1543  end function
1544 
1545  !> @ brief Calculate zero-order decay rate and constrain if necessary
1546  !!
1547  !! Function to calculate the zero-order decay rate from the user specified
1548  !! decay rate. If the decay rate is positive, then the decay rate must
1549  !! be constrained so that more mass is not removed than is available.
1550  !! Without this constraint, negative concentrations could result from
1551  !! zero-order decay.
1552  !<
1553  function get_zero_order_decay(decay_rate_usr, decay_rate_last, kiter, &
1554  cold, cnew, delt) result(decay_rate)
1555  ! -- dummy
1556  real(dp), intent(in) :: decay_rate_usr !< user-entered decay rate
1557  real(dp), intent(in) :: decay_rate_last !< decay rate used for last iteration
1558  integer(I4B), intent(in) :: kiter !< Picard iteration counter
1559  real(dp), intent(in) :: cold !< concentration at end of last time step
1560  real(dp), intent(in) :: cnew !< concentration at end of this time step
1561  real(dp), intent(in) :: delt !< length of time step
1562  ! -- return
1563  real(dp) :: decay_rate !< returned value for decay rate
1564  !
1565  ! -- Return user rate if production, otherwise constrain, if necessary
1566  if (decay_rate_usr < dzero) then
1567  !
1568  ! -- Production, no need to limit rate
1569  decay_rate = decay_rate_usr
1570  else
1571  !
1572  ! -- Need to ensure decay does not result in negative
1573  ! concentration, so reduce the rate if it would result in
1574  ! removing more mass than is in the cell.
1575  if (kiter == 1) then
1576  decay_rate = min(decay_rate_usr, cold / delt)
1577  else
1578  decay_rate = decay_rate_last
1579  if (cnew < dzero) then
1580  decay_rate = decay_rate_last + cnew / delt
1581  else if (cnew > cold) then
1582  decay_rate = decay_rate_last + cnew / delt
1583  end if
1584  decay_rate = min(decay_rate_usr, decay_rate)
1585  end if
1586  decay_rate = max(decay_rate, dzero)
1587  end if
1588  return
1589  end function get_zero_order_decay
1590 
1591 end module gwtmstmodule
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:664
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:67
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:78
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:36
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
– @ brief Mobile Storage and Transfer (MST) Module
Definition: gwt-mst.f90:10
subroutine mst_srb_term(isrb, volfracm, rhobm, vcell, tled, cnew, cold, swnew, swold, const1, const2, rate, hcofval, rhsval)
@ brief Calculate sorption terms
Definition: gwt-mst.f90:387
subroutine mst_cq_srb(this, nodes, cnew, cold, flowja)
@ brief Calculate sorption terms for package
Definition: gwt-mst.f90:713
subroutine addto_volfracim(this, volfracim)
@ brief Add volfrac values to volfracim
Definition: gwt-mst.f90:1420
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:1555
character(len=lenbudtxt), dimension(nbditems) budtxt
Definition: gwt-mst.f90:27
subroutine mst_cq(this, nodes, cnew, cold, flowja)
@ brief Calculate flows for package
Definition: gwt-mst.f90:566
subroutine mst_bd(this, isuppress_output, model_budget)
@ brief Calculate budget terms for package
Definition: gwt-mst.f90:876
subroutine allocate_arrays(this, nodes)
@ brief Allocate arrays for package
Definition: gwt-mst.f90:1052
real(dp) function get_freundlich_derivative(conc, kf, a)
@ brief Calculate sorption derivative using Freundlich
Definition: gwt-mst.f90:1509
subroutine mst_cq_dcy_srb(this, nodes, cnew, cold, flowja)
@ brief Calculate decay-sorption terms for package
Definition: gwt-mst.f90:773
real(dp) function get_freundlich_conc(conc, kf, a)
@ brief Calculate sorption concentration using Freundlich
Definition: gwt-mst.f90:1467
real(dp) function get_volfracm(this, node)
@ brief Return mobile domain volume fraction
Definition: gwt-mst.f90:1448
subroutine mst_ot_flow(this, icbcfl, icbcun)
@ brief Output flow terms for package
Definition: gwt-mst.f90:923
subroutine read_options(this)
@ brief Read options for package
Definition: gwt-mst.f90:1142
subroutine mst_ar(this, dis, ibound)
@ brief Allocate and read method for package
Definition: gwt-mst.f90:132
subroutine mst_fc_srb(this, nodes, cold, nja, matrix_sln, idxglo, rhs, cnew)
@ brief Fill sorption coefficient method for package
Definition: gwt-mst.f90:325
subroutine mst_cq_dcy(this, nodes, cnew, cold, flowja)
@ brief Calculate decay terms for package
Definition: gwt-mst.f90:655
subroutine read_data(this)
@ brief Read data for package
Definition: gwt-mst.f90:1219
subroutine allocate_scalars(this)
@ brief Allocate scalar variables for package
Definition: gwt-mst.f90:1025
subroutine mst_fc_dcy_srb(this, nodes, cold, nja, matrix_sln, idxglo, rhs, cnew, kiter)
@ brief Fill sorption-decay coefficient method for package
Definition: gwt-mst.f90:462
subroutine mst_da(this)
@ brief Deallocate
Definition: gwt-mst.f90:983
subroutine mst_fc_dcy(this, nodes, cold, cnew, nja, matrix_sln, idxglo, rhs, kiter)
@ brief Fill decay coefficient method for package
Definition: gwt-mst.f90:264
integer(i4b), parameter nbditems
Definition: gwt-mst.f90:26
subroutine mst_fc(this, nodes, cold, nja, matrix_sln, idxglo, cnew, rhs, kiter)
@ brief Fill coefficient method for package
Definition: gwt-mst.f90:170
subroutine mst_fc_sto(this, nodes, cold, nja, matrix_sln, idxglo, rhs)
@ brief Fill storage coefficient method for package
Definition: gwt-mst.f90:213
real(dp) function get_langmuir_derivative(conc, kl, sbar)
@ brief Calculate sorption derivative using Langmuir
Definition: gwt-mst.f90:1530
subroutine, public mst_cr(mstobj, name_model, inunit, iout, fmi)
@ brief Create a new package object
Definition: gwt-mst.f90:98
subroutine mst_cq_sto(this, nodes, cnew, cold, flowja)
@ brief Calculate storage terms for package
Definition: gwt-mst.f90:603
real(dp) function get_langmuir_conc(conc, kl, sbar)
@ brief Calculate sorption concentration using Langmuir
Definition: gwt-mst.f90:1488
This module defines variable data types.
Definition: kind.f90:8
This module contains the base numerical package 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
This module contains simulation variables.
Definition: SimVariables.f90:9
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
Definition: tdis.f90:29
Derived type for the Budget object.
Definition: Budget.f90:39
@ brief Mobile storage and transfer
Definition: gwt-mst.f90:37