MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
gwtmstmodule Module Reference

– @ brief Mobile Storage and Transfer (MST) Module More...

Data Types

type  gwtmsttype
 @ brief Mobile storage and transfer More...
 

Enumerations

enum  { decay_off = 0 , decay_first_order = 1 , decay_zero_order = 2 }
 Enumerator that defines the decay options. More...
 

Functions/Subroutines

subroutine, public mst_cr (mstobj, name_model, input_mempath, inunit, iout, fmi)
 @ brief Create a new package object More...
 
subroutine mst_ar (this, dis, ibound)
 @ brief Allocate and read method for package More...
 
subroutine mst_fc (this, nodes, cold, nja, matrix_sln, idxglo, cnew, rhs, kiter)
 @ brief Fill coefficient method for package More...
 
subroutine mst_fc_sto (this, nodes, cold, nja, matrix_sln, idxglo, rhs)
 @ brief Fill storage coefficient method for package More...
 
subroutine mst_fc_dcy (this, nodes, cold, cnew, nja, matrix_sln, idxglo, rhs, kiter)
 @ brief Fill decay coefficient method for package More...
 
subroutine mst_fc_srb (this, nodes, cold, nja, matrix_sln, idxglo, rhs, cnew)
 @ brief Fill sorption coefficient method for package More...
 
subroutine mst_fc_dcy_srb (this, nodes, cold, nja, matrix_sln, idxglo, rhs, cnew, kiter)
 @ brief Fill sorption-decay coefficient method for package More...
 
subroutine mst_cq (this, nodes, cnew, cold, flowja)
 @ brief Calculate flows for package More...
 
subroutine mst_cq_sto (this, nodes, cnew, cold, flowja)
 @ brief Calculate storage terms for package More...
 
subroutine mst_cq_dcy (this, nodes, cnew, cold, flowja)
 @ brief Calculate decay terms for package More...
 
subroutine mst_cq_srb (this, nodes, cnew, cold, flowja)
 @ brief Calculate sorption terms for package More...
 
subroutine mst_cq_dcy_srb (this, nodes, cnew, cold, flowja)
 @ brief Calculate decay-sorption terms for package More...
 
subroutine mst_calc_csrb (this, cnew)
 @ brief Calculate sorbed concentration More...
 
subroutine mst_bd (this, isuppress_output, model_budget)
 @ brief Calculate budget terms for package More...
 
subroutine mst_ot_flow (this, icbcfl, icbcun)
 @ brief Output flow terms for package More...
 
subroutine mst_ot_dv (this, idvsave)
 Save sorbate concentration array to binary file. More...
 
subroutine mst_da (this)
 @ brief Deallocate More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalar variables for package More...
 
subroutine allocate_arrays (this, nodes)
 @ brief Allocate arrays for package More...
 
subroutine source_options (this)
 @ brief Source options for package More...
 
subroutine log_options (this, found, sorbate_fname)
 Log user options to list file. More...
 
subroutine source_data (this)
 @ brief Source data for package More...
 
subroutine log_data (this, found)
 Log user data to list file. More...
 
subroutine addto_volfracim (this, volfracim)
 @ brief Add volfrac values to volfracim More...
 
real(dp) function get_volfracm (this, node)
 @ brief Return mobile domain volume fraction More...
 
real(dp) function get_zero_order_decay (decay_rate_usr, decay_rate_last, kiter, cold, cnew, delt)
 @ brief Calculate zero-order decay rate and constrain if necessary More...
 

Variables

integer(i4b), parameter nbditems = 4
 
character(len=lenbudtxt), dimension(nbditemsbudtxt
 

Detailed Description

The GwtMstModule contains the GwtMstType, which is the derived type responsible for adding the effects of

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

Enumeration Type Documentation

◆ anonymous enum

anonymous enum
Enumerator
decay_off 

Decay (or production) of mass inactive (default)

decay_first_order 

First-order decay.

decay_zero_order 

Zeroth-order decay.

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

Function/Subroutine Documentation

◆ addto_volfracim()

subroutine gwtmstmodule::addto_volfracim ( class(gwtmsttype this,
real(dp), dimension(:), intent(in)  volfracim 
)

Method to add immobile domain volume fracions, which are stored as a cumulative value in volfracim.

Parameters
thisGwtMstType object
[in]volfracimimmobile domain volume fraction that contributes to total immobile volume fraction

Definition at line 1407 of file gwt-mst.f90.

1408  ! -- dummy
1409  class(GwtMstType) :: this !< GwtMstType object
1410  real(DP), dimension(:), intent(in) :: volfracim !< immobile domain volume fraction that contributes to total immobile volume fraction
1411  ! -- local
1412  integer(I4B) :: n
1413  !
1414  ! -- Add to volfracim
1415  do n = 1, this%dis%nodes
1416  this%volfracim(n) = this%volfracim(n) + volfracim(n)
1417  end do
1418  !
1419  ! -- An immobile domain is adding a volume fraction, so update thetam
1420  ! accordingly.
1421  do n = 1, this%dis%nodes
1422  this%thetam(n) = this%get_volfracm(n) * this%porosity(n)
1423  end do

◆ allocate_arrays()

subroutine gwtmstmodule::allocate_arrays ( class(gwtmsttype this,
integer(i4b), intent(in)  nodes 
)

Method to allocate arrays for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes

Definition at line 1014 of file gwt-mst.f90.

1015  ! -- modules
1017  use constantsmodule, only: dzero
1018  ! -- dummy
1019  class(GwtMstType) :: this !< GwtMstType object
1020  integer(I4B), intent(in) :: nodes !< number of nodes
1021  ! -- local
1022  integer(I4B) :: n
1023  !
1024  ! -- Allocate
1025  ! -- sto
1026  call mem_allocate(this%porosity, nodes, 'POROSITY', this%memoryPath)
1027  call mem_allocate(this%thetam, nodes, 'THETAM', this%memoryPath)
1028  call mem_allocate(this%volfracim, nodes, 'VOLFRACIM', this%memoryPath)
1029  call mem_allocate(this%ratesto, nodes, 'RATESTO', this%memoryPath)
1030  !
1031  ! -- dcy
1032  if (this%idcy == decay_off) then
1033  call mem_allocate(this%ratedcy, 1, 'RATEDCY', this%memoryPath)
1034  call mem_allocate(this%decay, 1, 'DECAY', this%memoryPath)
1035  call mem_allocate(this%decaylast, 1, 'DECAYLAST', this%memoryPath)
1036  else
1037  call mem_allocate(this%ratedcy, this%dis%nodes, 'RATEDCY', this%memoryPath)
1038  call mem_allocate(this%decay, nodes, 'DECAY', this%memoryPath)
1039  call mem_allocate(this%decaylast, nodes, 'DECAYLAST', this%memoryPath)
1040  end if
1041  if (this%idcy /= decay_off .and. this%isrb /= sorption_off) then
1042  call mem_allocate(this%ratedcys, this%dis%nodes, 'RATEDCYS', &
1043  this%memoryPath)
1044  call mem_allocate(this%decayslast, this%dis%nodes, 'DECAYSLAST', &
1045  this%memoryPath)
1046  else
1047  call mem_allocate(this%ratedcys, 1, 'RATEDCYS', this%memoryPath)
1048  call mem_allocate(this%decayslast, 1, 'DECAYSLAST', this%memoryPath)
1049  end if
1050  call mem_allocate(this%decay_sorbed, 1, 'DECAY_SORBED', &
1051  this%memoryPath)
1052  !
1053  ! -- srb
1054  if (this%isrb == sorption_off) then
1055  call mem_allocate(this%bulk_density, 1, 'BULK_DENSITY', this%memoryPath)
1056  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
1057  call mem_allocate(this%distcoef, 1, 'DISTCOEF', this%memoryPath)
1058  call mem_allocate(this%ratesrb, 1, 'RATESRB', this%memoryPath)
1059  call mem_allocate(this%csrb, 1, 'CSRB', this%memoryPath)
1060  else
1061  call mem_allocate(this%bulk_density, nodes, 'BULK_DENSITY', this%memoryPath)
1062  call mem_allocate(this%distcoef, nodes, 'DISTCOEF', this%memoryPath)
1063  call mem_allocate(this%ratesrb, nodes, 'RATESRB', this%memoryPath)
1064  call mem_allocate(this%csrb, nodes, 'CSRB', this%memoryPath)
1065  if (this%isrb == sorption_linear) then
1066  call mem_allocate(this%sp2, 1, 'SP2', this%memoryPath)
1067  else
1068  call mem_allocate(this%sp2, nodes, 'SP2', this%memoryPath)
1069  end if
1070  end if
1071  !
1072  ! -- Initialize
1073  do n = 1, nodes
1074  this%porosity(n) = dzero
1075  this%thetam(n) = dzero
1076  this%volfracim(n) = dzero
1077  this%ratesto(n) = dzero
1078  end do
1079  do n = 1, size(this%decay)
1080  this%decay(n) = dzero
1081  this%ratedcy(n) = dzero
1082  this%decaylast(n) = dzero
1083  end do
1084  do n = 1, size(this%bulk_density)
1085  this%bulk_density(n) = dzero
1086  this%distcoef(n) = dzero
1087  this%ratesrb(n) = dzero
1088  this%csrb(n) = dzero
1089  end do
1090  do n = 1, size(this%sp2)
1091  this%sp2(n) = dzero
1092  end do
1093  do n = 1, size(this%ratedcys)
1094  this%ratedcys(n) = dzero
1095  this%decayslast(n) = dzero
1096  end do
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65

◆ allocate_scalars()

subroutine gwtmstmodule::allocate_scalars ( class(gwtmsttype this)

Method to allocate scalar variables for the package.

Parameters
thisGwtMstType object

Definition at line 990 of file gwt-mst.f90.

991  ! -- modules
993  ! -- dummy
994  class(GwtMstType) :: this !< GwtMstType object
995  !
996  ! -- Allocate scalars in NumericalPackageType
997  call this%NumericalPackageType%allocate_scalars()
998  !
999  ! -- Allocate
1000  call mem_allocate(this%isrb, 'ISRB', this%memoryPath)
1001  call mem_allocate(this%ioutsorbate, 'IOUTSORBATE', this%memoryPath)
1002  call mem_allocate(this%idcy, 'IDCY', this%memoryPath)
1003  !
1004  ! -- Initialize
1005  this%isrb = izero
1006  this%ioutsorbate = 0
1007  this%idcy = izero

◆ get_volfracm()

real(dp) function gwtmstmodule::get_volfracm ( class(gwtmsttype this,
integer(i4b), intent(in)  node 
)

Calculate and return the volume fraction of the aquifer that is mobile

Parameters
thisGwtMstType object
[in]nodenode number

Definition at line 1430 of file gwt-mst.f90.

1431  ! -- dummy
1432  class(GwtMstType) :: this !< GwtMstType object
1433  integer(I4B), intent(in) :: node !< node number
1434  ! -- return
1435  real(DP) :: volfracm
1436  !
1437  volfracm = done - this%volfracim(node)

◆ get_zero_order_decay()

real(dp) function gwtmstmodule::get_zero_order_decay ( real(dp), intent(in)  decay_rate_usr,
real(dp), intent(in)  decay_rate_last,
integer(i4b), intent(in)  kiter,
real(dp), intent(in)  cold,
real(dp), intent(in)  cnew,
real(dp), intent(in)  delt 
)

Function to calculate the zero-order decay rate from the user specified decay rate. If the decay rate is positive, then the decay rate must be constrained so that more mass is not removed than is available. Without this constraint, negative concentrations could result from zero-order decay.

Parameters
[in]decay_rate_usruser-entered decay rate
[in]decay_rate_lastdecay rate used for last iteration
[in]kiterPicard iteration counter
[in]coldconcentration at end of last time step
[in]cnewconcentration at end of this time step
[in]deltlength of time step
Returns
returned value for decay rate

Definition at line 1448 of file gwt-mst.f90.

1450  ! -- dummy
1451  real(DP), intent(in) :: decay_rate_usr !< user-entered decay rate
1452  real(DP), intent(in) :: decay_rate_last !< decay rate used for last iteration
1453  integer(I4B), intent(in) :: kiter !< Picard iteration counter
1454  real(DP), intent(in) :: cold !< concentration at end of last time step
1455  real(DP), intent(in) :: cnew !< concentration at end of this time step
1456  real(DP), intent(in) :: delt !< length of time step
1457  ! -- return
1458  real(DP) :: decay_rate !< returned value for decay rate
1459  !
1460  ! -- Return user rate if production, otherwise constrain, if necessary
1461  if (decay_rate_usr < dzero) then
1462  !
1463  ! -- Production, no need to limit rate
1464  decay_rate = decay_rate_usr
1465  else
1466  !
1467  ! -- Need to ensure decay does not result in negative
1468  ! concentration, so reduce the rate if it would result in
1469  ! removing more mass than is in the cell.
1470  if (kiter == 1) then
1471  decay_rate = min(decay_rate_usr, cold / delt)
1472  else
1473  decay_rate = decay_rate_last
1474  if (cnew < dzero) then
1475  decay_rate = decay_rate_last + cnew / delt
1476  else if (cnew > cold) then
1477  decay_rate = decay_rate_last + cnew / delt
1478  end if
1479  decay_rate = min(decay_rate_usr, decay_rate)
1480  end if
1481  decay_rate = max(decay_rate, dzero)
1482  end if
Here is the caller graph for this function:

◆ log_data()

subroutine gwtmstmodule::log_data ( class(gwtmsttype this,
type(gwtmstparamfoundtype), intent(in)  found 
)

Definition at line 1375 of file gwt-mst.f90.

1377  class(GwTMstType) :: this
1378  type(GwtMstParamFoundType), intent(in) :: found
1379 
1380  write (this%iout, '(1x,a)') 'PROCESSING GRIDDATA'
1381  if (found%porosity) then
1382  write (this%iout, '(4x,a)') 'MOBILE DOMAIN POROSITY set from input file'
1383  end if
1384  if (found%decay) then
1385  write (this%iout, '(4x,a)') 'DECAY RATE set from input file'
1386  end if
1387  if (found%decay_sorbed) then
1388  write (this%iout, '(4x,a)') 'DECAY SORBED RATE set from input file'
1389  end if
1390  if (found%bulk_density) then
1391  write (this%iout, '(4x,a)') 'BULK DENSITY set from input file'
1392  end if
1393  if (found%distcoef) then
1394  write (this%iout, '(4x,a)') 'DISTRIBUTION COEFFICIENT set from input file'
1395  end if
1396  if (found%sp2) then
1397  write (this%iout, '(4x,a)') 'SECOND SORPTION PARAM set from input file'
1398  end if
1399  write (this%iout, '(1x,a)') 'END PROCESSING GRIDDATA'

◆ log_options()

subroutine gwtmstmodule::log_options ( class(gwtmsttype this,
type(gwtmstparamfoundtype), intent(in)  found,
character(len=*), intent(in)  sorbate_fname 
)

Definition at line 1156 of file gwt-mst.f90.

1158  class(GwTMstType) :: this
1159  type(GwtMstParamFoundType), intent(in) :: found
1160  character(len=*), intent(in) :: sorbate_fname
1161  ! -- formats
1162  character(len=*), parameter :: fmtisvflow = &
1163  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
1164  &WHENEVER ICBCFL IS NOT ZERO.')"
1165  character(len=*), parameter :: fmtlinear = &
1166  &"(4x,'LINEAR SORPTION IS ACTIVE. ')"
1167  character(len=*), parameter :: fmtfreundlich = &
1168  &"(4x,'FREUNDLICH SORPTION IS ACTIVE. ')"
1169  character(len=*), parameter :: fmtlangmuir = &
1170  &"(4x,'LANGMUIR SORPTION IS ACTIVE. ')"
1171  character(len=*), parameter :: fmtidcy1 = &
1172  &"(4x,'FIRST-ORDER DECAY IS ACTIVE. ')"
1173  character(len=*), parameter :: fmtidcy2 = &
1174  &"(4x,'ZERO-ORDER DECAY IS ACTIVE. ')"
1175  character(len=*), parameter :: fmtfileout = &
1176  "(4x,'MST ',1x,a,1x,' WILL BE SAVED TO FILE: ',a,/4x,&
1177  &'OPENED ON UNIT: ',I7)"
1178 
1179  write (this%iout, '(1x,a)') 'PROCESSING MOBILE STORAGE AND TRANSFER OPTIONS'
1180  if (found%save_flows) then
1181  write (this%iout, fmtisvflow)
1182  end if
1183  if (found%order1_decay) then
1184  write (this%iout, fmtidcy1)
1185  end if
1186  if (found%order0_decay) then
1187  write (this%iout, fmtidcy2)
1188  end if
1189  if (found%sorption) then
1190  select case (this%isrb)
1191  case (sorption_linear)
1192  write (this%iout, fmtlinear)
1193  case (sorption_freund)
1194  write (this%iout, fmtfreundlich)
1195  case (sorption_lang)
1196  write (this%iout, fmtlangmuir)
1197  end select
1198  end if
1199  if (found%sorbatefile) then
1200  write (this%iout, fmtfileout) &
1201  'SORBATE', sorbate_fname, this%ioutsorbate
1202  end if
1203  write (this%iout, '(1x,a)') 'END OF MOBILE STORAGE AND TRANSFER OPTIONS'

◆ mst_ar()

subroutine gwtmstmodule::mst_ar ( class(gwtmsttype), intent(inout)  this,
class(disbasetype), intent(in), pointer  dis,
integer(i4b), dimension(:), pointer, contiguous  ibound 
)

Method to allocate and read static data for the package.

Parameters
[in,out]thisGwtMstType object
[in]dispointer to dis package
iboundpointer to GWT ibound array

Definition at line 142 of file gwt-mst.f90.

143  ! -- modules
144  ! -- dummy
145  class(GwtMstType), intent(inout) :: this !< GwtMstType object
146  class(DisBaseType), pointer, intent(in) :: dis !< pointer to dis package
147  integer(I4B), dimension(:), pointer, contiguous :: ibound !< pointer to GWT ibound array
148  ! -- local
149  ! -- formats
150  character(len=*), parameter :: fmtmst = &
151  "(1x,/1x,'MST -- MOBILE STORAGE AND TRANSFER PACKAGE, VERSION 1, &
152  &7/29/2020 INPUT READ FROM MEMPATH: ', A, //)"
153  !
154  ! --print a message identifying the mobile storage and transfer package.
155  write (this%iout, fmtmst) this%input_mempath
156  !
157  ! -- Source options
158  call this%source_options()
159  !
160  ! -- store pointers to arguments that were passed in
161  this%dis => dis
162  this%ibound => ibound
163  !
164  ! -- Allocate arrays
165  call this%allocate_arrays(dis%nodes)
166  !
167  ! -- source the data block
168  call this%source_data()
169  !
170  ! -- Create isotherm object if sorption is active
171  this%isotherm => create_isotherm(this%isrb, this%distcoef, this%sp2)
Here is the call graph for this function:

◆ mst_bd()

subroutine gwtmstmodule::mst_bd ( class(gwtmsttype this,
integer(i4b), intent(in)  isuppress_output,
type(budgettype), intent(inout)  model_budget 
)
Parameters
thisGwtMstType object
[in]isuppress_outputflag to suppress output
[in,out]model_budgetmodel budget object

Definition at line 810 of file gwt-mst.f90.

811  ! -- modules
812  use tdismodule, only: delt
814  ! -- dummy
815  class(GwtMstType) :: this !< GwtMstType object
816  integer(I4B), intent(in) :: isuppress_output !< flag to suppress output
817  type(BudgetType), intent(inout) :: model_budget !< model budget object
818  ! -- local
819  real(DP) :: rin
820  real(DP) :: rout
821  !
822  ! -- sto
823  call rate_accumulator(this%ratesto, rin, rout)
824  call model_budget%addentry(rin, rout, delt, budtxt(1), &
825  isuppress_output, rowlabel=this%packName)
826  !
827  ! -- dcy
828  if (this%idcy /= decay_off) then
829  call rate_accumulator(this%ratedcy, rin, rout)
830  call model_budget%addentry(rin, rout, delt, budtxt(2), &
831  isuppress_output, rowlabel=this%packName)
832  end if
833  !
834  ! -- srb
835  if (this%isrb /= sorption_off) then
836  call rate_accumulator(this%ratesrb, rin, rout)
837  call model_budget%addentry(rin, rout, delt, budtxt(3), &
838  isuppress_output, rowlabel=this%packName)
839  end if
840  !
841  ! -- srb dcy
842  if (this%isrb /= sorption_off .and. this%idcy /= decay_off) then
843  call rate_accumulator(this%ratedcys, rin, rout)
844  call model_budget%addentry(rin, rout, delt, budtxt(4), &
845  isuppress_output, rowlabel=this%packName)
846  end if
This module contains the BudgetModule.
Definition: Budget.f90:20
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:

◆ mst_calc_csrb()

subroutine gwtmstmodule::mst_calc_csrb ( class(gwtmsttype this,
real(dp), dimension(:), intent(in)  cnew 
)
Parameters
thisGwtMstType object
[in]cnewconcentration at end of this time step

Definition at line 789 of file gwt-mst.f90.

790  ! -- dummy
791  class(GwtMstType) :: this !< GwtMstType object
792  real(DP), intent(in), dimension(:) :: cnew !< concentration at end of this time step
793  ! -- local
794  integer(I4B) :: n
795  real(DP) :: csrb
796 
797  ! Calculate sorbed concentration
798  do n = 1, size(cnew)
799  csrb = dzero
800  if (this%ibound(n) > 0 .and. this%isrb /= sorption_off) then
801  csrb = this%isotherm%value(cnew, n)
802  end if
803  this%csrb(n) = csrb
804  end do
805 

◆ mst_cq()

subroutine gwtmstmodule::mst_cq ( class(gwtmsttype this,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cnew,
real(dp), dimension(nodes), intent(in)  cold,
real(dp), dimension(:), intent(inout), contiguous  flowja 
)

Method to calculate flows for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]cnewconcentration at end of this time step
[in]coldconcentration at end of last time step
[in,out]flowjaflow between two connected control volumes

Definition at line 482 of file gwt-mst.f90.

483  ! -- dummy
484  class(GwtMstType) :: this !< GwtMstType object
485  integer(I4B), intent(in) :: nodes !< number of nodes
486  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
487  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
488  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
489  !
490  ! - storage
491  call this%mst_cq_sto(nodes, cnew, cold, flowja)
492  !
493  ! -- decay
494  if (this%idcy /= decay_off) then
495  call this%mst_cq_dcy(nodes, cnew, cold, flowja)
496  end if
497  !
498  ! -- sorption
499  if (this%isrb /= sorption_off) then
500  call this%mst_cq_srb(nodes, cnew, cold, flowja)
501  end if
502  !
503  ! -- decay sorbed
504  if (this%isrb /= sorption_off .and. this%idcy /= decay_off) then
505  call this%mst_cq_dcy_srb(nodes, cnew, cold, flowja)
506  end if
507  !
508  ! -- calculate csrb
509  if (this%isrb /= sorption_off) then
510  call this%mst_calc_csrb(cnew)
511  end if

◆ mst_cq_dcy()

subroutine gwtmstmodule::mst_cq_dcy ( class(gwtmsttype this,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cnew,
real(dp), dimension(nodes), intent(in)  cold,
real(dp), dimension(:), intent(inout), contiguous  flowja 
)

Method to calculate decay terms for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]cnewconcentration at end of this time step
[in]coldconcentration at end of last time step
[in,out]flowjaflow between two connected control volumes

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

567  ! -- modules
568  use tdismodule, only: delt
569  ! -- dummy
570  class(GwtMstType) :: this !< GwtMstType object
571  integer(I4B), intent(in) :: nodes !< number of nodes
572  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
573  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
574  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
575  ! -- local
576  integer(I4B) :: n
577  integer(I4B) :: idiag
578  real(DP) :: rate
579  real(DP) :: swtpdt
580  real(DP) :: hhcof, rrhs
581  real(DP) :: vcell
582  real(DP) :: decay_rate
583  !
584  ! -- initialize
585  !
586  ! -- Calculate decay change
587  do n = 1, nodes
588  !
589  ! -- skip if transport inactive
590  this%ratedcy(n) = dzero
591  if (this%ibound(n) <= 0) cycle
592  !
593  ! -- calculate new and old water volumes
594  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
595  swtpdt = this%fmi%gwfsat(n)
596  !
597  ! -- calculate decay gains and losses
598  rate = dzero
599  hhcof = dzero
600  rrhs = dzero
601  if (this%idcy == decay_first_order) then
602  if (cnew(n) > dzero) then
603  hhcof = -this%decay(n) * vcell * swtpdt * this%thetam(n)
604  end if
605  elseif (this%idcy == decay_zero_order) then
606  decay_rate = get_zero_order_decay(this%decay(n), this%decaylast(n), &
607  0, cold(n), cnew(n), delt)
608  rrhs = decay_rate * vcell * swtpdt * this%thetam(n)
609  end if
610  rate = hhcof * cnew(n) - rrhs
611  this%ratedcy(n) = rate
612  idiag = this%dis%con%ia(n)
613  flowja(idiag) = flowja(idiag) + rate
614  !
615  end do
Here is the call graph for this function:

◆ mst_cq_dcy_srb()

subroutine gwtmstmodule::mst_cq_dcy_srb ( class(gwtmsttype this,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cnew,
real(dp), dimension(nodes), intent(in)  cold,
real(dp), dimension(:), intent(inout), contiguous  flowja 
)

Method to calculate decay-sorption terms for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]cnewconcentration at end of this time step
[in]coldconcentration at end of last time step
[in,out]flowjaflow between two connected control volumes

Definition at line 694 of file gwt-mst.f90.

695  ! -- modules
696  use tdismodule, only: delt
697  ! -- dummy
698  class(GwtMstType) :: this !< GwtMstType object
699  integer(I4B), intent(in) :: nodes !< number of nodes
700  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
701  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
702  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
703  ! -- local
704  integer(I4B) :: n
705  integer(I4B) :: idiag
706  real(DP) :: rate
707  real(DP) :: hhcof, rrhs
708  real(DP) :: vcell
709  real(DP) :: swnew
710  real(DP) :: distcoef
711  real(DP) :: volfracm
712  real(DP) :: rhobm
713  real(DP) :: term
714  real(DP) :: csrb
715  real(DP) :: csrbnew
716  real(DP) :: csrbold
717  real(DP) :: decay_rate
718  !
719  ! -- Calculate sorbed decay change
720  ! This routine will only be called if sorption and decay are active
721  do n = 1, nodes
722  !
723  ! -- initialize rates
724  this%ratedcys(n) = dzero
725  !
726  ! -- skip if transport inactive
727  if (this%ibound(n) <= 0) cycle
728  !
729  ! -- set variables
730  hhcof = dzero
731  rrhs = dzero
732  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
733  swnew = this%fmi%gwfsat(n)
734  distcoef = this%distcoef(n)
735  volfracm = this%get_volfracm(n)
736  rhobm = this%bulk_density(n)
737  term = this%decay_sorbed(n) * volfracm * rhobm * swnew * vcell
738  !
739  ! -- add sorbed mass decay rate terms to accumulators
740  select case (this%idcy)
741  case (decay_first_order)
742 
743  !
744  select case (this%isrb)
745  case (sorption_linear)
746  !
747  ! -- first order decay rate is a function of concentration, so add
748  ! to left hand side
749  if (cnew(n) > dzero) then
750  hhcof = -term * distcoef
751  end if
752  case (sorption_freund)
753  !
754  ! -- nonlinear Freundlich sorption, so add to RHS
755  csrb = this%isotherm%value(cnew, n)
756  rrhs = term * csrb
757  case (sorption_lang)
758  !
759  ! -- nonlinear Lanmuir sorption, so add to RHS
760  csrb = this%isotherm%value(cnew, n)
761  rrhs = term * csrb
762  end select
763  case (decay_zero_order)
764  !
765  ! -- Call function to get zero-order decay rate, which may be changed
766  ! from the user-specified rate to prevent negative concentrations
767  if (distcoef > dzero) then
768  csrbold = this%isotherm%value(cold, n)
769  csrbnew = this%isotherm%value(cnew, n)
770 
771  decay_rate = get_zero_order_decay(this%decay_sorbed(n), &
772  this%decayslast(n), &
773  0, csrbold, csrbnew, delt)
774  rrhs = decay_rate * volfracm * rhobm * swnew * vcell
775  end if
776  end select
777  !
778  ! -- calculate rate
779  rate = hhcof * cnew(n) - rrhs
780  this%ratedcys(n) = rate
781  idiag = this%dis%con%ia(n)
782  flowja(idiag) = flowja(idiag) + rate
783  !
784  end do
Here is the call graph for this function:

◆ mst_cq_srb()

subroutine gwtmstmodule::mst_cq_srb ( class(gwtmsttype this,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cnew,
real(dp), dimension(nodes), intent(in)  cold,
real(dp), dimension(:), intent(inout), contiguous  flowja 
)

Method to calculate sorption terms for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]cnewconcentration at end of this time step
[in]coldconcentration at end of last time step
[in,out]flowjaflow between two connected control volumes

Definition at line 622 of file gwt-mst.f90.

623  ! -- modules
624  use tdismodule, only: delt
625  ! -- dummy
626  class(GwtMstType) :: this !< GwtMstType object
627  integer(I4B), intent(in) :: nodes !< number of nodes
628  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
629  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
630  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
631  ! -- local
632  integer(I4B) :: n
633  integer(I4B) :: idiag
634  real(DP) :: rate
635  real(DP) :: tled
636  real(DP) :: vcell
637  real(DP) :: volfracm
638  real(DP) :: rhobm
639  real(DP) :: sat_new, sat_old
640  real(DP) :: contribution
641  !
642  ! -- initialize
643  tled = done / delt
644  !
645  ! -- Calculate sorption change
646  do n = 1, nodes
647  !
648  ! -- initialize rates
649  this%ratesrb(n) = dzero
650  rate = 0.0_dp
651  !
652  ! -- skip if transport inactive
653  if (this%ibound(n) <= 0) cycle
654  !
655  ! -- assign variables
656  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
657 
658  rhobm = this%bulk_density(n)
659  volfracm = this%get_volfracm(n)
660  sat_new = this%fmi%gwfsat(n)
661  sat_old = this%fmi%gwfsatold(n, delt)
662 
663  ! -- Matrix contribution for sorption term
664  contribution = -volfracm * rhobm * sat_new &
665  * this%isotherm%derivative(cnew, n) * cnew(n) * vcell * tled
666  rate = rate + contribution
667 
668  ! -- Right-hand side contribution due to linearization
669  ! -- Note: this contribution should cancel with the matrix contribution when the outer loop is converged
670  contribution = -volfracm * rhobm * sat_new * &
671  this%isotherm%derivative(cnew, n) * cnew(n) * vcell * tled
672  rate = rate - contribution
673 
674  contribution = volfracm * rhobm * sat_new * this%isotherm%value(cnew, n) &
675  * vcell * tled
676  rate = rate - contribution
677 
678  ! -- Right-hand side contribution from previous time step
679  contribution = -volfracm * rhobm * sat_old * this%isotherm%value(cold, n) &
680  * vcell * tled
681  rate = rate - contribution
682 
683  this%ratesrb(n) = rate
684  idiag = this%dis%con%ia(n)
685  flowja(idiag) = flowja(idiag) + rate
686  !
687  end do

◆ mst_cq_sto()

subroutine gwtmstmodule::mst_cq_sto ( class(gwtmsttype this,
integer(i4b), intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cnew,
real(dp), dimension(nodes), intent(in)  cold,
real(dp), dimension(:), intent(inout), contiguous  flowja 
)

Method to calculate storage terms for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]cnewconcentration at end of this time step
[in]coldconcentration at end of last time step
[in,out]flowjaflow between two connected control volumes

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

519  ! -- modules
520  use tdismodule, only: delt
521  ! -- dummy
522  class(GwtMstType) :: this !< GwtMstType object
523  integer(I4B), intent(in) :: nodes !< number of nodes
524  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
525  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
526  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< flow between two connected control volumes
527  ! -- local
528  integer(I4B) :: n
529  integer(I4B) :: idiag
530  real(DP) :: rate
531  real(DP) :: tled
532  real(DP) :: vnew, vold
533  real(DP) :: hhcof, rrhs
534  !
535  ! -- initialize
536  tled = done / delt
537  !
538  ! -- Calculate storage change
539  do n = 1, nodes
540  this%ratesto(n) = dzero
541  !
542  ! -- skip if transport inactive
543  if (this%ibound(n) <= 0) cycle
544  !
545  ! -- calculate new and old water volumes
546  vnew = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n)) * &
547  this%fmi%gwfsat(n) * this%thetam(n)
548  vold = vnew
549  if (this%fmi%igwfstrgss /= 0) vold = vold + this%fmi%gwfstrgss(n) * delt
550  if (this%fmi%igwfstrgsy /= 0) vold = vold + this%fmi%gwfstrgsy(n) * delt
551  !
552  ! -- calculate rate
553  hhcof = -vnew * tled
554  rrhs = -vold * tled * cold(n)
555  rate = hhcof * cnew(n) - rrhs
556  this%ratesto(n) = rate
557  idiag = this%dis%con%ia(n)
558  flowja(idiag) = flowja(idiag) + rate
559  end do

◆ mst_cr()

subroutine, public gwtmstmodule::mst_cr ( type(gwtmsttype), pointer  mstobj,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  input_mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
type(tspfmitype), intent(in), target  fmi 
)

Create a new MST object

Parameters
mstobjunallocated new mst object to create
[in]name_modelname of the model
[in]input_mempathmemory path of input
[in]inunitunit number of WEL package input file
[in]ioutunit number of model listing file
[in]fmifmi package for this GWT model

Definition at line 114 of file gwt-mst.f90.

115  ! -- dummy
116  type(GwtMstType), pointer :: mstobj !< unallocated new mst object to create
117  character(len=*), intent(in) :: name_model !< name of the model
118  character(len=*), intent(in) :: input_mempath !< memory path of input
119  integer(I4B), intent(in) :: inunit !< unit number of WEL package input file
120  integer(I4B), intent(in) :: iout !< unit number of model listing file
121  type(TspFmiType), intent(in), target :: fmi !< fmi package for this GWT model
122  !
123  ! -- Create the object
124  allocate (mstobj)
125  !
126  ! -- create name and memory path
127  call mstobj%set_names(1, name_model, 'MST', 'MST', input_mempath)
128  !
129  ! -- Allocate scalars
130  call mstobj%allocate_scalars()
131  !
132  ! -- Set variables
133  mstobj%inunit = inunit
134  mstobj%iout = iout
135  mstobj%fmi => fmi
Here is the caller graph for this function:

◆ mst_da()

subroutine gwtmstmodule::mst_da ( class(gwtmsttype this)

Method to deallocate memory for the package.

Parameters
thisGwtMstType object

Definition at line 944 of file gwt-mst.f90.

945  ! -- modules
947  ! -- dummy
948  class(GwtMstType) :: this !< GwtMstType object
949  !
950  ! -- Deallocate arrays if package was active
951  if (this%inunit > 0) then
952  call mem_deallocate(this%porosity)
953  call mem_deallocate(this%thetam)
954  call mem_deallocate(this%volfracim)
955  call mem_deallocate(this%ratesto)
956  call mem_deallocate(this%idcy)
957  call mem_deallocate(this%decay)
958  call mem_deallocate(this%decay_sorbed)
959  call mem_deallocate(this%ratedcy)
960  call mem_deallocate(this%decaylast)
961  call mem_deallocate(this%decayslast)
962  call mem_deallocate(this%isrb)
963  call mem_deallocate(this%ioutsorbate)
964  call mem_deallocate(this%bulk_density)
965  call mem_deallocate(this%distcoef)
966  call mem_deallocate(this%sp2)
967  call mem_deallocate(this%ratesrb)
968  call mem_deallocate(this%csrb)
969  call mem_deallocate(this%ratedcys)
970  this%ibound => null()
971  this%fmi => null()
972  end if
973  !
974  ! -- Scalars
975  !
976  ! -- Objects
977  if (associated(this%isotherm)) then
978  deallocate (this%isotherm)
979  nullify (this%isotherm)
980  end if
981  !
982  ! -- deallocate parent
983  call this%NumericalPackageType%da()

◆ mst_fc()

subroutine gwtmstmodule::mst_fc ( class(gwtmsttype this,
integer, intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cold,
integer(i4b), intent(in)  nja,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(nja), intent(in)  idxglo,
real(dp), dimension(nodes), intent(in)  cnew,
real(dp), dimension(nodes), intent(inout)  rhs,
integer(i4b), intent(in)  kiter 
)

Method to calculate and fill coefficients for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]coldconcentration at end of last time step
[in]njanumber of GWT connections
matrix_slnsolution matrix
[in]idxglomapping vector for model (local) to solution (global)
[in,out]rhsright-hand side vector for model
[in]cnewconcentration at end of this time step
[in]kitersolution outer iteration number

Definition at line 178 of file gwt-mst.f90.

180  ! -- modules
181  ! -- dummy
182  class(GwtMstType) :: this !< GwtMstType object
183  integer, intent(in) :: nodes !< number of nodes
184  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
185  integer(I4B), intent(in) :: nja !< number of GWT connections
186  class(MatrixBaseType), pointer :: matrix_sln !< solution matrix
187  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
188  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
189  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
190  integer(I4B), intent(in) :: kiter !< solution outer iteration number
191  !
192  ! -- storage contribution
193  call this%mst_fc_sto(nodes, cold, nja, matrix_sln, idxglo, rhs)
194  !
195  ! -- decay contribution
196  if (this%idcy /= decay_off) then
197  call this%mst_fc_dcy(nodes, cold, cnew, nja, matrix_sln, idxglo, &
198  rhs, kiter)
199  end if
200  !
201  ! -- sorption contribution
202  if (this%isrb /= sorption_off) then
203  call this%mst_fc_srb(nodes, cold, nja, matrix_sln, idxglo, rhs, cnew)
204  end if
205  !
206  ! -- decay sorbed contribution
207  if (this%isrb /= sorption_off .and. this%idcy /= decay_off) then
208  call this%mst_fc_dcy_srb(nodes, cold, nja, matrix_sln, idxglo, rhs, &
209  cnew, kiter)
210  end if

◆ mst_fc_dcy()

subroutine gwtmstmodule::mst_fc_dcy ( class(gwtmsttype this,
integer, intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cold,
real(dp), dimension(nodes), intent(in)  cnew,
integer(i4b), intent(in)  nja,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(nja), intent(in)  idxglo,
real(dp), dimension(nodes), intent(inout)  rhs,
integer(i4b), intent(in)  kiter 
)

Method to calculate and fill decay coefficients for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]coldconcentration at end of last time step
[in]cnewconcentration at end of this time step
[in]njanumber of GWT connections
matrix_slnsolution coefficient matrix
[in]idxglomapping vector for model (local) to solution (global)
[in,out]rhsright-hand side vector for model
[in]kitersolution outer iteration number

Definition at line 263 of file gwt-mst.f90.

265  ! -- modules
266  use tdismodule, only: delt
267  ! -- dummy
268  class(GwtMstType) :: this !< GwtMstType object
269  integer, intent(in) :: nodes !< number of nodes
270  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
271  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
272  integer(I4B), intent(in) :: nja !< number of GWT connections
273  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
274  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
275  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
276  integer(I4B), intent(in) :: kiter !< solution outer iteration number
277  ! -- local
278  integer(I4B) :: n, idiag
279  real(DP) :: hhcof, rrhs
280  real(DP) :: swtpdt
281  real(DP) :: vcell
282  real(DP) :: decay_rate
283  !
284  ! -- loop through and calculate decay contribution to hcof and rhs
285  do n = 1, this%dis%nodes
286  !
287  ! -- skip if transport inactive
288  if (this%ibound(n) <= 0) cycle
289  !
290  ! -- calculate new and old water volumes
291  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
292  swtpdt = this%fmi%gwfsat(n)
293  !
294  ! -- add decay rate terms to accumulators
295  idiag = this%dis%con%ia(n)
296  select case (this%idcy)
297  case (decay_first_order)
298  !
299  ! -- first order decay rate is a function of concentration, so add
300  ! to left hand side
301  if (cnew(n) > dzero) then
302  hhcof = -this%decay(n) * vcell * swtpdt * this%thetam(n)
303  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
304  end if
305  case (decay_zero_order)
306  !
307  ! -- Call function to get zero-order decay rate, which may be changed
308  ! from the user-specified rate to prevent negative concentrations
309  decay_rate = get_zero_order_decay(this%decay(n), this%decaylast(n), &
310  kiter, cold(n), cnew(n), delt)
311  this%decaylast(n) = decay_rate
312  rrhs = decay_rate * vcell * swtpdt * this%thetam(n)
313  rhs(n) = rhs(n) + rrhs
314  end select
315  !
316  end do
Here is the call graph for this function:

◆ mst_fc_dcy_srb()

subroutine gwtmstmodule::mst_fc_dcy_srb ( class(gwtmsttype this,
integer, intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cold,
integer(i4b), intent(in)  nja,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(nja), intent(in)  idxglo,
real(dp), dimension(nodes), intent(inout)  rhs,
real(dp), dimension(nodes), intent(in)  cnew,
integer(i4b), intent(in)  kiter 
)

Method to calculate and fill sorption-decay coefficients for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]coldconcentration at end of last time step
[in]njanumber of GWT connections
matrix_slnsolution coefficient matrix
[in]idxglomapping vector for model (local) to solution (global)
[in,out]rhsright-hand side vector for model
[in]cnewconcentration at end of this time step
[in]kitersolution outer iteration number

Definition at line 388 of file gwt-mst.f90.

390  ! -- modules
391  use tdismodule, only: delt
392  ! -- dummy
393  class(GwtMstType) :: this !< GwtMstType object
394  integer, intent(in) :: nodes !< number of nodes
395  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
396  integer(I4B), intent(in) :: nja !< number of GWT connections
397  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
398  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
399  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
400  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
401  integer(I4B), intent(in) :: kiter !< solution outer iteration number
402  ! -- local
403  integer(I4B) :: n, idiag
404  real(DP) :: hhcof, rrhs
405  real(DP) :: vcell
406  real(DP) :: swnew
407  real(DP) :: distcoef
408  real(DP) :: volfracm
409  real(DP) :: rhobm
410  real(DP) :: term
411  real(DP) :: csrb
412  real(DP) :: decay_rate
413  real(DP) :: csrbold
414  real(DP) :: csrbnew
415  !
416  ! -- loop through and calculate sorption contribution to hcof and rhs
417  do n = 1, this%dis%nodes
418  !
419  ! -- skip if transport inactive
420  if (this%ibound(n) <= 0) cycle
421  !
422  ! -- set variables
423  hhcof = dzero
424  rrhs = dzero
425  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
426  swnew = this%fmi%gwfsat(n)
427  distcoef = this%distcoef(n)
428  idiag = this%dis%con%ia(n)
429  volfracm = this%get_volfracm(n)
430  rhobm = this%bulk_density(n)
431  term = this%decay_sorbed(n) * volfracm * rhobm * swnew * vcell
432  !
433  ! -- add sorbed mass decay rate terms to accumulators
434  select case (this%idcy)
435  case (decay_first_order)
436  select case (this%isrb)
437  case (sorption_linear)
438  !
439  ! -- first order decay rate is a function of concentration, so add
440  ! to left hand side
441  if (cnew(n) > dzero) then
442  hhcof = -term * distcoef
443  end if
444  case (sorption_freund)
445  !
446  ! -- nonlinear Freundlich sorption, so add to RHS
447  csrb = this%isotherm%value(cnew, n)
448  rrhs = term * csrb
449  case (sorption_lang)
450  !
451  ! -- nonlinear Lanmuir sorption, so add to RHS
452  csrb = this%isotherm%value(cnew, n)
453  rrhs = term * csrb
454  end select
455  case (decay_zero_order)
456  !
457  ! -- call function to get zero-order decay rate, which may be changed
458  ! from the user-specified rate to prevent negative concentrations
459  if (distcoef > dzero) then
460  csrbold = this%isotherm%value(cold, n)
461  csrbnew = this%isotherm%value(cnew, n)
462  !
463  decay_rate = get_zero_order_decay(this%decay_sorbed(n), &
464  this%decayslast(n), &
465  kiter, csrbold, csrbnew, delt)
466  this%decayslast(n) = decay_rate
467  rrhs = decay_rate * volfracm * rhobm * swnew * vcell
468  end if
469  end select
470  !
471  ! -- Add hhcof to diagonal and rrhs to right-hand side
472  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
473  rhs(n) = rhs(n) + rrhs
474  !
475  end do
Here is the call graph for this function:

◆ mst_fc_srb()

subroutine gwtmstmodule::mst_fc_srb ( class(gwtmsttype this,
integer, intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cold,
integer(i4b), intent(in)  nja,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(nja), intent(in)  idxglo,
real(dp), dimension(nodes), intent(inout)  rhs,
real(dp), dimension(nodes), intent(in)  cnew 
)

Method to calculate and fill sorption coefficients for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]coldconcentration at end of last time step
[in]njanumber of GWT connections
matrix_slnsolution coefficient matrix
[in]idxglomapping vector for model (local) to solution (global)
[in,out]rhsright-hand side vector for model
[in]cnewconcentration at end of this time step

Definition at line 323 of file gwt-mst.f90.

325  ! -- modules
326  use tdismodule, only: delt
327  ! -- dummy
328  class(GwtMstType) :: this !< GwtMstType object
329  integer, intent(in) :: nodes !< number of nodes
330  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
331  integer(I4B), intent(in) :: nja !< number of GWT connections
332  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
333  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
334  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
335  real(DP), intent(in), dimension(nodes) :: cnew !< concentration at end of this time step
336  ! -- local
337  integer(I4B) :: n, idiag
338  real(DP) :: tled
339  real(DP) :: hhcof, rrhs
340  real(DP) :: vcell
341  real(DP) :: volfracm
342  real(DP) :: rhobm
343  real(DP) :: sat_new, sat_old
344  !
345  ! -- set variables
346  tled = done / delt
347  !
348  ! -- loop through and calculate sorption contribution to hcof and rhs
349  do n = 1, this%dis%nodes
350  !
351  ! -- skip if transport inactive
352  if (this%ibound(n) <= 0) cycle
353  !
354  ! -- assign variables
355  vcell = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n))
356  volfracm = this%get_volfracm(n)
357  rhobm = this%bulk_density(n)
358  sat_new = this%fmi%gwfsat(n)
359  sat_old = this%fmi%gwfsatold(n, delt)
360 
361  ! -- Matrix contribution for sorption term
362  hhcof = -volfracm * rhobm * sat_new * this%isotherm%derivative(cnew, n) &
363  * vcell * tled
364  idiag = this%dis%con%ia(n)
365  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
366 
367  ! -- Right-hand side contribution due to linearization
368  rrhs = -volfracm * rhobm * sat_new * this%isotherm%derivative(cnew, n) &
369  * cnew(n) * vcell * tled
370  rhs(n) = rhs(n) + rrhs
371 
372  rrhs = volfracm * rhobm * sat_new * this%isotherm%value(cnew, n) * &
373  vcell * tled
374  rhs(n) = rhs(n) + rrhs
375 
376  ! -- Right-hand side contribution from previous time step
377  rrhs = -volfracm * rhobm * sat_old * this%isotherm%value(cold, n) * &
378  vcell * tled
379  rhs(n) = rhs(n) + rrhs
380 
381  end do

◆ mst_fc_sto()

subroutine gwtmstmodule::mst_fc_sto ( class(gwtmsttype this,
integer, intent(in)  nodes,
real(dp), dimension(nodes), intent(in)  cold,
integer(i4b), intent(in)  nja,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(nja), intent(in)  idxglo,
real(dp), dimension(nodes), intent(inout)  rhs 
)

Method to calculate and fill storage coefficients for the package.

Parameters
thisGwtMstType object
[in]nodesnumber of nodes
[in]coldconcentration at end of last time step
[in]njanumber of GWT connections
matrix_slnsolution coefficient matrix
[in]idxglomapping vector for model (local) to solution (global)
[in,out]rhsright-hand side vector for model

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

218  ! -- modules
219  use tdismodule, only: delt
220  ! -- dummy
221  class(GwtMstType) :: this !< GwtMstType object
222  integer, intent(in) :: nodes !< number of nodes
223  real(DP), intent(in), dimension(nodes) :: cold !< concentration at end of last time step
224  integer(I4B), intent(in) :: nja !< number of GWT connections
225  class(MatrixBaseType), pointer :: matrix_sln !< solution coefficient matrix
226  integer(I4B), intent(in), dimension(nja) :: idxglo !< mapping vector for model (local) to solution (global)
227  real(DP), intent(inout), dimension(nodes) :: rhs !< right-hand side vector for model
228  ! -- local
229  integer(I4B) :: n, idiag
230  real(DP) :: tled
231  real(DP) :: hhcof, rrhs
232  real(DP) :: vnew, vold
233  !
234  ! -- set variables
235  tled = done / delt
236  !
237  ! -- loop through and calculate storage contribution to hcof and rhs
238  do n = 1, this%dis%nodes
239  !
240  ! -- skip if transport inactive
241  if (this%ibound(n) <= 0) cycle
242  !
243  ! -- calculate new and old water volumes
244  vnew = this%dis%area(n) * (this%dis%top(n) - this%dis%bot(n)) * &
245  this%fmi%gwfsat(n) * this%thetam(n)
246  vold = vnew
247  if (this%fmi%igwfstrgss /= 0) vold = vold + this%fmi%gwfstrgss(n) * delt
248  if (this%fmi%igwfstrgsy /= 0) vold = vold + this%fmi%gwfstrgsy(n) * delt
249  !
250  ! -- add terms to diagonal and rhs accumulators
251  hhcof = -vnew * tled
252  rrhs = -vold * tled * cold(n)
253  idiag = this%dis%con%ia(n)
254  call matrix_sln%add_value_pos(idxglo(idiag), hhcof)
255  rhs(n) = rhs(n) + rrhs
256  end do

◆ mst_ot_dv()

subroutine gwtmstmodule::mst_ot_dv ( class(gwtmsttype this,
integer(i4b), intent(in)  idvsave 
)

Definition at line 906 of file gwt-mst.f90.

907  ! -- dummy
908  class(GwtMstType) :: this
909  integer(I4B), intent(in) :: idvsave
910  ! -- local
911  character(len=1) :: cdatafmp = ' ', editdesc = ' '
912  integer(I4B) :: ibinun
913  integer(I4B) :: iprint
914  integer(I4B) :: nvaluesp
915  integer(I4B) :: nwidthp
916  real(DP) :: dinact
917 
918  ! Set unit number for sorbate output
919  if (this%ioutsorbate /= 0) then
920  ibinun = 1
921  else
922  ibinun = 0
923  end if
924  if (idvsave == 0) ibinun = 0
925 
926  ! save sorbate concentration array
927  if (ibinun /= 0) then
928  iprint = 0
929  dinact = dhnoflo
930  if (this%ioutsorbate /= 0) then
931  ibinun = this%ioutsorbate
932  call this%dis%record_array(this%csrb, this%iout, iprint, ibinun, &
933  ' SORBATE', cdatafmp, nvaluesp, &
934  nwidthp, editdesc, dinact)
935  end if
936  end if
937 

◆ mst_ot_flow()

subroutine gwtmstmodule::mst_ot_flow ( class(gwtmsttype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  icbcun 
)

Method to output terms for the package.

Parameters
thisGwtMstType object
[in]icbcflflag and unit number for cell-by-cell output
[in]icbcunflag indication if cell-by-cell data should be saved

Definition at line 853 of file gwt-mst.f90.

854  ! -- dummy
855  class(GwtMstType) :: this !< GwtMstType object
856  integer(I4B), intent(in) :: icbcfl !< flag and unit number for cell-by-cell output
857  integer(I4B), intent(in) :: icbcun !< flag indication if cell-by-cell data should be saved
858  ! -- local
859  integer(I4B) :: ibinun
860  integer(I4B) :: iprint, nvaluesp, nwidthp
861  character(len=1) :: cdatafmp = ' ', editdesc = ' '
862  real(DP) :: dinact
863  !
864  ! -- Set unit number for binary output
865  if (this%ipakcb < 0) then
866  ibinun = icbcun
867  elseif (this%ipakcb == 0) then
868  ibinun = 0
869  else
870  ibinun = this%ipakcb
871  end if
872  if (icbcfl == 0) ibinun = 0
873  !
874  ! -- Record the storage rate if requested
875  if (ibinun /= 0) then
876  iprint = 0
877  dinact = dzero
878  !
879  ! -- sto
880  call this%dis%record_array(this%ratesto, this%iout, iprint, -ibinun, &
881  budtxt(1), cdatafmp, nvaluesp, &
882  nwidthp, editdesc, dinact)
883  !
884  ! -- dcy
885  if (this%idcy /= decay_off) &
886  call this%dis%record_array(this%ratedcy, this%iout, iprint, -ibinun, &
887  budtxt(2), cdatafmp, nvaluesp, &
888  nwidthp, editdesc, dinact)
889  !
890  ! -- srb
891  if (this%isrb /= sorption_off) &
892  call this%dis%record_array(this%ratesrb, this%iout, iprint, -ibinun, &
893  budtxt(3), cdatafmp, nvaluesp, &
894  nwidthp, editdesc, dinact)
895  !
896  ! -- dcy srb
897  if (this%isrb /= sorption_off .and. this%idcy /= decay_off) &
898  call this%dis%record_array(this%ratedcys, this%iout, iprint, -ibinun, &
899  budtxt(4), cdatafmp, nvaluesp, &
900  nwidthp, editdesc, dinact)
901  end if

◆ source_data()

subroutine gwtmstmodule::source_data ( class(gwtmsttype this)

Method to source data for the package.

Definition at line 1210 of file gwt-mst.f90.

1211  ! -- modules
1215  ! -- dummy
1216  class(GwtMsttype) :: this
1217  ! -- locals
1218  character(len=LINELENGTH) :: errmsg
1219  type(GwtMstParamFoundType) :: found
1220  integer(I4B) :: n, asize
1221  integer(I4B), dimension(:), pointer, contiguous :: map
1222  !
1223  ! -- set map to convert user input data into reduced data
1224  map => null()
1225  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
1226  !
1227  ! -- reallocate
1228  if (this%isrb == sorption_off) then
1229  call get_isize('BULK_DENSITY', this%input_mempath, asize)
1230  if (asize > 0) &
1231  call mem_reallocate(this%bulk_density, this%dis%nodes, &
1232  'BULK_DENSITY', this%memoryPath)
1233  call get_isize('DISTCOEF', this%input_mempath, asize)
1234  if (asize > 0) &
1235  call mem_reallocate(this%distcoef, this%dis%nodes, 'DISTCOEF', &
1236  this%memoryPath)
1237  end if
1238  if (this%idcy == decay_off) then
1239  call get_isize('DECAY', this%input_mempath, asize)
1240  if (asize > 0) &
1241  call mem_reallocate(this%decay, this%dis%nodes, 'DECAY', this%memoryPath)
1242  end if
1243  call get_isize('DECAY_SORBED', this%input_mempath, asize)
1244  if (asize > 0) then
1245  call mem_reallocate(this%decay_sorbed, this%dis%nodes, &
1246  'DECAY_SORBED', this%memoryPath)
1247  end if
1248  if (this%isrb == sorption_off .or. this%isrb == sorption_linear) then
1249  call get_isize('SP2', this%input_mempath, asize)
1250  if (asize > 0) &
1251  call mem_reallocate(this%sp2, this%dis%nodes, 'SP2', this%memoryPath)
1252  end if
1253  !
1254  ! -- update defaults with memory sourced values
1255  call mem_set_value(this%porosity, 'POROSITY', this%input_mempath, map, &
1256  found%porosity)
1257  call mem_set_value(this%decay, 'DECAY', this%input_mempath, map, &
1258  found%decay)
1259  call mem_set_value(this%decay_sorbed, 'DECAY_SORBED', this%input_mempath, &
1260  map, found%decay_sorbed)
1261  call mem_set_value(this%bulk_density, 'BULK_DENSITY', this%input_mempath, &
1262  map, found%bulk_density)
1263  call mem_set_value(this%distcoef, 'DISTCOEF', this%input_mempath, map, &
1264  found%distcoef)
1265  call mem_set_value(this%sp2, 'SP2', this%input_mempath, map, &
1266  found%sp2)
1267 
1268  ! -- log options
1269  if (this%iout > 0) then
1270  call this%log_data(found)
1271  end if
1272 
1273  ! -- Check for required porosity
1274  if (.not. found%porosity) then
1275  write (errmsg, '(a)') 'POROSITY not specified in GRIDDATA block.'
1276  call store_error(errmsg)
1277  end if
1278 
1279  ! -- Check for required sorption variables
1280  if (this%isrb == sorption_linear .or. this%isrb == sorption_freund .or. &
1281  this%isrb == sorption_lang) then
1282  if (.not. found%bulk_density) then
1283  write (errmsg, '(a)') 'Sorption is active but BULK_DENSITY &
1284  &not specified. BULK_DENSITY must be specified in GRIDDATA block.'
1285  call store_error(errmsg)
1286  end if
1287  if (.not. found%distcoef) then
1288  write (errmsg, '(a)') 'Sorption is active but distribution &
1289  &coefficient not specified. DISTCOEF must be specified in &
1290  &GRIDDATA block.'
1291  call store_error(errmsg)
1292  end if
1293  if (this%isrb == sorption_freund .or. this%isrb == sorption_lang) then
1294  if (.not. found%sp2) then
1295  write (errmsg, '(a)') 'Freundlich or langmuir sorption is active &
1296  &but SP2 not specified. SP2 must be specified in &
1297  &GRIDDATA block.'
1298  call store_error(errmsg)
1299  end if
1300  end if
1301  else
1302  if (found%bulk_density) then
1303  write (warnmsg, '(a)') 'Sorption is not active but &
1304  &BULK_DENSITY was specified. BULK_DENSITY will have no affect on &
1305  &simulation results.'
1306  call store_warning(warnmsg)
1307  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1308  end if
1309  if (found%distcoef) then
1310  write (warnmsg, '(a)') 'Sorption is not active but &
1311  &distribution coefficient was specified. DISTCOEF will have &
1312  &no affect on simulation results.'
1313  call store_warning(warnmsg)
1314  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1315  end if
1316  if (found%sp2) then
1317  write (warnmsg, '(a)') 'Sorption is not active but &
1318  &SP2 was specified. SP2 will have &
1319  &no affect on simulation results.'
1320  call store_warning(warnmsg)
1321  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1322  end if
1323  end if
1324 
1325  ! -- Check for required decay/production rate coefficients
1326  if (this%idcy /= decay_off) then
1327  if (.not. found%decay) then
1328  write (errmsg, '(a)') 'First or zero order decay is &
1329  &active but the first rate coefficient is not specified. DECAY &
1330  &must be specified in GRIDDATA block.'
1331  call store_error(errmsg)
1332  end if
1333  if (.not. found%decay_sorbed) then
1334  !
1335  ! -- If DECAY_SORBED not specified and sorption is active, then
1336  ! terminate with an error
1337  if (this%isrb == sorption_linear .or. this%isrb == sorption_freund .or. &
1338  this%isrb == sorption_lang) then
1339  write (errmsg, '(a)') 'DECAY_SORBED not provided in GRIDDATA &
1340  &block but decay and sorption are active. Specify DECAY_SORBED &
1341  &in GRIDDATA block.'
1342  call store_error(errmsg)
1343  end if
1344  end if
1345  else
1346  if (found%decay) then
1347  write (warnmsg, '(a)') 'First- or zero-order decay &
1348  &is not active but decay was specified. DECAY will &
1349  &have no affect on simulation results.'
1350  call store_warning(warnmsg)
1351  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1352  end if
1353  if (found%decay_sorbed) then
1354  write (warnmsg, '(a)') 'First- or zero-order decay &
1355  &is not active but DECAY_SORBED was specified. &
1356  &DECAY_SORBED will have no affect on simulation results.'
1357  call store_warning(warnmsg)
1358  write (this%iout, '(1x,a)') 'WARNING. '//warnmsg
1359  end if
1360  end if
1361 
1362  ! -- terminate if errors
1363  if (count_errors() > 0) then
1364  call store_error_filename(this%input_fname)
1365  end if
1366 
1367  ! -- initialize thetam from porosity
1368  do n = 1, size(this%porosity)
1369  this%thetam(n) = this%porosity(n)
1370  end do
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
Here is the call graph for this function:

◆ source_options()

subroutine gwtmstmodule::source_options ( class(gwtmsttype this)

Method to source options for the package.

Definition at line 1103 of file gwt-mst.f90.

1104  ! -- modules
1105  use constantsmodule, only: lenvarname
1106  use openspecmodule, only: access, form
1107  use inputoutputmodule, only: getunit, openfile
1110  ! -- dummy
1111  class(GwtMstType) :: this
1112  ! -- locals
1113  type(GwtMstParamFoundType) :: found
1114  character(len=LENVARNAME), dimension(3) :: sorption_method = &
1115  &[character(len=LENVARNAME) :: 'LINEAR', 'FREUNDLICH', 'LANGMUIR']
1116  character(len=linelength) :: fname
1117  !
1118  ! -- update defaults with memory sourced values
1119  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', this%input_mempath, &
1120  found%save_flows)
1121  call mem_set_value(this%idcy, 'ORDER1_DECAY', this%input_mempath, &
1122  found%order1_decay)
1123  call mem_set_value(this%idcy, 'ORDER0_DECAY', this%input_mempath, &
1124  found%order0_decay)
1125  call mem_set_value(this%isrb, 'SORPTION', this%input_mempath, &
1126  sorption_method, found%sorption)
1127  call mem_set_value(fname, 'SORBATEFILE', this%input_mempath, &
1128  found%sorbatefile)
1129 
1130  ! -- found side effects
1131  if (found%save_flows) this%ipakcb = -1
1132  if (found%order1_decay) this%idcy = decay_first_order
1133  if (found%order0_decay) this%idcy = decay_zero_order
1134  if (found%sorption) then
1135  if (this%isrb == izero) then
1136  call store_error('Unknown sorption type was specified. &
1137  &Sorption must be specified as LINEAR, &
1138  &FREUNDLICH, or LANGMUIR.')
1139  call store_error_filename(this%input_fname)
1140  end if
1141  end if
1142  if (found%sorbatefile) then
1143  this%ioutsorbate = getunit()
1144  call openfile(this%ioutsorbate, this%iout, fname, 'DATA(BINARY)', &
1145  form, access, 'REPLACE', mode_opt=mnormal)
1146  end if
1147  !
1148  ! -- log options
1149  if (this%iout > 0) then
1150  call this%log_options(found, fname)
1151  end if
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
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) gwtmstmodule::budtxt

Definition at line 31 of file gwt-mst.f90.

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

◆ nbditems

integer(i4b), parameter gwtmstmodule::nbditems = 4

Definition at line 30 of file gwt-mst.f90.

30  integer(I4B), parameter :: NBDITEMS = 4