MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
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, mempath, 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_calc_csrb (this, cim)
 @ brief Calculate immobile sorbed concentration 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 dependent variables. More...
 
subroutine output_immobile_concentration (this, idvsave, idvprint)
 @ brief Output immobile domain aqueous concentration. More...
 
subroutine output_immobile_sorbate_concentration (this, idvsave, idvprint)
 @ brief Output immobile domain sorbate 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 source_options (this)
 @ brief Source options for package More...
 
subroutine log_options (this, found, cim6_fname, budget_fname, budgetcsv_fname, sorbate_fname)
 Log user options to list file. More...
 
subroutine read_options (this)
 @ brief Read options for package More...
 
subroutine source_data (this)
 @ brief Source data for package More...
 
subroutine log_data (this, found)
 Log user data to list file. More...
 
subroutine ist_read_dimensions (this)
 @ brief Read dimensions for package More...
 
real(dp) function get_thetaim (this, node)
 @ brief Return thetaim More...
 
subroutine get_ddterm (thetaim, vcell, delt, swtpdt, volfracim, rhobim, kdnew, kdold, 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 1517 of file gwt-ist.f90.

1518  ! -- modules
1519  ! -- dummy
1520  real(DP), dimension(:, :), intent(inout) :: budterm !<
1521  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1522  real(DP), intent(in) :: cimnew !< immobile domain concenration at the end of this time step
1523  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1524  real(DP), intent(in) :: cnew !< mobile domain concentration at the end of this time step
1525  integer(I4B), intent(in) :: idcy !< order of decay rate (0:none, 1:first, 2:zero)
1526  ! -- local
1527  real(DP) :: rate
1528  integer(I4B) :: i
1529  !
1530  ! -- calculate STORAGE-AQUEOUS
1531  i = 1
1532  rate = -ddterm(1) * cimnew + ddterm(2) * cimold
1533  if (rate > dzero) then
1534  budterm(1, i) = budterm(1, i) + rate
1535  else
1536  budterm(2, i) = budterm(2, i) - rate
1537  end if
1538  !
1539  ! -- calculate STORAGE-SORBED
1540  i = 2
1541  rate = -ddterm(3) * cimnew + ddterm(4) * cimold
1542  if (rate > dzero) then
1543  budterm(1, i) = budterm(1, i) + rate
1544  else
1545  budterm(2, i) = budterm(2, i) - rate
1546  end if
1547  !
1548  ! -- calculate DECAY-AQUEOUS
1549  i = 3
1550  rate = dzero
1551  if (idcy == 1) then
1552  rate = -ddterm(5) * cimnew
1553  else if (idcy == 2) then
1554  rate = -ddterm(7)
1555  else
1556  rate = dzero
1557  end if
1558  if (rate > dzero) then
1559  budterm(1, i) = budterm(1, i) + rate
1560  else
1561  budterm(2, i) = budterm(2, i) - rate
1562  end if
1563  !
1564  ! -- calculate DECAY-SORBED
1565  i = 4
1566  if (idcy == 1) then
1567  rate = -ddterm(6) * cimnew
1568  else if (idcy == 2) then
1569  rate = -ddterm(8)
1570  else
1571  rate = dzero
1572  end if
1573  if (rate > dzero) then
1574  budterm(1, i) = budterm(1, i) + rate
1575  else
1576  budterm(2, i) = budterm(2, i) - rate
1577  end if
1578  !
1579  ! -- calculate MOBILE-DOMAIN
1580  i = 5
1581  rate = ddterm(9) * cnew - ddterm(9) * cimnew
1582  if (rate > dzero) then
1583  budterm(1, i) = budterm(1, i) + rate
1584  else
1585  budterm(2, i) = budterm(2, i) - rate
1586  end if
1587  !
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 854 of file gwt-ist.f90.

855  ! -- modules
858  ! -- dummy
859  class(GwtIstType) :: this !< GwtIstType object
860  ! -- local
861  !
862  ! -- call standard BndType allocate scalars
863  call this%BndType%allocate_scalars()
864  !
865  ! -- Allocate
866  call mem_allocate(this%icimout, 'ICIMOUT', this%memoryPath)
867  call mem_allocate(this%ibudgetout, 'IBUDGETOUT', this%memoryPath)
868  call mem_allocate(this%ibudcsv, 'IBUDCSV', this%memoryPath)
869  call mem_allocate(this%ioutsorbate, 'IOUTSORBATE', this%memoryPath)
870  call mem_allocate(this%isrb, 'ISRB', this%memoryPath)
871  call mem_allocate(this%idcy, 'IDCY', this%memoryPath)
872  call mem_allocate(this%kiter, 'KITER', this%memoryPath)
873  !
874  ! -- Initialize
875  this%lstfmt = ''
876  this%icimout = 0
877  this%ibudgetout = 0
878  this%ibudcsv = 0
879  this%ioutsorbate = 0
880  this%isrb = 0
881  this%idcy = 0
882  this%kiter = 0
883  !
884  ! -- Create the ocd object, which is used to manage printing and saving
885  ! of the immobile domain concentrations
886  call ocd_cr(this%ocd)
Output control data module.
subroutine, public ocd_cr(ocdobj)
@ brief Create a new output control data type.
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 1495 of file gwt-ist.f90.

1496  ! -- dummy
1497  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1498  real(DP), intent(in) :: f !< the f term used to calculate the immobile domain concentration
1499  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1500  real(DP), intent(in) :: cnew !< concentration of the mobile domain at the end of the time step
1501  ! -- result
1502  real(DP) :: cimnew !< calculated concentration of the immobile domain
1503  ! -- local
1504  !
1505  ! -- calculate ddconc
1506  cimnew = (ddterm(2) + ddterm(4)) * cimold + ddterm(9) * cnew - ddterm(7) &
1507  - ddterm(8)
1508  cimnew = cimnew / f
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)  kdnew,
real(dp), intent(in)  kdold,
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]kdneweffective distribution coefficient for new time
[in]kdoldeffective distribution coefficient for old time
[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 1422 of file gwt-ist.f90.

1426  ! -- dummy
1427  real(DP), intent(in) :: thetaim !< immobile domain porosity
1428  real(DP), intent(in) :: vcell !< volume of cell
1429  real(DP), intent(in) :: delt !< length of time step
1430  real(DP), intent(in) :: swtpdt !< cell saturation at end of time step
1431  real(DP), intent(in) :: volfracim !< volume fraction of immobile domain
1432  real(DP), intent(in) :: rhobim !< bulk density for the immobile domain (fim * rhob)
1433  real(DP), intent(in) :: kdnew !< effective distribution coefficient for new time
1434  real(DP), intent(in) :: kdold !< effective distribution coefficient for old time
1435  real(DP), intent(in) :: lambda1im !< first-order decay rate in aqueous phase
1436  real(DP), intent(in) :: lambda2im !< first-order decay rate in sorbed phase
1437  real(DP), intent(in) :: gamma1im !< zero-order decay rate in aqueous phase
1438  real(DP), intent(in) :: gamma2im !< zero-order decay rate in sorbed phase
1439  real(DP), intent(in) :: zetaim !< transfer coefficient between mobile and immobile domains
1440  real(DP), dimension(:), intent(inout) :: ddterm !< nine terms comprising the balance equation of the immobile domain
1441  real(DP), intent(inout) :: f !< the f term used to calculate the immobile domain concentration
1442  ! -- local
1443  real(DP) :: tled
1444  !
1445  ! -- initialize
1446  tled = done / delt
1447  !
1448  ! -- Calculate terms. These terms correspond to the concentration
1449  ! coefficients in equation 7-4 of the GWT model report. However,
1450  ! an updated equation is presented as 9-9 in the supplemental technical
1451  ! information guide (mf6suptechinfo.pdf)
1452  ddterm(1) = thetaim * vcell * tled
1453  ddterm(2) = thetaim * vcell * tled
1454  ddterm(3) = volfracim * rhobim * vcell * kdnew * tled
1455  ddterm(4) = volfracim * rhobim * vcell * kdold * tled
1456  ddterm(5) = thetaim * lambda1im * vcell
1457  ddterm(6) = lambda2im * volfracim * rhobim * kdnew * vcell
1458  ddterm(7) = thetaim * gamma1im * vcell
1459  ddterm(8) = gamma2im * volfracim * rhobim * vcell
1460  ddterm(9) = vcell * swtpdt * zetaim
1461  !
1462  ! -- calculate denominator term, f
1463  f = ddterm(1) + ddterm(3) + ddterm(5) + ddterm(6) + ddterm(9)
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 1472 of file gwt-ist.f90.

1473  ! -- dummy
1474  real(DP), dimension(:), intent(in) :: ddterm !< terms comprising the balance equation of the immobile domain
1475  real(DP), intent(in) :: f !< the f term used to calculate the immobile domain concentration
1476  real(DP), intent(in) :: cimold !< immobile domain concentration at end of last time step
1477  real(DP), intent(inout) :: hcof !< calculated contribution for the a-matrix diagonal position
1478  real(DP), intent(inout) :: rhs !< calculated contribution for the solution right-hand side
1479  !
1480  ! -- calculate hcof
1481  hcof = ddterm(9)**2 / f - ddterm(9)
1482  !
1483  ! -- calculate rhs, and switch the sign because this term needs to
1484  ! be moved to the left hand side
1485  rhs = (ddterm(2) + ddterm(4)) * cimold - ddterm(7) - ddterm(8)
1486  rhs = rhs * ddterm(9) / f
1487  rhs = -rhs
Here is the caller graph for this function:

◆ get_thetaim()

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

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

Parameters
thisGwtIstType object
[in]nodenode number

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

1404  ! -- modules
1405  ! -- dummy
1406  class(GwtIstType) :: this !< GwtIstType object
1407  integer(I4B), intent(in) :: node !< node number
1408  ! -- return
1409  real(DP) :: thetaim
1410  !
1411  thetaim = this%volfrac(node) * this%porosity(node)

◆ 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 242 of file gwt-ist.f90.

243  ! -- modules
245  ! -- dummy variables
246  class(GwtIstType) :: this !< GwtIstType object
247  ! -- local variables
248  integer(I4B) :: n
249  !
250  ! -- Call parent advance
251  call this%BndType%bnd_ad()
252  !
253  ! -- set independent kiter counter to zero
254  this%kiter = 0
255  !
256  ! -- copy cimnew into cimold or vice versa if this is a repeat of
257  ! a failed time step
258  if (ifailedstepretry == 0) then
259  do n = 1, this%dis%nodes
260  this%cimold(n) = this%cimnew(n)
261  end do
262  else
263  do n = 1, this%dis%nodes
264  this%cimnew(n) = this%cimold(n)
265  end do
266  end if
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 894 of file gwt-ist.f90.

895  ! -- modules
897  ! -- dummy
898  class(GwtIstType), intent(inout) :: this !< GwtIstType object
899  ! -- local
900  integer(I4B) :: n
901  !
902  ! -- call standard BndType allocate scalars
903  ! nbound and maxbound are 0 in order to keep memory footprint low
904  call this%BndType%allocate_arrays()
905  !
906  ! -- allocate ist arrays of size nodes
907  call mem_allocate(this%strg, this%dis%nodes, 'STRG', this%memoryPath)
908  call mem_allocate(this%cim, this%dis%nodes, 'CIM', this%memoryPath)
909  call mem_allocate(this%cimnew, this%dis%nodes, 'CIMNEW', this%memoryPath)
910  call mem_allocate(this%cimold, this%dis%nodes, 'CIMOLD', this%memoryPath)
911  call mem_allocate(this%porosity, this%dis%nodes, 'POROSITY', this%memoryPath)
912  call mem_allocate(this%zetaim, this%dis%nodes, 'ZETAIM', this%memoryPath)
913  call mem_allocate(this%volfrac, this%dis%nodes, 'VOLFRAC', this%memoryPath)
914  if (this%isrb == 0) then
915  call mem_allocate(this%bulk_density, 1, 'BULK_DENSITY', this%memoryPath)
916  call mem_allocate(this%distcoef, 1, 'DISTCOEF', this%memoryPath)
917  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
918  call mem_allocate(this%cimsrb, 1, 'SORBATE', this%memoryPath)
919  else
920  call mem_allocate(this%bulk_density, this%dis%nodes, 'BULK_DENSITY', &
921  this%memoryPath)
922  call mem_allocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
923  this%memoryPath)
924  call mem_allocate(this%cimsrb, this%dis%nodes, 'SORBATE', &
925  this%memoryPath)
926  if (this%isrb == 1) then
927  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
928  else
929  call mem_allocate(this%sp2, this%dis%nodes, 'SP2', this%memoryPath)
930  end if
931  end if
932  if (this%idcy == 0) then
933  call mem_allocate(this%decay, 1, 'DECAY', this%memoryPath)
934  call mem_allocate(this%decaylast, 1, 'DECAYLAST', this%memoryPath)
935  else
936  call mem_allocate(this%decay, this%dis%nodes, 'DECAY', this%memoryPath)
937  call mem_allocate(this%decaylast, this%dis%nodes, 'DECAYLAST', &
938  this%memoryPath)
939  end if
940  if (this%isrb == 0 .and. this%idcy == 0) then
941  call mem_allocate(this%decayslast, 1, 'DECAYSLAST', this%memoryPath)
942  else
943  call mem_allocate(this%decayslast, this%dis%nodes, 'DECAYSLAST', &
944  this%memoryPath)
945  end if
946  call mem_allocate(this%decay_sorbed, 1, 'DECAY_SORBED', this%memoryPath)
947  !
948  ! -- initialize
949  do n = 1, this%dis%nodes
950  this%strg(n) = dzero
951  this%cim(n) = dzero
952  this%cimnew(n) = dzero
953  this%cimold(n) = dzero
954  this%porosity(n) = dzero
955  this%zetaim(n) = dzero
956  this%volfrac(n) = dzero
957  end do
958  do n = 1, size(this%bulk_density)
959  this%bulk_density(n) = dzero
960  this%distcoef(n) = dzero
961  this%cimsrb(n) = dzero
962  end do
963  do n = 1, size(this%sp2)
964  this%sp2(n) = dzero
965  end do
966  do n = 1, size(this%decay)
967  this%decay(n) = dzero
968  this%decaylast(n) = dzero
969  end do
970  do n = 1, size(this%decayslast)
971  this%decayslast(n) = dzero
972  end do
973  !
974  ! -- Set pointers
975  this%ocd%dis => this%dis

◆ 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 166 of file gwt-ist.f90.

167  ! -- modules
169  use budgetmodule, only: budget_cr
170  ! -- dummy
171  class(GwtIstType), intent(inout) :: this !< GwtIstType object
172  ! -- local
173  integer(I4B) :: n
174  !
175  ! -- Allocate arrays
176  call this%ist_allocate_arrays()
177  !
178  ! -- Now that arrays are allocated, check in the cimnew array to
179  ! the output control manager for subsequent printing/saving
180  call this%ocd%init_dbl('CIM', this%cimnew, this%dis, 'PRINT LAST ', &
181  'COLUMNS 10 WIDTH 11 DIGITS 4 GENERAL ', &
182  this%iout, dhnoflo)
183 
184  ! -- apply user override if provided
185  if (this%lstfmt /= '') then
186  call this%ocd%set_option(trim(this%lstfmt)//" ", 0, this%iout)
187  end if
188  !
189  ! -- source the data block
190  call this%source_data()
191  !
192  ! -- set cimnew to the cim start values read from input
193  do n = 1, this%dis%nodes
194  this%cimnew(n) = this%cim(n)
195  end do
196  !
197  ! -- add volfrac to the volfracim accumulator in mst package
198  call this%mst%addto_volfracim(this%volfrac)
199  !
200  ! -- setup the immobile domain budget
201  call budget_cr(this%budget, this%memoryPath)
202  call this%budget%budget_df(nbditems, 'MASS', 'M', bdzone=this%packName)
203  call this%budget%set_ibudcsv(this%ibudcsv)
204  !
205  ! -- Perform a check to ensure that sorption and decay are set
206  ! consistently between the MST and IST packages.
207  if (this%idcy /= this%mst%idcy) then
208  call store_error('DECAY must be activated consistently between the &
209  &MST and IST Packages. Activate or deactivate DECAY for both &
210  &Packages.')
211  end if
212  if (this%isrb /= this%mst%isrb) then
213  call store_error('SORPTION must be activated consistently between the &
214  &MST and IST Packages. Activate or deactivate SORPTION for both &
215  &Packages. If activated, the same type of sorption (LINEAR, &
216  &FREUNDLICH, or LANGMUIR) must be specified in the options block of &
217  &both the MST and IST Packages.')
218  end if
219  if (count_errors() > 0) then
220  call store_error_filename(this%input_fname)
221  end if
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
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
Here is the call graph for this function:

◆ ist_bd()

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

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 599 of file gwt-ist.f90.

600  ! -- modules
601  use tdismodule, only: delt
603  ! -- dummy
604  class(GwtIstType) :: this !< GwtIstType object
605  type(BudgetType), intent(inout) :: model_budget !< model budget object
606  ! -- local
607  real(DP) :: ratin
608  real(DP) :: ratout
609  integer(I4B) :: isuppress_output
610  isuppress_output = 0
611  call rate_accumulator(this%strg(:), ratin, ratout)
612  call model_budget%addentry(ratin, ratout, delt, this%text, &
613  isuppress_output, this%packName)
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
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_calc_csrb()

subroutine gwtistmodule::ist_calc_csrb ( class(gwtisttype this,
real(dp), dimension(:), intent(in)  cim 
)
Parameters
thisGwtMstType object
[in]cimimmobile domain aqueous concentration at end of this time step

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

567  ! -- dummy
568  class(GwtIstType) :: this !< GwtMstType object
569  real(DP), intent(in), dimension(:) :: cim !< immobile domain aqueous concentration at end of this time step
570  ! -- local
571  integer(I4B) :: n
572  real(DP) :: distcoef
573  real(DP) :: csrb
574 
575  ! Calculate sorbed concentration
576  do n = 1, size(cim)
577  csrb = dzero
578  if (this%ibound(n) > 0) then
579  distcoef = this%distcoef(n)
580  if (this%isrb == 1) then
581  csrb = cim(n) * distcoef
582  else if (this%isrb == 2) then
583  csrb = get_freundlich_conc(cim(n), distcoef, this%sp2(n))
584  else if (this%isrb == 3) then
585  csrb = get_langmuir_conc(cim(n), distcoef, this%sp2(n))
586  end if
587  end if
588  this%cimsrb(n) = csrb
589  end do
590 
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 
)
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 409 of file gwt-ist.f90.

410  ! modules
411  use tdismodule, only: delt
412  use constantsmodule, only: dzero
413  ! dummy
414  class(GwtIstType), intent(inout) :: this !< GwtIstType object
415  real(DP), dimension(:), intent(in) :: x !< current dependent-variable value
416  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
417  integer(I4B), optional, intent(in) :: iadv !< flag that indicates if this is an advance package
418  ! local
419  integer(I4B) :: idiag
420  integer(I4B) :: n
421  real(DP) :: rate
422  real(DP) :: swt, swtpdt
423  real(DP) :: hhcof, rrhs
424  real(DP) :: vcell
425  real(DP) :: thetaim
426  real(DP) :: zetaim
427  real(DP) :: kdnew
428  real(DP) :: kdold
429  real(DP) :: volfracim
430  real(DP) :: rhobim
431  real(DP) :: lambda1im
432  real(DP) :: lambda2im
433  real(DP) :: gamma1im
434  real(DP) :: gamma2im
435  real(DP) :: cimnew
436  real(DP) :: cimold
437  real(DP) :: f
438  real(DP) :: cimsrbold
439  real(DP) :: cimsrbnew
440  real(DP), dimension(9) :: ddterm
441  ! formats
442 
443  ! initialize
444  this%budterm(:, :) = dzero
445 
446  ! Calculate immobile domain transfer rate
447  do n = 1, this%dis%nodes
448 
449  ! skip if transport inactive
450  rate = dzero
451  cimnew = dzero
452  if (this%ibound(n) > 0) then
453 
454  ! calculate new and old water volumes
455  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
456  swtpdt = this%fmi%gwfsat(n)
457  swt = this%fmi%gwfsatold(n, delt)
458  thetaim = this%get_thetaim(n)
459 
460  ! set exchange coefficient
461  zetaim = this%zetaim(n)
462 
463  ! Calculate exchange with immobile domain
464  rate = dzero
465  hhcof = dzero
466  rrhs = dzero
467  kdnew = dzero
468  kdold = dzero
469  volfracim = dzero
470  rhobim = dzero
471  lambda1im = dzero
472  lambda2im = dzero
473  gamma1im = dzero
474  gamma2im = dzero
475  if (this%idcy == 1) lambda1im = this%decay(n)
476  if (this%idcy == 2) then
477  gamma1im = get_zero_order_decay(this%decay(n), this%decaylast(n), 0, &
478  this%cimold(n), this%cimnew(n), delt)
479  end if
480 
481  ! setup sorption variables
482  if (this%isrb > 0) then
483 
484  ! initialize sorption variables
485  volfracim = this%volfrac(n)
486  rhobim = this%bulk_density(n)
487 
488  ! set isotherm dependent sorption variables
489  select case (this%isrb)
490  case (1)
491  ! linear
492  kdnew = this%distcoef(n)
493  kdold = this%distcoef(n)
494  cimsrbnew = this%cimnew(n) * kdnew
495  cimsrbold = this%cimold(n) * kdold
496  case (2)
497  ! freundlich
498  kdnew = get_freundlich_kd(this%cimnew(n), this%distcoef(n), &
499  this%sp2(n))
500  kdold = get_freundlich_kd(this%cimold(n), this%distcoef(n), &
501  this%sp2(n))
502  cimsrbnew = get_freundlich_conc(this%cimnew(n), this%distcoef(n), &
503  this%sp2(n))
504  cimsrbold = get_freundlich_conc(this%cimold(n), this%distcoef(n), &
505  this%sp2(n))
506  case (3)
507  ! langmuir
508  kdnew = get_langmuir_kd(this%cimnew(n), this%distcoef(n), &
509  this%sp2(n))
510  kdold = get_langmuir_kd(this%cimold(n), this%distcoef(n), &
511  this%sp2(n))
512  cimsrbnew = get_langmuir_conc(this%cimnew(n), this%distcoef(n), &
513  this%sp2(n))
514  cimsrbold = get_langmuir_conc(this%cimold(n), this%distcoef(n), &
515  this%sp2(n))
516  end select
517 
518  ! set decay of sorbed solute parameters
519  if (this%idcy == 1) then
520  lambda2im = this%decay_sorbed(n)
521  else if (this%idcy == 2) then
522  gamma2im = get_zero_order_decay(this%decay_sorbed(n), &
523  this%decayslast(n), &
524  0, cimsrbold, &
525  cimsrbnew, delt)
526  end if
527  end if
528 
529  ! calculate the terms and then get the hcof and rhs contributions
530  call get_ddterm(thetaim, vcell, delt, swtpdt, &
531  volfracim, rhobim, kdnew, kdold, lambda1im, lambda2im, &
532  gamma1im, gamma2im, zetaim, ddterm, f)
533  cimold = this%cimold(n)
534  call get_hcofrhs(ddterm, f, cimold, hhcof, rrhs)
535 
536  ! calculate rate from hcof and rhs
537  rate = hhcof * x(n) - rrhs
538 
539  ! calculate immobile domain concentration
540  cimnew = get_ddconc(ddterm, f, cimold, x(n))
541 
542  ! accumulate the budget terms
543  call accumulate_budterm(this%budterm, ddterm, cimnew, cimold, x(n), &
544  this%idcy)
545  end if
546 
547  ! store rate and add to flowja
548  this%strg(n) = rate
549  idiag = this%dis%con%ia(n)
550  flowja(idiag) = flowja(idiag) + rate
551 
552  ! store immobile domain concentration
553  this%cimnew(n) = cimnew
554 
555  end do
556 
557  ! calculate csrb
558  if (this%isrb /= 0) then
559  call this%ist_calc_csrb(this%cimnew)
560  end if
561 
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
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,
character(len=*), intent(in)  mempath,
type(tspfmitype), intent(in), pointer  fmi,
type(gwtmsttype), intent(in), 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 118 of file gwt-ist.f90.

120  ! -- dummy
121  class(BndType), pointer :: packobj !< BndType pointer that will point to new IST Package
122  integer(I4B), intent(in) :: id !< name of the model
123  integer(I4B), intent(in) :: ibcnum !< consecutive package number
124  integer(I4B), intent(in) :: inunit !< unit number of package input file
125  integer(I4B), intent(in) :: iout !< unit number of model listing file
126  character(len=*), intent(in) :: namemodel !< name of the model
127  character(len=*), intent(in) :: pakname !< name of the package
128  character(len=*), intent(in) :: mempath
129  type(TspFmiType), pointer, intent(in) :: fmi
130  type(GwtMstType), pointer, intent(in) :: mst
131  ! -- local
132  type(GwtIstType), pointer :: istobj
133  !
134  ! -- allocate the object and assign values to object variables
135  allocate (istobj)
136  packobj => istobj
137  !
138  ! -- create name and memory path
139  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
140  packobj%text = text
141  !
142  ! -- allocate scalars
143  call packobj%allocate_scalars()
144  !
145  ! -- initialize package
146  call packobj%pack_initialize()
147  !
148  ! -- store values
149  packobj%inunit = inunit
150  packobj%iout = iout
151  packobj%id = id
152  packobj%ibcnum = ibcnum
153  packobj%ncolbnd = 1
154  packobj%iscloc = 1
155  !
156  ! -- Point IST specific variables
157  istobj%fmi => fmi
158  istobj%mst => mst
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 803 of file gwt-ist.f90.

804  ! -- modules
806  ! -- dummy
807  class(GwtIstType) :: this !< GwtIstType object
808  !
809  ! -- Deallocate arrays if package was active
810  if (this%inunit > 0) then
811  call mem_deallocate(this%icimout)
812  call mem_deallocate(this%ibudgetout)
813  call mem_deallocate(this%ibudcsv)
814  call mem_deallocate(this%ioutsorbate)
815  call mem_deallocate(this%idcy)
816  call mem_deallocate(this%isrb)
817  call mem_deallocate(this%kiter)
818  call mem_deallocate(this%cim)
819  call mem_deallocate(this%cimnew)
820  call mem_deallocate(this%cimold)
821  call mem_deallocate(this%cimsrb)
822  call mem_deallocate(this%zetaim)
823  call mem_deallocate(this%porosity)
824  call mem_deallocate(this%volfrac)
825  call mem_deallocate(this%bulk_density)
826  call mem_deallocate(this%distcoef)
827  call mem_deallocate(this%sp2)
828  call mem_deallocate(this%decay)
829  call mem_deallocate(this%decaylast)
830  call mem_deallocate(this%decayslast)
831  call mem_deallocate(this%decay_sorbed)
832  call mem_deallocate(this%strg)
833  this%fmi => null()
834  this%mst => null()
835  end if
836  !
837  ! -- Scalars
838  !
839  ! -- Objects
840  call this%budget%budget_da()
841  deallocate (this%budget)
842  call this%ocd%ocd_da()
843  deallocate (this%ocd)
844  !
845  ! -- deallocate parent
846  call this%BndType%bnd_da()

◆ 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 
)
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 271 of file gwt-ist.f90.

272  ! modules
273  use tdismodule, only: delt
274  ! dummy
275  class(GwtIstType) :: this !< GwtIstType object
276  real(DP), dimension(:), intent(inout) :: rhs !< right-hand side vector for model
277  integer(I4B), dimension(:), intent(in) :: ia !< solution CRS row pointers
278  integer(I4B), dimension(:), intent(in) :: idxglo !< mapping vector for model (local) to solution (global)
279  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
280  ! local
281  integer(I4B) :: n, idiag
282  real(DP) :: tled
283  real(DP) :: hhcof, rrhs
284  real(DP) :: swt, swtpdt
285  real(DP) :: vcell
286  real(DP) :: thetaim
287  real(DP) :: zetaim
288  real(DP) :: kdnew
289  real(DP) :: kdold
290  real(DP) :: volfracim
291  real(DP) :: rhobim
292  real(DP) :: lambda1im
293  real(DP) :: lambda2im
294  real(DP) :: gamma1im
295  real(DP) :: gamma2im
296  real(DP) :: cimold
297  real(DP) :: f
298  real(DP) :: cimsrbold
299  real(DP) :: cimsrbnew
300  real(DP), dimension(9) :: ddterm
301 
302  ! set variables
303  tled = done / delt
304  this%kiter = this%kiter + 1
305 
306  ! loop through each node and calculate immobile domain contribution
307  ! to hcof and rhs
308  do n = 1, this%dis%nodes
309 
310  ! skip if transport inactive
311  if (this%ibound(n) <= 0) cycle
312 
313  ! calculate new and old water volumes
314  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
315  swtpdt = this%fmi%gwfsat(n)
316  swt = this%fmi%gwfsatold(n, delt)
317  thetaim = this%get_thetaim(n)
318  idiag = ia(n)
319 
320  ! set exchange coefficient
321  zetaim = this%zetaim(n)
322 
323  ! Add dual domain mass transfer contributions to rhs and hcof
324  ! dcimsrbdc = DZERO
325  kdnew = dzero
326  kdold = dzero
327  volfracim = dzero
328  rhobim = dzero
329  lambda1im = dzero
330  lambda2im = dzero
331  gamma1im = dzero
332  gamma2im = dzero
333 
334  ! set variables for decay of aqueous solute
335  if (this%idcy == 1) lambda1im = this%decay(n)
336  if (this%idcy == 2) then
337  gamma1im = get_zero_order_decay(this%decay(n), this%decaylast(n), &
338  this%kiter, this%cimold(n), &
339  this%cimnew(n), delt)
340  this%decaylast(n) = gamma1im
341  end if
342 
343  ! setup sorption variables
344  if (this%isrb > 0) then
345 
346  ! initialize sorption variables
347  volfracim = this%volfrac(n)
348  rhobim = this%bulk_density(n)
349 
350  ! set isotherm dependent sorption variables
351  select case (this%isrb)
352  case (1)
353  ! linear
354  kdnew = this%distcoef(n)
355  kdold = this%distcoef(n)
356  cimsrbnew = this%cimnew(n) * kdnew
357  cimsrbold = this%cimold(n) * kdold
358  case (2)
359  ! freundlich
360  kdnew = get_freundlich_kd(this%cimnew(n), this%distcoef(n), &
361  this%sp2(n))
362  kdold = get_freundlich_kd(this%cimold(n), this%distcoef(n), &
363  this%sp2(n))
364  cimsrbnew = get_freundlich_conc(this%cimnew(n), this%distcoef(n), &
365  this%sp2(n))
366  cimsrbold = get_freundlich_conc(this%cimold(n), this%distcoef(n), &
367  this%sp2(n))
368  case (3)
369  ! langmuir
370  kdnew = get_langmuir_kd(this%cimnew(n), this%distcoef(n), &
371  this%sp2(n))
372  kdold = get_langmuir_kd(this%cimold(n), this%distcoef(n), &
373  this%sp2(n))
374  cimsrbnew = get_langmuir_conc(this%cimnew(n), this%distcoef(n), &
375  this%sp2(n))
376  cimsrbold = get_langmuir_conc(this%cimold(n), this%distcoef(n), &
377  this%sp2(n))
378  end select
379 
380  ! set decay of sorbed solute parameters
381  if (this%idcy == 1) then
382  lambda2im = this%decay_sorbed(n)
383  else if (this%idcy == 2) then
384  gamma2im = get_zero_order_decay(this%decay_sorbed(n), &
385  this%decayslast(n), &
386  this%kiter, cimsrbold, &
387  cimsrbnew, delt)
388  this%decayslast(n) = gamma2im
389  end if
390  end if
391 
392  ! calculate dual domain terms and get the hcof and rhs contributions
393  call get_ddterm(thetaim, vcell, delt, swtpdt, &
394  volfracim, rhobim, kdnew, kdold, lambda1im, lambda2im, &
395  gamma1im, gamma2im, zetaim, ddterm, f)
396  cimold = this%cimold(n)
397  call get_hcofrhs(ddterm, f, cimold, hhcof, rrhs)
398 
399  ! update solution accumulators
400  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
401  rhs(n) = rhs(n) + rrhs
402 
403  end do
404 
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 
)
private

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 772 of file gwt-ist.f90.

773  ! -- modules
774  use tdismodule, only: delt, totim
775  ! -- dummy variables
776  class(GwtIstType) :: this !< GwtIstType object
777  integer(I4B), intent(in) :: kstp !< time step number
778  integer(I4B), intent(in) :: kper !< period number
779  integer(I4B), intent(in) :: iout !< flag and unit number for the model listing file
780  integer(I4B), intent(in) :: ibudfl !< flag indicating budget should be written
781  ! -- local
782  integer(I4B) :: isuppress_output = 0
783  !
784  ! -- Fill budget terms
785  call this%budget%reset()
786  call this%budget%addentry(this%budterm, delt, budtxt, isuppress_output)
787  !
788  ! -- Write budget to list file
789  call this%budget%finalize_step(delt)
790  if (ibudfl /= 0) then
791  call this%budget%budget_ot(kstp, kper, iout)
792  end if
793  !
794  ! -- Write budget csv
795  call this%budget%writecsv(totim)
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 685 of file gwt-ist.f90.

686  ! dummy variables
687  class(GwtIstType) :: this !< BndType object
688  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
689  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
690 
691  call this%output_immobile_concentration(idvsave, idvprint)
692  call this%output_immobile_sorbate_concentration(idvsave, idvprint)
693 

◆ 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 622 of file gwt-ist.f90.

623  ! -- modules
624  use constantsmodule, only: dzero
625  ! -- dummy
626  class(GwtIstType) :: this !< GwtIstType object
627  integer(I4B), intent(in) :: icbcfl !< flag for cell-by-cell output
628  integer(I4B), intent(in) :: ibudfl !< flag indication if cell-by-cell data should be saved
629  integer(I4B), intent(in) :: icbcun !< unit number for cell-by-cell output
630  integer(I4B), dimension(:), optional, intent(in) :: imap !< mapping vector
631  ! -- local
632  integer(I4B) :: n
633  integer(I4B) :: ibinun
634  integer(I4B) :: nbound
635  integer(I4B) :: naux
636  real(DP) :: rate
637  real(DP), dimension(0) :: auxrow
638  !
639  ! -- Set unit number for binary output
640  if (this%ipakcb < 0) then
641  ibinun = icbcun
642  elseif (this%ipakcb == 0) then
643  ibinun = 0
644  else
645  ibinun = this%ipakcb
646  end if
647  if (icbcfl == 0) ibinun = 0
648  !
649  ! -- Record the storage rate if requested
650  !
651  ! -- If cell-by-cell flows will be saved as a list, write header.
652  if (ibinun /= 0) then
653  nbound = this%dis%nodes
654  naux = 0
655  call this%dis%record_srcdst_list_header(this%text, this%name_model, &
656  this%name_model, this%name_model, &
657  this%packName, naux, this%auxname, &
658  ibinun, nbound, this%iout)
659  end if
660  !
661  ! -- Calculate immobile domain rhs and hcof
662  do n = 1, this%dis%nodes
663  !
664  ! -- skip if transport inactive
665  rate = dzero
666  if (this%ibound(n) > 0) then
667  !
668  ! -- set rate from this%strg
669  rate = this%strg(n)
670  end if
671  !
672  ! -- If saving cell-by-cell flows in list, write flow
673  if (ibinun /= 0) then
674  call this%dis%record_mf6_list_entry(ibinun, n, n, rate, &
675  naux, auxrow, &
676  olconv=.true., &
677  olconv2=.true.)
678  end if
679  !
680  end do

◆ 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 1391 of file gwt-ist.f90.

1392  ! -- dummy
1393  class(GwtIstType), intent(inout) :: this !< GwtIstType object
1394  ! -- local
1395  ! -- format

◆ 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 229 of file gwt-ist.f90.

230  ! -- dummy
231  class(GwtIstType), intent(inout) :: this !< GwtIstType object
232  ! -- local
233  ! -- format

◆ log_data()

subroutine gwtistmodule::log_data ( class(gwtisttype), intent(inout)  this,
type(gwtistparamfoundtype), intent(in)  found 
)

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

1352  class(GwTIstType), intent(inout) :: this
1353  type(GwtIstParamFoundType), intent(in) :: found
1354 
1355  write (this%iout, '(1x,a)') 'PROCESSING GRIDDATA'
1356  if (found%porosity) then
1357  write (this%iout, '(4x,a)') 'MOBILE DOMAIN POROSITY set from input file'
1358  end if
1359  if (found%volfrac) then
1360  write (this%iout, '(4x,a)') 'VOLFRAC set from input file'
1361  end if
1362  if (found%zetaim) then
1363  write (this%iout, '(4x,a)') 'ZETAIM set from input file'
1364  end if
1365  if (found%cim) then
1366  write (this%iout, '(4x,a)') 'CIM set from input file'
1367  end if
1368  if (found%decay) then
1369  write (this%iout, '(4x,a)') 'DECAY RATE set from input file'
1370  end if
1371  if (found%decay_sorbed) then
1372  write (this%iout, '(4x,a)') 'DECAY SORBED RATE set from input file'
1373  end if
1374  if (found%bulk_density) then
1375  write (this%iout, '(4x,a)') 'BULK DENSITY set from input file'
1376  end if
1377  if (found%distcoef) then
1378  write (this%iout, '(4x,a)') 'DISTRIBUTION COEFFICIENT set from input file'
1379  end if
1380  if (found%sp2) then
1381  write (this%iout, '(4x,a)') 'SECOND SORPTION PARAM set from input file'
1382  end if
1383  write (this%iout, '(1x,a)') 'END PROCESSING GRIDDATA'

◆ log_options()

subroutine gwtistmodule::log_options ( class(gwtisttype), intent(inout)  this,
type(gwtistparamfoundtype), intent(in)  found,
character(len=*), intent(in)  cim6_fname,
character(len=*), intent(in)  budget_fname,
character(len=*), intent(in)  budgetcsv_fname,
character(len=*), intent(in)  sorbate_fname 
)

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

1083  class(GwTIstType), intent(inout) :: this
1084  type(GwtIstParamFoundType), intent(in) :: found
1085  character(len=*), intent(in) :: cim6_fname
1086  character(len=*), intent(in) :: budget_fname
1087  character(len=*), intent(in) :: budgetcsv_fname
1088  character(len=*), intent(in) :: sorbate_fname
1089  ! -- formats
1090  character(len=*), parameter :: fmtisvflow = &
1091  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
1092  &WHENEVER ICBCFL IS NOT ZERO.')"
1093  character(len=*), parameter :: fmtlinear = &
1094  &"(4x,'LINEAR SORPTION IS SELECTED. ')"
1095  character(len=*), parameter :: fmtfreundlich = &
1096  &"(4x,'FREUNDLICH SORPTION IS ACTIVE. ')"
1097  character(len=*), parameter :: fmtlangmuir = &
1098  &"(4x,'LANGMUIR SORPTION IS ACTIVE. ')"
1099  character(len=*), parameter :: fmtidcy1 = &
1100  &"(4x,'FIRST-ORDER DECAY IS ACTIVE. ')"
1101  character(len=*), parameter :: fmtidcy2 = &
1102  &"(4x,'ZERO-ORDER DECAY IS ACTIVE. ')"
1103  character(len=*), parameter :: fmtistbin = &
1104  "(4x, 'IST ', 1x, a, 1x, ' WILL BE SAVED TO FILE: ', a, &
1105  &/4x, 'OPENED ON UNIT: ', I0)"
1106 
1107  write (this%iout, '(1x,a)') 'PROCESSING IMMOBILE STORAGE AND TRANSFER &
1108  &OPTIONS'
1109  if (found%save_flows) then
1110  write (this%iout, fmtisvflow)
1111  end if
1112  if (found%budgetfile) then
1113  write (this%iout, fmtistbin) 'BUDGET', trim(adjustl(budget_fname)), &
1114  this%ibudgetout
1115  end if
1116  if (found%budgetcsvfile) then
1117  write (this%iout, fmtistbin) 'BUDGET CSV', trim(adjustl(budgetcsv_fname)), &
1118  this%ibudcsv
1119  end if
1120  if (found%sorption) then
1121  select case (this%isrb)
1122  case (1)
1123  write (this%iout, fmtlinear)
1124  case (2)
1125  write (this%iout, fmtfreundlich)
1126  case (3)
1127  write (this%iout, fmtlangmuir)
1128  end select
1129  end if
1130  if (found%order1_decay) then
1131  write (this%iout, fmtidcy1)
1132  end if
1133  if (found%order0_decay) then
1134  write (this%iout, fmtidcy2)
1135  end if
1136  if (found%sorbatefile) then
1137  write (this%iout, fmtistbin) &
1138  'SORBATE', sorbate_fname, this%ioutsorbate
1139  end if
1140  write (this%iout, '(1x,a)') 'END OF IMMOBILE STORAGE AND TRANSFER &
1141  &OPTIONS'

◆ output_immobile_concentration()

subroutine gwtistmodule::output_immobile_concentration ( class(gwtisttype this,
integer(i4b), intent(in)  idvsave,
integer(i4b), intent(in)  idvprint 
)
private
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 698 of file gwt-ist.f90.

699  ! modules
700  use tdismodule, only: kstp, endofperiod
701  ! dummy variables
702  class(GwtIstType) :: this !< BndType object
703  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
704  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
705  ! local
706  integer(I4B) :: ipflg
707  integer(I4B) :: ibinun
708  !
709  ! Save cim to a binary file. ibinun is a flag where 1 indicates that
710  ! cim should be written to a binary file if a binary file is open for it.
711  ipflg = 0
712  ibinun = 1
713  if (idvsave == 0) ibinun = 0
714  if (ibinun /= 0) then
715  call this%ocd%ocd_ot(ipflg, kstp, endofperiod, this%iout, &
716  iprint_opt=0, isav_opt=ibinun)
717  end if
718  !
719  ! Print immobile domain concentrations to listing file
720  if (idvprint /= 0) then
721  call this%ocd%ocd_ot(ipflg, kstp, endofperiod, this%iout, &
722  iprint_opt=idvprint, isav_opt=0)
723  end if
724 
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

◆ output_immobile_sorbate_concentration()

subroutine gwtistmodule::output_immobile_sorbate_concentration ( 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 729 of file gwt-ist.f90.

730  ! modules
731  ! dummy
732  class(GwtIstType) :: this !< BndType object
733  integer(I4B), intent(in) :: idvsave !< flag and unit number for dependent-variable output
734  integer(I4B), intent(in) :: idvprint !< flag indicating if dependent-variable should be written to the model listing file
735  ! local
736  character(len=1) :: cdatafmp = ' ', editdesc = ' '
737  integer(I4B) :: ibinun
738  integer(I4B) :: iprint, nvaluesp, nwidthp
739  real(DP) :: dinact
740 
741  ! Save cimsrb to a binary file. ibinun is a flag where 1 indicates that
742  ! cim should be written to a binary file if a binary file is open for it.
743  ! Set unit number for sorbate output
744  if (this%ioutsorbate /= 0) then
745  ibinun = 1
746  else
747  ibinun = 0
748  end if
749  if (idvsave == 0) ibinun = 0
750 
751  ! save sorbate concentration array
752  if (ibinun /= 0) then
753  iprint = 0
754  dinact = dhnoflo
755  if (this%ioutsorbate /= 0) then
756  ibinun = this%ioutsorbate
757  call this%dis%record_array(this%cimsrb, this%iout, iprint, ibinun, &
758  ' SORBATE', cdatafmp, nvaluesp, &
759  nwidthp, editdesc, dinact)
760  end if
761  end if
762 

◆ read_options()

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

Read options for boundary packages.

Parameters
[in,out]thisGwtIstType object

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

1150  ! -- modules
1151  ! -- dummy
1152  class(GwtIstType), intent(inout) :: this !< GwtIstType object
1153 
1154  ! -- source options
1155  call this%source_options()

◆ source_data()

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

Method to source data for the package.

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

1163  ! -- modules
1164  use constantsmodule, only: linelength
1165  use simvariablesmodule, only: errmsg, warnmsg
1171  ! -- dummy
1172  class(GwtIsttype) :: this
1173  ! -- locals
1174  !character(len=LINELENGTH) :: errmsg
1175  type(GwtIstParamFoundType) :: found
1176  integer(I4B) :: asize
1177  integer(I4B), dimension(:), pointer, contiguous :: map
1178  !
1179  ! -- set map to convert user input data into reduced data
1180  map => null()
1181  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1182  !
1183  ! -- reallocate
1184  if (this%isrb == 0) then
1185  call get_isize('BULK_DENSITY', this%input_mempath, asize)
1186  if (asize > 0) &
1187  call mem_reallocate(this%bulk_density, this%dis%nodes, &
1188  'BULK_DENSITY', this%memoryPath)
1189  call get_isize('DISTCOEF', this%input_mempath, asize)
1190  if (asize > 0) &
1191  call mem_reallocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
1192  this%memoryPath)
1193  end if
1194  if (this%idcy == 0) then
1195  call get_isize('DECAY', this%input_mempath, asize)
1196  if (asize > 0) &
1197  call mem_reallocate(this%decay, this%dis%nodes, 'DECAY', this%memoryPath)
1198  end if
1199  call get_isize('DECAY_SORBED', this%input_mempath, asize)
1200  if (asize > 0) then
1201  call mem_reallocate(this%decay_sorbed, this%dis%nodes, &
1202  'DECAY_SORBED', this%memoryPath)
1203  end if
1204  if (this%isrb < 2) then
1205  call get_isize('SP2', this%input_mempath, asize)
1206  if (asize > 0) &
1207  call mem_reallocate(this%sp2, this%dis%nodes, 'SP2', this%memoryPath)
1208  end if
1209  !
1210  ! -- update defaults with memory sourced values
1211  call mem_set_value(this%porosity, 'POROSITY', this%input_mempath, map, &
1212  found%porosity)
1213  call mem_set_value(this%volfrac, 'VOLFRAC', this%input_mempath, map, &
1214  found%volfrac)
1215  call mem_set_value(this%zetaim, 'ZETAIM', this%input_mempath, map, &
1216  found%zetaim)
1217  call mem_set_value(this%cim, 'CIM', this%input_mempath, map, &
1218  found%cim)
1219  call mem_set_value(this%decay, 'DECAY', this%input_mempath, map, &
1220  found%decay)
1221  call mem_set_value(this%decay_sorbed, 'DECAY_SORBED', this%input_mempath, &
1222  map, found%decay_sorbed)
1223  call mem_set_value(this%bulk_density, 'BULK_DENSITY', this%input_mempath, &
1224  map, found%bulk_density)
1225  call mem_set_value(this%distcoef, 'DISTCOEF', this%input_mempath, map, &
1226  found%distcoef)
1227  call mem_set_value(this%sp2, 'SP2', this%input_mempath, map, &
1228  found%sp2)
1229 
1230  ! -- log options
1231  if (this%iout > 0) then
1232  call this%log_data(found)
1233  end if
1234 
1235  ! -- Check for required sorption variables
1236  if (this%isrb > 0) then
1237  if (.not. found%bulk_density) then
1238  write (errmsg, '(a)') 'Sorption is active but BULK_DENSITY &
1239  &not specified. BULK_DENSITY must be specified in GRIDDATA block.'
1240  call store_error(errmsg)
1241  end if
1242  if (.not. found%distcoef) then
1243  write (errmsg, '(a)') 'Sorption is active but distribution &
1244  &coefficient not specified. DISTCOEF must be specified in &
1245  &GRIDDATA block.'
1246  call store_error(errmsg)
1247  end if
1248  if (this%isrb > 1) then
1249  if (.not. found%sp2) then
1250  write (errmsg, '(a)') 'Freundlich or langmuir sorption is active &
1251  &but SP2 not specified. SP2 must be specified in &
1252  &GRIDDATA block.'
1253  call store_error(errmsg)
1254  end if
1255  end if
1256  else
1257  if (found%bulk_density) then
1258  write (warnmsg, '(a)') 'Sorption is not active but &
1259  &BULK_DENSITY was specified. BULK_DENSITY will have no affect on &
1260  &simulation results.'
1261  call store_warning(warnmsg)
1262  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1263  end if
1264  if (found%distcoef) then
1265  write (warnmsg, '(a)') 'Sorption is not active but &
1266  &distribution coefficient was specified. DISTCOEF will have &
1267  &no affect on simulation results.'
1268  call store_warning(warnmsg)
1269  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1270  end if
1271  if (found%sp2) then
1272  write (warnmsg, '(a)') 'Sorption is not active but &
1273  &SP2 was specified. SP2 will have &
1274  &no affect on simulation results.'
1275  call store_warning(warnmsg)
1276  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1277  end if
1278  end if
1279 
1280  ! -- Check for required decay/production rate coefficients
1281  if (this%idcy > 0) then
1282  if (.not. found%decay) then
1283  write (errmsg, '(a)') 'First or zero order decay is &
1284  &active but the first rate coefficient is not specified. DECAY &
1285  &must be specified in GRIDDATA block.'
1286  call store_error(errmsg)
1287  end if
1288  if (.not. found%decay_sorbed) then
1289  !
1290  ! -- If DECAY_SORBED not specified and sorption is active, then
1291  ! terminate with an error
1292  write (errmsg, '(a)') 'DECAY_SORBED not provided in GRIDDATA &
1293  &block but decay and sorption are active. Specify DECAY_SORBED &
1294  &in GRIDDATA block.'
1295  call store_error(errmsg)
1296  end if
1297  else
1298  if (found%decay) then
1299  write (warnmsg, '(a)') 'First- or zero-order decay &
1300  &is not active but decay was specified. DECAY will &
1301  &have no affect on simulation results.'
1302  call store_warning(warnmsg)
1303  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1304  end if
1305  if (found%decay_sorbed) then
1306  write (warnmsg, '(a)') 'First- or zero-order decay &
1307  &is not active but DECAY_SORBED was specified. &
1308  &DECAY_SORBED will have no affect on simulation results.'
1309  call store_warning(warnmsg)
1310  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1311  end if
1312  end if
1313 
1314  ! -- Check for required dual domain arrays or warn if they are specified
1315  ! but won't be used.
1316  if (.not. found%cim) then
1317  write (this%iout, '(1x,a)') 'Warning. Dual domain is active but &
1318  &initial immobile domain concentration was not specified. &
1319  &Setting CIM to zero.'
1320  end if
1321  if (.not. found%zetaim) then
1322  write (errmsg, '(a)') 'Dual domain is active but dual &
1323  &domain mass transfer rate (ZETAIM) was not specified. ZETAIM &
1324  &must be specified in GRIDDATA block.'
1325  call store_error(errmsg)
1326  end if
1327  if (.not. found%porosity) then
1328  write (errmsg, '(a)') 'Dual domain is active but &
1329  &immobile domain POROSITY was not specified. POROSITY &
1330  &must be specified in GRIDDATA block.'
1331  call store_error(errmsg)
1332  end if
1333  if (.not. found%volfrac) then
1334  write (errmsg, '(a)') 'Dual domain is active but &
1335  &immobile domain VOLFRAC was not specified. VOLFRAC &
1336  &must be specified in GRIDDATA block. This is a new &
1337  &requirement for MODFLOW versions later than version &
1338  &6.4.1.'
1339  call store_error(errmsg)
1340  end if
1341 
1342  ! -- terminate if errors
1343  if (count_errors() > 0) then
1344  call store_error_filename(this%input_fname)
1345  end if
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
character(len=maxcharlen) errmsg
error message string
character(len=maxcharlen) warnmsg
warning message string
Here is the call graph for this function:

◆ source_options()

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

Method to source options for the package.

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

983  ! -- modules
986  use openspecmodule, only: access, form
990  ! -- dummy
991  class(GwtIstType), intent(inout) :: this
992  ! -- locals
993  character(len=LINELENGTH) :: prnfmt
994  integer(I4B), pointer :: columns, width, digits
995  type(GwtIstParamFoundType) :: found
996  character(len=LENVARNAME), dimension(3) :: sorption_method = &
997  &[character(len=LENVARNAME) :: 'LINEAR', 'FREUNDLICH', 'LANGMUIR']
998  character(len=linelength) :: sorbate_fname, cim6_fname, budget_fname, &
999  budgetcsv_fname
1000  allocate (columns)
1001  allocate (width)
1002  allocate (digits)
1003  !
1004  ! -- update defaults with memory sourced values
1005  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', this%input_mempath, &
1006  found%save_flows)
1007  call mem_set_value(budget_fname, 'BUDGETFILE', this%input_mempath, &
1008  found%budgetfile)
1009  call mem_set_value(budgetcsv_fname, 'BUDGETCSVFILE', this%input_mempath, &
1010  found%budgetcsvfile)
1011  call mem_set_value(this%isrb, 'SORPTION', this%input_mempath, &
1012  sorption_method, found%sorption)
1013  call mem_set_value(this%idcy, 'ORDER1_DECAY', this%input_mempath, &
1014  found%order1_decay)
1015  call mem_set_value(this%idcy, 'ORDER0_DECAY', this%input_mempath, &
1016  found%order0_decay)
1017  call mem_set_value(cim6_fname, 'CIMFILE', this%input_mempath, &
1018  found%cimfile)
1019  call mem_set_value(sorbate_fname, 'SORBATEFILE', this%input_mempath, &
1020  found%sorbatefile)
1021  call mem_set_value(columns, 'COLUMNS', this%input_mempath, &
1022  found%columns)
1023  call mem_set_value(width, 'WIDTH', this%input_mempath, &
1024  found%width)
1025  call mem_set_value(digits, 'DIGITS', this%input_mempath, &
1026  found%digits)
1027  call mem_set_value(prnfmt, 'FORMAT', this%input_mempath, &
1028  found%format)
1029 
1030  ! -- found side effects
1031  if (found%save_flows) this%ipakcb = -1
1032  if (found%budgetfile) then
1033  call assign_iounit(this%ibudgetout, this%inunit, "BUDGET fileout")
1034  call openfile(this%ibudgetout, this%iout, budget_fname, 'DATA(BINARY)', &
1035  form, access, 'REPLACE', mode_opt=mnormal)
1036  end if
1037  if (found%budgetcsvfile) then
1038  call assign_iounit(this%ibudcsv, this%inunit, "BUDGETCSV fileout")
1039  call openfile(this%ibudcsv, this%iout, budgetcsv_fname, 'CSV', &
1040  filstat_opt='REPLACE')
1041  end if
1042  if (found%sorption) then
1043  if (this%isrb == 0) then
1044  call store_error('Unknown sorption type was specified. &
1045  &Sorption must be specified as LINEAR, &
1046  &FREUNDLICH, or LANGMUIR.')
1047  call store_error_filename(this%input_fname)
1048  end if
1049  end if
1050  if (found%order1_decay) this%idcy = 1
1051  if (found%order0_decay) this%idcy = 2
1052  if (found%cimfile) then
1053  call this%ocd%set_option('FILEOUT '//trim(cim6_fname)//" ", 0, &
1054  this%iout)
1055  end if
1056  if (found%columns .and. found%width .and. &
1057  found%digits .and. found%format) then
1058  write (this%lstfmt, '(a,i0,a,i0,a,i0,a)') 'PRINT_FORMAT COLUMNS ', &
1059  columns, ' WIDTH ', width, ' DIGITS ', digits, ' '//trim(prnfmt)
1060  end if
1061  if (found%sorbatefile) then
1062  this%ioutsorbate = getunit()
1063  call openfile(this%ioutsorbate, this%iout, sorbate_fname, &
1064  'DATA(BINARY)', form, access, 'REPLACE', mode_opt=mnormal)
1065  end if
1066  !
1067  ! -- log options
1068  if (this%iout > 0) then
1069  call this%log_options(found, cim6_fname, budget_fname, &
1070  budgetcsv_fname, sorbate_fname)
1071  end if
1072 
1073  deallocate (columns)
1074  deallocate (width)
1075  deallocate (digits)
@ mnormal
normal output mode
Definition: Constants.f90:206
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
subroutine, public assign_iounit(iounit, errunit, description)
@ brief assign io unit number
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
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 38 of file gwt-ist.f90.

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

◆ ftype

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

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

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

◆ nbditems

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

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

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

◆ text

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

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

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