MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwtistmodule Module Reference

– @ brief Immobile Storage and Transfer (IST) Module More...

Data Types

type  gwtisttype
 @ brief Immobile storage and transfer More...
 

Functions/Subroutines

subroutine, public ist_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, mst)
 @ brief Create a new package object More...
 
subroutine ist_ar (this)
 @ brief Allocate and read method for package More...
 
subroutine ist_rp (this)
 @ brief Read and prepare method for package More...
 
subroutine ist_ad (this)
 @ brief Advance the ist package More...
 
subroutine ist_fc (this, rhs, ia, idxglo, matrix_sln)
 @ brief Fill coefficient method for package More...
 
subroutine ist_cq (this, x, flowja, iadv)
 @ brief Calculate package flows. More...
 
subroutine ist_bd (this, model_budget)
 @ brief Add package flows to model budget. More...
 
subroutine ist_ot_model_flows (this, icbcfl, ibudfl, icbcun, imap)
 @ brief Output model flow terms. More...
 
subroutine ist_ot_dv (this, idvsave, idvprint)
 @ brief Output immobile domain concentration. More...
 
subroutine ist_ot_bdsummary (this, kstp, kper, iout, ibudfl)
 @ brief Output IST package budget summary. More...
 
subroutine ist_da (this)
 @ brief Deallocate package memory More...
 
subroutine allocate_scalars (this)
 @ brief Allocate package scalars More...
 
subroutine ist_allocate_arrays (this)
 @ brief Allocate package arrays More...
 
subroutine read_options (this)
 @ brief Read options for package More...
 
subroutine ist_read_dimensions (this)
 @ brief Read dimensions for package More...
 
subroutine read_data (this)
 @ brief Read data for package More...
 
real(dp) function get_thetaim (this, node)
 @ brief Return thetaim More...
 
subroutine get_ddterm (thetaim, vcell, delt, swtpdt, volfracim, rhobim, kd, lambda1im, lambda2im, gamma1im, gamma2im, zetaim, ddterm, f)
 @ brief Calculate immobile domain equation terms More...
 
subroutine get_hcofrhs (ddterm, f, cimold, hcof, rhs)
 @ brief Calculate the hcof and rhs terms for immobile domain More...
 
real(dp) function get_ddconc (ddterm, f, cimold, cnew)
 @ brief Calculate the immobile domain concentration More...
 
subroutine accumulate_budterm (budterm, ddterm, cimnew, cimold, cnew, idcy)
 @ brief Calculate the immobile domain budget terms More...
 

Variables

character(len=lenftype) ftype = 'IST'
 
character(len=lenpackagename) text = ' IMMOBILE DOMAIN'
 
integer(i4b), parameter nbditems = 5
 
character(len=lenbudtxt), dimension(nbditemsbudtxt
 

Detailed Description

The GwtIstModule is contains the GwtIstType, which is the derived type responsible for adding the effects of an immobile domain. In addition to representing transfer of mass between the mobile and immobile domain, the IST Package also represents the following processes within the immobile domain

  1. Changes in dissolved solute mass
  2. Decay of dissolved solute mass
  3. Sorption
  4. Decay of sorbed solute mass

Function/Subroutine Documentation

◆ accumulate_budterm()

subroutine gwtistmodule::accumulate_budterm ( real(dp), dimension(:, :), intent(inout)  budterm,
real(dp), dimension(:), intent(in)  ddterm,
real(dp), intent(in)  cimnew,
real(dp), intent(in)  cimold,
real(dp), intent(in)  cnew,
integer(i4b), intent(in)  idcy 
)
private

This subroutine calculates and accumulates the immobile domain budget terms into the budterm accumulator

Parameters
[in]ddtermterms comprising the balance equation of the immobile domain
[in]cimnewimmobile domain concenration at the end of this time step
[in]cimoldimmobile domain concentration at end of last time step
[in]cnewmobile domain concentration at the end of this time step
[in]idcyorder of decay rate (0:none, 1:first, 2:zero)

Definition at line 1273 of file gwt-ist.f90.

1274  ! -- modules
1275  ! -- dummy
1276  real(DP), dimension(:, :), intent(inout) :: budterm !<
1277  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1278  real(DP), intent(in) :: cimnew !< immobile domain concenration at the end of this time step
1279  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1280  real(DP), intent(in) :: cnew !< mobile domain concentration at the end of this time step
1281  integer(I4B), intent(in) :: idcy !< order of decay rate (0:none, 1:first, 2:zero)
1282  ! -- local
1283  real(DP) :: rate
1284  integer(I4B) :: i
1285  !
1286  ! -- calculate STORAGE-AQUEOUS
1287  i = 1
1288  rate = -ddterm(1) * cimnew + ddterm(2) * cimold
1289  if (rate > dzero) then
1290  budterm(1, i) = budterm(1, i) + rate
1291  else
1292  budterm(2, i) = budterm(2, i) - rate
1293  end if
1294  !
1295  ! -- calculate STORAGE-SORBED
1296  i = 2
1297  rate = -ddterm(3) * cimnew + ddterm(4) * cimold
1298  if (rate > dzero) then
1299  budterm(1, i) = budterm(1, i) + rate
1300  else
1301  budterm(2, i) = budterm(2, i) - rate
1302  end if
1303  !
1304  ! -- calculate DECAY-AQUEOUS
1305  i = 3
1306  rate = dzero
1307  if (idcy == 1) then
1308  rate = -ddterm(5) * cimnew
1309  else if (idcy == 2) then
1310  rate = -ddterm(7)
1311  else
1312  rate = dzero
1313  end if
1314  if (rate > dzero) then
1315  budterm(1, i) = budterm(1, i) + rate
1316  else
1317  budterm(2, i) = budterm(2, i) - rate
1318  end if
1319  !
1320  ! -- calculate DECAY-SORBED
1321  i = 4
1322  if (idcy == 1) then
1323  rate = -ddterm(6) * cimnew
1324  else if (idcy == 2) then
1325  rate = -ddterm(8)
1326  else
1327  rate = dzero
1328  end if
1329  if (rate > dzero) then
1330  budterm(1, i) = budterm(1, i) + rate
1331  else
1332  budterm(2, i) = budterm(2, i) - rate
1333  end if
1334  !
1335  ! -- calculate MOBILE-DOMAIN
1336  i = 5
1337  rate = ddterm(9) * cnew - ddterm(9) * cimnew
1338  if (rate > dzero) then
1339  budterm(1, i) = budterm(1, i) + rate
1340  else
1341  budterm(2, i) = budterm(2, i) - rate
1342  end if
1343  !
1344  !
1345  ! -- Return
1346  return
Here is the caller graph for this function:

◆ allocate_scalars()

subroutine gwtistmodule::allocate_scalars ( class(gwtisttype this)

Allocate and initialize package scalars.

Parameters
thisGwtIstType object

Definition at line 703 of file gwt-ist.f90.

704  ! -- modules
707  ! -- dummy
708  class(GwtIstType) :: this !< GwtIstType object
709  ! -- local
710  !
711  ! -- call standard BndType allocate scalars
712  call this%BndType%allocate_scalars()
713  !
714  ! -- Allocate
715  call mem_allocate(this%icimout, 'ICIMOUT', this%memoryPath)
716  call mem_allocate(this%ibudgetout, 'IBUDGETOUT', this%memoryPath)
717  call mem_allocate(this%ibudcsv, 'IBUDCSV', this%memoryPath)
718  call mem_allocate(this%isrb, 'ISRB', this%memoryPath)
719  call mem_allocate(this%idcy, 'IDCY', this%memoryPath)
720  call mem_allocate(this%kiter, 'KITER', this%memoryPath)
721  !
722  ! -- Initialize
723  this%icimout = 0
724  this%ibudgetout = 0
725  this%ibudcsv = 0
726  this%isrb = 0
727  this%idcy = 0
728  this%kiter = 0
729  !
730  ! -- Create the ocd object, which is used to manage printing and saving
731  ! of the immobile domain concentrations
732  call ocd_cr(this%ocd)
733  !
734  ! -- Return
735  return
This module contains the OutputControlDataModule.
subroutine, public ocd_cr(ocdobj)
@ brief Create OutputControlDataType
Here is the call graph for this function:

◆ get_ddconc()

real(dp) function gwtistmodule::get_ddconc ( real(dp), dimension(:), intent(in)  ddterm,
real(dp), intent(in)  f,
real(dp), intent(in)  cimold,
real(dp), intent(in)  cnew 
)
private

This function calculates the concentration of the immobile domain.

Parameters
[in]ddtermterms comprising the balance equation of the immobile domain
[in]fthe f term used to calculate the immobile domain concentration
[in]cimoldimmobile domain concentration at end of last time step
[in]cnewconcentration of the mobile domain at the end of the time step
Returns
calculated concentration of the immobile domain

Definition at line 1248 of file gwt-ist.f90.

1249  ! -- dummy
1250  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1251  real(DP), intent(in) :: f !< the f term used to calculate the immobile domain concentration
1252  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1253  real(DP), intent(in) :: cnew !< concentration of the mobile domain at the end of the time step
1254  ! -- result
1255  real(DP) :: cimnew !< calculated concentration of the immobile domain
1256  ! -- local
1257  !
1258  ! -- calculate ddconc
1259  cimnew = (ddterm(2) + ddterm(4)) * cimold + ddterm(9) * cnew - ddterm(7) &
1260  - ddterm(8)
1261  cimnew = cimnew / f
1262  !
1263  ! -- Return
1264  return
Here is the caller graph for this function:

◆ get_ddterm()

subroutine gwtistmodule::get_ddterm ( real(dp), intent(in)  thetaim,
real(dp), intent(in)  vcell,
real(dp), intent(in)  delt,
real(dp), intent(in)  swtpdt,
real(dp), intent(in)  volfracim,
real(dp), intent(in)  rhobim,
real(dp), intent(in)  kd,
real(dp), intent(in)  lambda1im,
real(dp), intent(in)  lambda2im,
real(dp), intent(in)  gamma1im,
real(dp), intent(in)  gamma2im,
real(dp), intent(in)  zetaim,
real(dp), dimension(:), intent(inout)  ddterm,
real(dp), intent(inout)  f 
)
private

This subroutine calculates the immobile domain (or dual domain) terms. The resulting ddterm and f terms are described in the GWT model report. The terms are concentration coefficients used in the balance equation for the immobile domain.

Parameters
[in]thetaimimmobile domain porosity
[in]vcellvolume of cell
[in]deltlength of time step
[in]swtpdtcell saturation at end of time step
[in]volfracimvolume fraction of immobile domain
[in]rhobimbulk density for the immobile domain (fim * rhob)
[in]kddistribution coefficient for linear isotherm
[in]lambda1imfirst-order decay rate in aqueous phase
[in]lambda2imfirst-order decay rate in sorbed phase
[in]gamma1imzero-order decay rate in aqueous phase
[in]gamma2imzero-order decay rate in sorbed phase
[in]zetaimtransfer coefficient between mobile and immobile domains
[in,out]ddtermnine terms comprising the balance equation of the immobile domain
[in,out]fthe f term used to calculate the immobile domain concentration

Definition at line 1171 of file gwt-ist.f90.

1174  ! -- dummy
1175  real(DP), intent(in) :: thetaim !< immobile domain porosity
1176  real(DP), intent(in) :: vcell !< volume of cell
1177  real(DP), intent(in) :: delt !< length of time step
1178  real(DP), intent(in) :: swtpdt !< cell saturation at end of time step
1179  real(DP), intent(in) :: volfracim !< volume fraction of immobile domain
1180  real(DP), intent(in) :: rhobim !< bulk density for the immobile domain (fim * rhob)
1181  real(DP), intent(in) :: kd !< distribution coefficient for linear isotherm
1182  real(DP), intent(in) :: lambda1im !< first-order decay rate in aqueous phase
1183  real(DP), intent(in) :: lambda2im !< first-order decay rate in sorbed phase
1184  real(DP), intent(in) :: gamma1im !< zero-order decay rate in aqueous phase
1185  real(DP), intent(in) :: gamma2im !< zero-order decay rate in sorbed phase
1186  real(DP), intent(in) :: zetaim !< transfer coefficient between mobile and immobile domains
1187  real(DP), dimension(:), intent(inout) :: ddterm !< nine terms comprising the balance equation of the immobile domain
1188  real(DP), intent(inout) :: f !< the f term used to calculate the immobile domain concentration
1189  ! -- local
1190  real(DP) :: tled
1191  !
1192  ! -- initialize
1193  tled = done / delt
1194  !
1195  ! -- Calculate terms. These terms correspond to the concentration
1196  ! coefficients in equation 7-4 of the GWT model report. However,
1197  ! an updated equation is presented as 9-9 in the supplemental technical
1198  ! information guide (mf6suptechinfo.pdf)
1199  ddterm(1) = thetaim * vcell * tled
1200  ddterm(2) = thetaim * vcell * tled
1201  ddterm(3) = volfracim * rhobim * vcell * kd * tled
1202  ddterm(4) = volfracim * rhobim * vcell * kd * tled
1203  ddterm(5) = thetaim * lambda1im * vcell
1204  ddterm(6) = lambda2im * volfracim * rhobim * kd * vcell
1205  ddterm(7) = thetaim * gamma1im * vcell
1206  ddterm(8) = gamma2im * volfracim * rhobim * vcell
1207  ddterm(9) = vcell * swtpdt * zetaim
1208  !
1209  ! -- calculate denominator term, f
1210  f = ddterm(1) + ddterm(3) + ddterm(5) + ddterm(6) + ddterm(9)
1211  !
1212  ! -- Return
1213  return
Here is the caller graph for this function:

◆ get_hcofrhs()

subroutine gwtistmodule::get_hcofrhs ( real(dp), dimension(:), intent(in)  ddterm,
real(dp), intent(in)  f,
real(dp), intent(in)  cimold,
real(dp), intent(inout)  hcof,
real(dp), intent(inout)  rhs 
)
private

This subroutine calculates the hcof and rhs terms that must be added to the solution system of equations

Parameters
[in]ddtermterms comprising the balance equation of the immobile domain
[in]fthe f term used to calculate the immobile domain concentration
[in]cimoldimmobile domain concentration at end of last time step
[in,out]hcofcalculated contribution for the a-matrix diagonal position
[in,out]rhscalculated contribution for the solution right-hand side

Definition at line 1222 of file gwt-ist.f90.

1223  ! -- dummy
1224  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1225  real(DP), intent(in) :: f !< the f term used to calculate the immobile domain concentration
1226  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1227  real(DP), intent(inout) :: hcof !< calculated contribution for the a-matrix diagonal position
1228  real(DP), intent(inout) :: rhs !< calculated contribution for the solution right-hand side
1229  !
1230  ! -- calculate hcof
1231  hcof = ddterm(9)**2 / f - ddterm(9)
1232  !
1233  ! -- calculate rhs, and switch the sign because this term needs to
1234  ! be moved to the left hand side
1235  rhs = (ddterm(2) + ddterm(4)) * cimold - ddterm(7) - ddterm(8)
1236  rhs = rhs * ddterm(9) / f
1237  rhs = -rhs
1238  !
1239  ! -- Return
1240  return
Here is the caller graph for this function:

◆ get_thetaim()

real(dp) function gwtistmodule::get_thetaim ( class(gwtisttype this,
integer(i4b), intent(in)  node 
)

Calculate and return thetaim, volume of immobile voids per aquifer volume

Parameters
thisGwtIstType object
[in]nodenode number

Definition at line 1149 of file gwt-ist.f90.

1150  ! -- modules
1151  ! -- dummy
1152  class(GwtIstType) :: this !< GwtIstType object
1153  integer(I4B), intent(in) :: node !< node number
1154  ! -- return
1155  real(DP) :: thetaim
1156  !
1157  thetaim = this%volfrac(node) * this%porosity(node)
1158  !
1159  ! -- Return
1160  return

◆ ist_ad()

subroutine gwtistmodule::ist_ad ( class(gwtisttype this)
private

Advance the IST Package and handle the adaptive time stepping feature by copying from new to old or old to new accordingly

Parameters
thisGwtIstType object

Definition at line 233 of file gwt-ist.f90.

234  ! -- modules
236  ! -- dummy variables
237  class(GwtIstType) :: this !< GwtIstType object
238  ! -- local variables
239  integer(I4B) :: n
240  !
241  ! -- Call parent advance
242  call this%BndType%bnd_ad()
243  !
244  ! -- set independent kiter counter to zero
245  this%kiter = 0
246  !
247  ! -- copy cimnew into cimold or vice versa if this is a repeat of
248  ! a failed time step
249  if (ifailedstepretry == 0) then
250  do n = 1, this%dis%nodes
251  this%cimold(n) = this%cimnew(n)
252  end do
253  else
254  do n = 1, this%dis%nodes
255  this%cimnew(n) = this%cimold(n)
256  end do
257  end if
258  !
259  return
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) ifailedstepretry
current retry for this time step

◆ ist_allocate_arrays()

subroutine gwtistmodule::ist_allocate_arrays ( class(gwtisttype), intent(inout)  this)

Allocate and initialize package arrays.

Parameters
[in,out]thisGwtIstType object

Definition at line 743 of file gwt-ist.f90.

744  ! -- modules
746  ! -- dummy
747  class(GwtIstType), intent(inout) :: this !< GwtIstType object
748  ! -- local
749  integer(I4B) :: n
750  !
751  ! -- call standard BndType allocate scalars
752  ! nbound and maxbound are 0 in order to keep memory footprint low
753  call this%BndType%allocate_arrays()
754  !
755  ! -- allocate ist arrays of size nodes
756  call mem_allocate(this%strg, this%dis%nodes, 'STRG', this%memoryPath)
757  call mem_allocate(this%cim, this%dis%nodes, 'CIM', this%memoryPath)
758  call mem_allocate(this%cimnew, this%dis%nodes, 'CIMNEW', this%memoryPath)
759  call mem_allocate(this%cimold, this%dis%nodes, 'CIMOLD', this%memoryPath)
760  call mem_allocate(this%porosity, this%dis%nodes, 'POROSITY', this%memoryPath)
761  call mem_allocate(this%zetaim, this%dis%nodes, 'ZETAIM', this%memoryPath)
762  call mem_allocate(this%volfrac, this%dis%nodes, 'VOLFRAC', this%memoryPath)
763  if (this%isrb == 0) then
764  call mem_allocate(this%bulk_density, 1, 'BULK_DENSITY', this%memoryPath)
765  call mem_allocate(this%distcoef, 1, 'DISTCOEF', this%memoryPath)
766  else
767  call mem_allocate(this%bulk_density, this%dis%nodes, 'BULK_DENSITY', &
768  this%memoryPath)
769  call mem_allocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
770  this%memoryPath)
771  end if
772  if (this%idcy == 0) then
773  call mem_allocate(this%decay, 1, 'DECAY', this%memoryPath)
774  call mem_allocate(this%decaylast, 1, 'DECAYLAST', this%memoryPath)
775  else
776  call mem_allocate(this%decay, this%dis%nodes, 'DECAY', this%memoryPath)
777  call mem_allocate(this%decaylast, this%dis%nodes, 'DECAYLAST', &
778  this%memoryPath)
779  end if
780  if (this%isrb == 0 .and. this%idcy == 0) then
781  call mem_allocate(this%decayslast, 1, 'DECAYSLAST', this%memoryPath)
782  else
783  call mem_allocate(this%decayslast, this%dis%nodes, 'DECAYSLAST', &
784  this%memoryPath)
785  end if
786  call mem_allocate(this%decay_sorbed, 1, 'DECAY_SORBED', this%memoryPath)
787  !
788  ! -- initialize
789  do n = 1, this%dis%nodes
790  this%strg(n) = dzero
791  this%cim(n) = dzero
792  this%cimnew(n) = dzero
793  this%cimold(n) = dzero
794  this%porosity(n) = dzero
795  this%zetaim(n) = dzero
796  this%volfrac(n) = dzero
797  end do
798  do n = 1, size(this%decay)
799  this%decay(n) = dzero
800  this%decaylast(n) = dzero
801  end do
802  do n = 1, size(this%decayslast)
803  this%decayslast(n) = dzero
804  end do
805  !
806  ! -- Set pointers
807  this%ocd%dis => this%dis
808  !
809  ! -- return
810  return

◆ ist_ar()

subroutine gwtistmodule::ist_ar ( class(gwtisttype), intent(inout)  this)
private

Method to allocate and read static data for the package.

Parameters
[in,out]thisGwtIstType object

Definition at line 157 of file gwt-ist.f90.

158  ! -- modules
160  use budgetmodule, only: budget_cr
161  ! -- dummy
162  class(GwtIstType), intent(inout) :: this !< GwtIstType object
163  ! -- local
164  integer(I4B) :: n
165  !
166  ! -- Allocate arrays
167  call this%ist_allocate_arrays()
168  !
169  ! -- Now that arrays are allocated, check in the cimnew array to
170  ! the output control manager for subsequent printing/saving
171  call this%ocd%init_dbl('CIM', this%cimnew, this%dis, 'PRINT LAST ', &
172  'COLUMNS 10 WIDTH 11 DIGITS 4 GENERAL ', &
173  this%iout, dhnoflo)
174  !
175  ! -- read the data block
176  call this%read_data()
177  !
178  ! -- set cimnew to the cim start values read from input
179  do n = 1, this%dis%nodes
180  this%cimnew(n) = this%cim(n)
181  end do
182  !
183  ! -- add volfrac to the volfracim accumulator in mst package
184  call this%mst%addto_volfracim(this%volfrac)
185  !
186  ! -- setup the immobile domain budget
187  call budget_cr(this%budget, this%memoryPath)
188  call this%budget%budget_df(nbditems, 'MASS', 'M', bdzone=this%packName)
189  call this%budget%set_ibudcsv(this%ibudcsv)
190  !
191  ! -- Perform a check to ensure that sorption and decay are set
192  ! consistently between the MST and IST packages.
193  if (this%idcy /= this%mst%idcy) then
194  call store_error('DECAY must be activated consistently between the &
195  &MST and IST Packages. Activate or deactivate DECAY for both &
196  &Packages.')
197  end if
198  if (this%isrb /= this%mst%isrb) then
199  call store_error('Sorption is active for the IST Package but it is not &
200  &compatible with the sorption option selected for the MST Package. &
201  &If sorption is active for the IST Package, then SORPTION LINEAR must &
202  &be specified in the options block of the MST Package.')
203  end if
204  if (count_errors() > 0) then
205  call this%parser%StoreErrorUnit()
206  end if
207  !
208  ! -- Return
209  return
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
This module contains simulation methods.
Definition: Sim.f90:10
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
Here is the call graph for this function:

◆ ist_bd()

subroutine gwtistmodule::ist_bd ( class(gwtisttype this,
type(budgettype), intent(inout)  model_budget 
)

Add the flow between IST package and the model (ratin and ratout) to the model budget.

Parameters
thisGwtIstType object
[in,out]model_budgetmodel budget object

Definition at line 494 of file gwt-ist.f90.

495  ! -- modules
496  use tdismodule, only: delt
498  ! -- dummy
499  class(GwtIstType) :: this !< GwtIstType object
500  type(BudgetType), intent(inout) :: model_budget !< model budget object
501  ! -- local
502  real(DP) :: ratin
503  real(DP) :: ratout
504  integer(I4B) :: isuppress_output
505  isuppress_output = 0
506  call rate_accumulator(this%strg(:), ratin, ratout)
507  call model_budget%addentry(ratin, ratout, delt, this%text, &
508  isuppress_output, this%packName)
509  return
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:664
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
Here is the call graph for this function:

◆ ist_cq()

subroutine gwtistmodule::ist_cq ( class(gwtisttype), intent(inout)  this,
real(dp), dimension(:), intent(in)  x,
real(dp), dimension(:), intent(inout), contiguous  flowja,
integer(i4b), intent(in), optional  iadv 
)

Calculate the flow between connected package control volumes.

Parameters
[in,out]thisGwtIstType object
[in]xcurrent dependent-variable value
[in,out]flowjaflow between two connected control volumes
[in]iadvflag that indicates if this is an advance package

Definition at line 374 of file gwt-ist.f90.

375  ! -- modules
376  use tdismodule, only: delt
377  use constantsmodule, only: dzero
378  ! -- dummy
379  class(GwtIstType), intent(inout) :: this !< GwtIstType object
380  real(DP), dimension(:), intent(in) :: x !< current dependent-variable value
381  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
382  integer(I4B), optional, intent(in) :: iadv !< flag that indicates if this is an advance package
383  ! -- local
384  integer(I4B) :: idiag
385  integer(I4B) :: n
386  real(DP) :: rate
387  real(DP) :: swt, swtpdt
388  real(DP) :: hhcof, rrhs
389  real(DP) :: vcell
390  real(DP) :: thetaim
391  real(DP) :: zetaim
392  real(DP) :: kd
393  real(DP) :: volfracim
394  real(DP) :: rhobim
395  real(DP) :: lambda1im
396  real(DP) :: lambda2im
397  real(DP) :: gamma1im
398  real(DP) :: gamma2im
399  real(DP) :: cimnew
400  real(DP) :: cimold
401  real(DP) :: f
402  real(DP) :: cimsrbold
403  real(DP) :: cimsrbnew
404  real(DP), dimension(9) :: ddterm
405  ! -- formats
406  !
407  ! -- initialize
408  this%budterm(:, :) = dzero
409  !
410  ! -- Calculate immobile domain transfer rate
411  do n = 1, this%dis%nodes
412  !
413  ! -- skip if transport inactive
414  rate = dzero
415  cimnew = dzero
416  if (this%ibound(n) > 0) then
417  !
418  ! -- calculate new and old water volumes
419  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
420  swtpdt = this%fmi%gwfsat(n)
421  swt = this%fmi%gwfsatold(n, delt)
422  thetaim = this%get_thetaim(n)
423  !
424  ! -- set exchange coefficient
425  zetaim = this%zetaim(n)
426  !
427  ! -- Calculate exchange with immobile domain
428  rate = dzero
429  hhcof = dzero
430  rrhs = dzero
431  kd = dzero
432  volfracim = dzero
433  rhobim = dzero
434  lambda1im = dzero
435  lambda2im = dzero
436  gamma1im = dzero
437  gamma2im = dzero
438  if (this%idcy == 1) lambda1im = this%decay(n)
439  if (this%idcy == 2) then
440  gamma1im = get_zero_order_decay(this%decay(n), this%decaylast(n), 0, &
441  this%cimold(n), this%cimnew(n), delt)
442  end if
443  if (this%isrb > 0) then
444  kd = this%distcoef(n)
445  volfracim = this%volfrac(n)
446  rhobim = this%bulk_density(n)
447  if (this%idcy == 1) lambda2im = this%decay_sorbed(n)
448  if (this%idcy == 2) then
449  cimsrbold = this%cimold(n) * kd
450  cimsrbnew = this%cimnew(n) * kd
451  gamma2im = get_zero_order_decay(this%decay_sorbed(n), &
452  this%decayslast(n), &
453  0, cimsrbold, &
454  cimsrbnew, delt)
455  end if
456  end if
457  !
458  ! -- calculate the terms and then get the hcof and rhs contributions
459  call get_ddterm(thetaim, vcell, delt, swtpdt, &
460  volfracim, rhobim, kd, lambda1im, lambda2im, &
461  gamma1im, gamma2im, zetaim, ddterm, f)
462  cimold = this%cimold(n)
463  call get_hcofrhs(ddterm, f, cimold, hhcof, rrhs)
464  !
465  ! -- calculate rate from hcof and rhs
466  rate = hhcof * x(n) - rrhs
467  !
468  ! -- calculate immobile domain concentration
469  cimnew = get_ddconc(ddterm, f, cimold, x(n))
470  !
471  ! -- accumulate the budget terms
472  call accumulate_budterm(this%budterm, ddterm, cimnew, cimold, x(n), &
473  this%idcy)
474  end if
475  !
476  ! -- store rate and add to flowja
477  this%strg(n) = rate
478  idiag = this%dis%con%ia(n)
479  flowja(idiag) = flowja(idiag) + rate
480  !
481  ! -- store immobile domain concentration
482  this%cimnew(n) = cimnew
483  !
484  end do
485  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
Here is the call graph for this function:

◆ ist_create()

subroutine, public gwtistmodule::ist_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=*), intent(in)  pakname,
type(tspfmitype), pointer  fmi,
type(gwtmsttype), pointer  mst 
)

Create a new IST object

Parameters
packobjBndType pointer that will point to new IST Package
[in]idname of the model
[in]ibcnumconsecutive package number
[in]inunitunit number of package input file
[in]ioutunit number of model listing file
[in]namemodelname of the model
[in]paknamename of the package

Definition at line 107 of file gwt-ist.f90.

109  ! -- dummy
110  class(BndType), pointer :: packobj !< BndType pointer that will point to new IST Package
111  integer(I4B), intent(in) :: id !< name of the model
112  integer(I4B), intent(in) :: ibcnum !< consecutive package number
113  integer(I4B), intent(in) :: inunit !< unit number of package input file
114  integer(I4B), intent(in) :: iout !< unit number of model listing file
115  character(len=*), intent(in) :: namemodel !< name of the model
116  character(len=*), intent(in) :: pakname !< name of the package
117  ! -- local
118  type(GwtIstType), pointer :: istobj
119  type(TspFmiType), pointer :: fmi
120  type(GwtMstType), pointer :: mst
121  !
122  ! -- allocate the object and assign values to object variables
123  allocate (istobj)
124  packobj => istobj
125  !
126  ! -- create name and memory path
127  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
128  packobj%text = text
129  !
130  ! -- allocate scalars
131  call packobj%allocate_scalars()
132  !
133  ! -- initialize package
134  call packobj%pack_initialize()
135  !
136  ! -- store values
137  packobj%inunit = inunit
138  packobj%iout = iout
139  packobj%id = id
140  packobj%ibcnum = ibcnum
141  packobj%ncolbnd = 1
142  packobj%iscloc = 1
143  !
144  ! -- Point IST specific variables
145  istobj%fmi => fmi
146  istobj%mst => mst
147  !
148  ! -- return
149  return
Here is the caller graph for this function:

◆ ist_da()

subroutine gwtistmodule::ist_da ( class(gwtisttype this)

Deallocate package scalars and arrays.

Parameters
thisGwtIstType object

Definition at line 652 of file gwt-ist.f90.

653  ! -- modules
655  ! -- dummy
656  class(GwtIstType) :: this !< GwtIstType object
657  !
658  ! -- Deallocate arrays if package was active
659  if (this%inunit > 0) then
660  call mem_deallocate(this%icimout)
661  call mem_deallocate(this%ibudgetout)
662  call mem_deallocate(this%ibudcsv)
663  call mem_deallocate(this%idcy)
664  call mem_deallocate(this%isrb)
665  call mem_deallocate(this%kiter)
666  call mem_deallocate(this%cim)
667  call mem_deallocate(this%cimnew)
668  call mem_deallocate(this%cimold)
669  call mem_deallocate(this%zetaim)
670  call mem_deallocate(this%porosity)
671  call mem_deallocate(this%volfrac)
672  call mem_deallocate(this%bulk_density)
673  call mem_deallocate(this%distcoef)
674  call mem_deallocate(this%decay)
675  call mem_deallocate(this%decaylast)
676  call mem_deallocate(this%decayslast)
677  call mem_deallocate(this%decay_sorbed)
678  call mem_deallocate(this%strg)
679  this%fmi => null()
680  this%mst => null()
681  end if
682  !
683  ! -- Scalars
684  !
685  ! -- Objects
686  call this%budget%budget_da()
687  deallocate (this%budget)
688  call this%ocd%ocd_da()
689  deallocate (this%ocd)
690  !
691  ! -- deallocate parent
692  call this%BndType%bnd_da()
693  !
694  ! -- Return
695  return

◆ ist_fc()

subroutine gwtistmodule::ist_fc ( class(gwtisttype this,
real(dp), dimension(:), intent(inout)  rhs,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), dimension(:), intent(in)  idxglo,
class(matrixbasetype), pointer  matrix_sln 
)

Method to calculate and fill coefficients for the package.

Parameters
thisGwtIstType object
[in,out]rhsright-hand side vector for model
[in]iasolution CRS row pointers
[in]idxglomapping vector for model (local) to solution (global)
matrix_slnsolution coefficient matrix

Definition at line 267 of file gwt-ist.f90.

268  ! -- modules
269  use tdismodule, only: delt
270  ! -- dummy
271  class(GwtIstType) :: this !< GwtIstType object
272  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
273  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
274  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
275  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
276  ! -- local
277  integer(I4B) :: n, idiag
278  real(DP) :: tled
279  real(DP) :: hhcof, rrhs
280  real(DP) :: swt, swtpdt
281  real(DP) :: vcell
282  real(DP) :: thetaim
283  real(DP) :: zetaim
284  real(DP) :: kd
285  real(DP) :: volfracim
286  real(DP) :: rhobim
287  real(DP) :: lambda1im
288  real(DP) :: lambda2im
289  real(DP) :: gamma1im
290  real(DP) :: gamma2im
291  real(DP) :: cimold
292  real(DP) :: f
293  real(DP) :: cimsrbold
294  real(DP) :: cimsrbnew
295  real(DP), dimension(9) :: ddterm
296  !
297  ! -- set variables
298  tled = done / delt
299  this%kiter = this%kiter + 1
300  !
301  ! -- loop through and calculate immobile domain contribution to hcof and rhs
302  do n = 1, this%dis%nodes
303  !
304  ! -- skip if transport inactive
305  if (this%ibound(n) <= 0) cycle
306  !
307  ! -- calculate new and old water volumes
308  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
309  swtpdt = this%fmi%gwfsat(n)
310  swt = this%fmi%gwfsatold(n, delt)
311  thetaim = this%get_thetaim(n)
312  idiag = ia(n)
313  !
314  ! -- set exchange coefficient
315  zetaim = this%zetaim(n)
316  !
317  ! -- Add dual domain mass transfer contributions to rhs and hcof
318  kd = dzero
319  volfracim = dzero
320  rhobim = dzero
321  lambda1im = dzero
322  lambda2im = dzero
323  gamma1im = dzero
324  gamma2im = dzero
325  !
326  ! -- setup decay variables
327  if (this%idcy == 1) lambda1im = this%decay(n)
328  if (this%idcy == 2) then
329  gamma1im = get_zero_order_decay(this%decay(n), this%decaylast(n), &
330  this%kiter, this%cimold(n), &
331  this%cimnew(n), delt)
332  this%decaylast(n) = gamma1im
333  end if
334  !
335  ! -- setup sorption variables
336  if (this%isrb > 0) then
337  kd = this%distcoef(n)
338  volfracim = this%volfrac(n)
339  rhobim = this%bulk_density(n)
340  if (this%idcy == 1) lambda2im = this%decay_sorbed(n)
341  if (this%idcy == 2) then
342  cimsrbold = this%cimold(n) * kd
343  cimsrbnew = this%cimnew(n) * kd
344  gamma2im = get_zero_order_decay(this%decay_sorbed(n), &
345  this%decayslast(n), &
346  this%kiter, cimsrbold, &
347  cimsrbnew, delt)
348  this%decayslast(n) = gamma2im
349  end if
350  end if
351  !
352  ! -- calculate the terms and then get the hcof and rhs contributions
353  call get_ddterm(thetaim, vcell, delt, swtpdt, &
354  volfracim, rhobim, kd, lambda1im, lambda2im, &
355  gamma1im, gamma2im, zetaim, ddterm, f)
356  cimold = this%cimold(n)
357  call get_hcofrhs(ddterm, f, cimold, hhcof, rrhs)
358  !
359  ! -- update solution accumulators
360  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
361  rhs(n) = rhs(n) + rrhs
362  !
363  end do
364  !
365  ! -- Return
366  return
Here is the call graph for this function:

◆ ist_ot_bdsummary()

subroutine gwtistmodule::ist_ot_bdsummary ( class(gwtisttype this,
integer(i4b), intent(in)  kstp,
integer(i4b), intent(in)  kper,
integer(i4b), intent(in)  iout,
integer(i4b), intent(in)  ibudfl 
)

Output advanced boundary package budget summary. This method only needs to be overridden for advanced packages that save budget summaries to the model listing file.

Parameters
thisGwtIstType object
[in]kstptime step number
[in]kperperiod number
[in]ioutflag and unit number for the model listing file
[in]ibudflflag indicating budget should be written

Definition at line 620 of file gwt-ist.f90.

621  ! -- modules
622  use tdismodule, only: delt, totim
623  ! -- dummy variables
624  class(GwtIstType) :: this !< GwtIstType object
625  integer(I4B), intent(in) :: kstp !< time step number
626  integer(I4B), intent(in) :: kper !< period number
627  integer(I4B), intent(in) :: iout !< flag and unit number for the model listing file
628  integer(I4B), intent(in) :: ibudfl !< flag indicating budget should be written
629  ! -- local
630  integer(I4B) :: isuppress_output = 0
631  !
632  ! -- Fill budget terms
633  call this%budget%reset()
634  call this%budget%addentry(this%budterm, delt, budtxt, isuppress_output)
635  !
636  ! -- Write budget to list file
637  call this%budget%finalize_step(delt)
638  if (ibudfl /= 0) then
639  call this%budget%budget_ot(kstp, kper, iout)
640  end if
641  !
642  ! -- Write budget csv
643  call this%budget%writecsv(totim)
644  return
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32

◆ ist_ot_dv()

subroutine gwtistmodule::ist_ot_dv ( class(gwtisttype this,
integer(i4b), intent(in)  idvsave,
integer(i4b), intent(in)  idvprint 
)
Parameters
thisBndType object
[in]idvsaveflag and unit number for dependent-variable output
[in]idvprintflag indicating if dependent-variable should be written to the model listing file

Definition at line 584 of file gwt-ist.f90.

585  ! -- modules
586  use tdismodule, only: kstp, endofperiod
587  ! -- dummy variables
588  class(GwtIstType) :: this !< BndType object
589  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
590  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
591  ! -- local
592  integer(I4B) :: ipflg
593  integer(I4B) :: ibinun
594  !
595  ! -- Save cim to a binary file. ibinun is a flag where 1 indicates that
596  ! cim should be written to a binary file if a binary file is open
597  ! for it.
598  ipflg = 0
599  ibinun = 1
600  if (idvsave == 0) ibinun = 0
601  if (ibinun /= 0) then
602  call this%ocd%ocd_ot(ipflg, kstp, endofperiod, this%iout, &
603  iprint_opt=0, isav_opt=ibinun)
604  end if
605  !
606  ! -- Print immobile domain concentrations to listing file
607  if (idvprint /= 0) then
608  call this%ocd%ocd_ot(ipflg, kstp, endofperiod, this%iout, &
609  iprint_opt=idvprint, isav_opt=0)
610  end if
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24

◆ ist_ot_model_flows()

subroutine gwtistmodule::ist_ot_model_flows ( class(gwtisttype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(in)  icbcun,
integer(i4b), dimension(:), intent(in), optional  imap 
)

Output flow terms between the IST package and model to a binary file and/or print flows to the model listing file.

Parameters
thisGwtIstType object
[in]icbcflflag for cell-by-cell output
[in]ibudflflag indication if cell-by-cell data should be saved
[in]icbcununit number for cell-by-cell output
[in]imapmapping vector

Definition at line 518 of file gwt-ist.f90.

519  ! -- modules
520  use constantsmodule, only: dzero
521  ! -- dummy
522  class(GwtIstType) :: this !< GwtIstType object
523  integer(I4B), intent(in) :: icbcfl !< flag for cell-by-cell output
524  integer(I4B), intent(in) :: ibudfl !< flag indication if cell-by-cell data should be saved
525  integer(I4B), intent(in) :: icbcun !< unit number for cell-by-cell output
526  integer(I4B), dimension(:), optional, intent(in) :: imap !< mapping vector
527  ! -- loca
528  integer(I4B) :: n
529  integer(I4B) :: ibinun
530  integer(I4B) :: nbound
531  integer(I4B) :: naux
532  real(DP) :: rate
533  !
534  ! -- Set unit number for binary output
535  if (this%ipakcb < 0) then
536  ibinun = icbcun
537  elseif (this%ipakcb == 0) then
538  ibinun = 0
539  else
540  ibinun = this%ipakcb
541  end if
542  if (icbcfl == 0) ibinun = 0
543  !
544  ! -- Record the storage rate if requested
545  !
546  ! -- If cell-by-cell flows will be saved as a list, write header.
547  if (ibinun /= 0) then
548  nbound = this%dis%nodes
549  naux = 0
550  call this%dis%record_srcdst_list_header(this%text, this%name_model, &
551  this%name_model, this%name_model, &
552  this%packName, naux, this%auxname, &
553  ibinun, nbound, this%iout)
554  end if
555  !
556  ! -- Calculate immobile domain rhs and hcof
557  do n = 1, this%dis%nodes
558  !
559  ! -- skip if transport inactive
560  rate = dzero
561  if (this%ibound(n) > 0) then
562  !
563  ! -- set rate from this%strg
564  rate = this%strg(n)
565  end if
566  !
567  ! -- If saving cell-by-cell flows in list, write flow
568  if (ibinun /= 0) then
569  call this%dis%record_mf6_list_entry(ibinun, n, n, rate, &
570  naux, this%auxvar(:, n), &
571  olconv=.true., &
572  olconv2=.true.)
573  end if
574  !
575  end do
576  !
577  ! -- Return
578  return

◆ ist_read_dimensions()

subroutine gwtistmodule::ist_read_dimensions ( class(gwtisttype), intent(inout)  this)

Read dimensions for package.

Parameters
[in,out]thisGwtIstType object

Definition at line 922 of file gwt-ist.f90.

923  ! -- dummy
924  class(GwtIstType), intent(inout) :: this !< GwtIstType object
925  ! -- local
926  ! -- format
927  !
928  ! -- return
929  return

◆ ist_rp()

subroutine gwtistmodule::ist_rp ( class(gwtisttype), intent(inout)  this)

Method to read and prepare package data

Parameters
[in,out]thisGwtIstType object

Definition at line 217 of file gwt-ist.f90.

218  ! -- dummy
219  class(GwtIstType), intent(inout) :: this !< GwtIstType object
220  ! -- local
221  ! -- format
222  !
223  ! -- return
224  return

◆ read_data()

subroutine gwtistmodule::read_data ( class(gwtisttype this)
private

Read data for package.

Parameters
thisGwtIstType object

Definition at line 937 of file gwt-ist.f90.

938  ! -- modules
939  use constantsmodule, only: linelength
942  ! -- dummy
943  class(GwtIstType) :: this !< GwtIstType object
944  ! -- local
945  character(len=LINELENGTH) :: errmsg, keyword
946  character(len=:), allocatable :: line
947  integer(I4B) :: istart, istop, lloc, ierr
948  logical :: isfound, endOfBlock
949  logical, dimension(8) :: lname
950  character(len=24), dimension(8) :: aname
951  ! -- formats
952  ! -- data
953  data aname(1)/' BULK DENSITY'/
954  data aname(2)/'DISTRIBUTION COEFFICIENT'/
955  data aname(3)/' DECAY RATE'/
956  data aname(4)/' DECAY SORBED RATE'/
957  data aname(5)/' INITIAL IMMOBILE CONC'/
958  data aname(6)/' FIRST ORDER TRANS RATE'/
959  data aname(7)/'IMMOBILE DOMAIN POROSITY'/
960  data aname(8)/'IMMOBILE VOLUME FRACTION'/
961  !
962  ! -- initialize
963  isfound = .false.
964  lname(:) = .false.
965  !
966  ! -- get griddata block
967  call this%parser%GetBlock('GRIDDATA', isfound, ierr)
968  if (isfound) then
969  write (this%iout, '(1x,a)') 'PROCESSING GRIDDATA'
970  do
971  call this%parser%GetNextLine(endofblock)
972  if (endofblock) exit
973  call this%parser%GetStringCaps(keyword)
974  call this%parser%GetRemainingLine(line)
975  lloc = 1
976  select case (keyword)
977  case ('BULK_DENSITY')
978  if (this%isrb == 0) &
979  call mem_reallocate(this%bulk_density, this%dis%nodes, &
980  'BULK_DENSITY', trim(this%memoryPath))
981  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
982  this%parser%iuactive, &
983  this%bulk_density, aname(1))
984  lname(1) = .true.
985  case ('DISTCOEF')
986  if (this%isrb == 0) &
987  call mem_reallocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
988  trim(this%memoryPath))
989  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
990  this%parser%iuactive, this%distcoef, &
991  aname(2))
992  lname(2) = .true.
993  case ('DECAY')
994  if (this%idcy == 0) &
995  call mem_reallocate(this%decay, this%dis%nodes, 'DECAY', &
996  trim(this%memoryPath))
997  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
998  this%parser%iuactive, this%decay, &
999  aname(3))
1000  lname(3) = .true.
1001  case ('DECAY_SORBED')
1002  call mem_reallocate(this%decay_sorbed, this%dis%nodes, &
1003  'DECAY_SORBED', trim(this%memoryPath))
1004  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1005  this%parser%iuactive, &
1006  this%decay_sorbed, aname(4))
1007  lname(4) = .true.
1008  case ('CIM')
1009  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1010  this%parser%iuactive, this%cim, &
1011  aname(5))
1012  lname(5) = .true.
1013  case ('ZETAIM')
1014  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1015  this%parser%iuactive, this%zetaim, &
1016  aname(6))
1017  lname(6) = .true.
1018  case ('POROSITY')
1019  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1020  this%parser%iuactive, this%porosity, &
1021  aname(7))
1022  lname(7) = .true.
1023  case ('VOLFRAC')
1024  call this%dis%read_grid_array(line, lloc, istart, istop, this%iout, &
1025  this%parser%iuactive, this%volfrac, &
1026  aname(8))
1027  lname(8) = .true.
1028  case ('THETAIM')
1029  write (errmsg, '(a)') &
1030  'THETAIM is no longer supported. See Chapter 9 in &
1031  &mf6suptechinfo.pdf for revised parameterization of mobile and &
1032  &immobile domain simulations.'
1033  call store_error(errmsg)
1034  call this%parser%StoreErrorUnit()
1035  case default
1036  write (errmsg, '(a,a)') 'Unknown GRIDDATA tag: ', trim(keyword)
1037  call store_error(errmsg)
1038  call this%parser%StoreErrorUnit()
1039  end select
1040  end do
1041  write (this%iout, '(1x,a)') 'END PROCESSING GRIDDATA'
1042  else
1043  write (errmsg, '(a)') 'Required GRIDDATA block not found.'
1044  call store_error(errmsg)
1045  call this%parser%StoreErrorUnit()
1046  end if
1047  !
1048  ! -- Check for required sorption variables
1049  if (this%isrb > 0) then
1050  if (.not. lname(1)) then
1051  write (errmsg, '(a)') 'Sorption is active but BULK_DENSITY &
1052  &not specified. BULK_DENSITY must be specified in griddata block.'
1053  call store_error(errmsg)
1054  end if
1055  if (.not. lname(2)) then
1056  write (errmsg, '(a)') 'Sorption is active but distribution &
1057  &coefficient not specified. DISTCOEF must be specified in &
1058  &GRIDDATA block.'
1059  call store_error(errmsg)
1060  end if
1061  else
1062  if (lname(1)) then
1063  write (this%iout, '(1x,a)') 'Warning. Sorption is not active but &
1064  &BULK_DENSITY was specified. BULK_DENSITY will have no affect on &
1065  &simulation results.'
1066  end if
1067  if (lname(2)) then
1068  write (this%iout, '(1x,a)') 'Warning. Sorption is not active but &
1069  &distribution coefficient was specified. DISTCOEF will have &
1070  &no affect on simulation results.'
1071  end if
1072  end if
1073  !
1074  ! -- Check for required decay/production rate coefficients
1075  if (this%idcy > 0) then
1076  if (.not. lname(3)) then
1077  write (errmsg, '(a)') 'First- or zero-order decay is &
1078  &active but the first rate coefficient was not specified. &
1079  &Decay must be specified in GRIDDATA block.'
1080  call store_error(errmsg)
1081  end if
1082  if (.not. lname(4)) then
1083  !
1084  ! -- If DECAY_SORBED not specified and sorption is active, then set
1085  ! decay_sorbed equal to decay
1086  if (this%isrb > 0) then
1087  write (errmsg, '(a)') 'DECAY_SORBED not provided in GRIDDATA &
1088  &block but decay and sorption are active. Specify DECAY_SORBED &
1089  &in GRIDDATA block.'
1090  call store_error(errmsg)
1091  end if
1092  end if
1093  else
1094  if (lname(3)) then
1095  write (this%iout, '(1x,a)') 'Warning. First- or zero-order decay &
1096  &is not active but DECAY was specified. DECAY will &
1097  &have no affect on simulation results.'
1098  end if
1099  if (lname(4)) then
1100  write (this%iout, '(1x,a)') 'Warning. First- or zero-order decay &
1101  &is not active but DECAY_SORBED was specified. &
1102  &DECAY_SORBED will have no affect on simulation &
1103  &results.'
1104  end if
1105  end if
1106  !
1107  ! -- Check for required dual domain arrays or warn if they are specified
1108  ! but won't be used.
1109  if (.not. lname(5)) then
1110  write (this%iout, '(1x,a)') 'Warning. Dual domain is active but &
1111  &initial immobile domain concentration was not specified. &
1112  &Setting CIM to zero.'
1113  end if
1114  if (.not. lname(6)) then
1115  write (errmsg, '(a)') 'Dual domain is active but dual &
1116  &domain mass transfer rate (ZETAIM) was not specified. ZETAIM &
1117  &must be specified in GRIDDATA block.'
1118  call store_error(errmsg)
1119  end if
1120  if (.not. lname(7)) then
1121  write (errmsg, '(a)') 'Dual domain is active but &
1122  &immobile domain POROSITY was not specified. POROSITY &
1123  &must be specified in GRIDDATA block.'
1124  call store_error(errmsg)
1125  end if
1126  if (.not. lname(8)) then
1127  write (errmsg, '(a)') 'Dual domain is active but &
1128  &immobile domain VOLFRAC was not specified. VOLFRAC &
1129  &must be specified in GRIDDATA block. This is a new &
1130  &requirement for MODFLOW versions later than version &
1131  &6.4.1.'
1132  call store_error(errmsg)
1133  end if
1134  !
1135  ! -- terminate if errors
1136  if (count_errors() > 0) then
1137  call this%parser%StoreErrorUnit()
1138  end if
1139  !
1140  ! -- Return
1141  return
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
Here is the call graph for this function:

◆ read_options()

subroutine gwtistmodule::read_options ( class(gwtisttype), intent(inout)  this)

Read options for boundary packages.

Parameters
[in,out]thisGwtIstType object

Definition at line 818 of file gwt-ist.f90.

819  ! -- modules
821  use simmodule, only: store_error
822  use openspecmodule, only: access, form
824  ! -- dummy
825  class(GwtIstType), intent(inout) :: this !< GwtIstType object
826  ! -- local
827  character(len=LINELENGTH) :: errmsg, keyword
828  character(len=LINELENGTH) :: fname
829  character(len=:), allocatable :: keyword2
830  integer(I4B) :: ierr
831  logical :: isfound, endOfBlock
832  logical :: found
833  ! -- formats
834  character(len=*), parameter :: fmtisvflow = &
835  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
836  &WHENEVER ICBCFL IS NOT ZERO.')"
837  character(len=*), parameter :: fmtisrb = &
838  &"(4x,'LINEAR SORPTION IS SELECTED. ')"
839  character(len=*), parameter :: fmtidcy1 = &
840  &"(4x,'FIRST-ORDER DECAY IS ACTIVE. ')"
841  character(len=*), parameter :: fmtidcy2 = &
842  &"(4x,'ZERO-ORDER DECAY IS ACTIVE. ')"
843  character(len=*), parameter :: fmtistbin = &
844  "(4x, 'IST ', 1x, a, 1x, ' WILL BE SAVED TO FILE: ', a, &
845  &/4x, 'OPENED ON UNIT: ', I0)"
846  !
847  ! -- get options block
848  call this%parser%GetBlock('OPTIONS', isfound, ierr, blockrequired=.false., &
849  supportopenclose=.true.)
850  !
851  ! -- parse options block if detected
852  if (isfound) then
853  write (this%iout, '(1x,a)') 'PROCESSING IMMOBILE STORAGE AND TRANSFER &
854  &OPTIONS'
855  do
856  call this%parser%GetNextLine(endofblock)
857  if (endofblock) exit
858  call this%parser%GetStringCaps(keyword)
859  select case (keyword)
860  case ('SAVE_FLOWS')
861  this%ipakcb = -1
862  write (this%iout, fmtisvflow)
863  case ('CIM')
864  call this%parser%GetRemainingLine(keyword2)
865  call this%ocd%set_option(keyword2, this%inunit, this%iout)
866  case ('BUDGET')
867  call this%parser%GetStringCaps(keyword)
868  if (keyword == 'FILEOUT') then
869  call this%parser%GetString(fname)
870  this%ibudgetout = getunit()
871  call openfile(this%ibudgetout, this%iout, fname, 'DATA(BINARY)', &
872  form, access, 'REPLACE', mode_opt=mnormal)
873  write (this%iout, fmtistbin) 'BUDGET', trim(adjustl(fname)), &
874  this%ibudgetout
875  found = .true.
876  else
877  call store_error('Optional BUDGET keyword must &
878  &be followed by FILEOUT')
879  end if
880  case ('BUDGETCSV')
881  call this%parser%GetStringCaps(keyword)
882  if (keyword == 'FILEOUT') then
883  call this%parser%GetString(fname)
884  this%ibudcsv = getunit()
885  call openfile(this%ibudcsv, this%iout, fname, 'CSV', &
886  filstat_opt='REPLACE')
887  write (this%iout, fmtistbin) 'BUDGET CSV', trim(adjustl(fname)), &
888  this%ibudcsv
889  else
890  call store_error('Optional BUDGETCSV keyword must be followed by &
891  &FILEOUT')
892  end if
893  case ('SORBTION', 'SORPTION')
894  this%isrb = 1
895  write (this%iout, fmtisrb)
896  case ('FIRST_ORDER_DECAY')
897  this%idcy = 1
898  write (this%iout, fmtidcy1)
899  case ('ZERO_ORDER_DECAY')
900  this%idcy = 2
901  write (this%iout, fmtidcy2)
902  case default
903  write (errmsg, '(a,a)') 'Unknown IST option: ', &
904  trim(keyword)
905  call store_error(errmsg)
906  call this%parser%StoreErrorUnit()
907  end select
908  end do
909  write (this%iout, '(1x,a)') 'END OF IMMOBILE STORAGE AND TRANSFER &
910  &OPTIONS'
911  end if
912  !
913  ! -- Return
914  return
@ mnormal
normal output mode
Definition: Constants.f90:205
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
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function:

Variable Documentation

◆ budtxt

character(len=lenbudtxt), dimension(nbditems) gwtistmodule::budtxt
private

Definition at line 35 of file gwt-ist.f90.

35  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt

◆ ftype

character(len=lenftype) gwtistmodule::ftype = 'IST'
private

Definition at line 32 of file gwt-ist.f90.

32  character(len=LENFTYPE) :: ftype = 'IST'

◆ nbditems

integer(i4b), parameter gwtistmodule::nbditems = 5
private

Definition at line 34 of file gwt-ist.f90.

34  integer(I4B), parameter :: NBDITEMS = 5

◆ text

character(len=lenpackagename) gwtistmodule::text = ' IMMOBILE DOMAIN'
private

Definition at line 33 of file gwt-ist.f90.

33  character(len=LENPACKAGENAME) :: text = ' IMMOBILE DOMAIN'