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

This module contains the storage package methods. More...

Data Types

type  gwfstotype
 

Functions/Subroutines

subroutine, public sto_cr (stoobj, name_model, mempath, inunit, iout)
 @ brief Create a new package object More...
 
subroutine sto_ar (this, dis, ibound)
 @ brief Allocate and read method for package More...
 
subroutine sto_rp (this)
 @ brief Read and prepare method for package More...
 
subroutine sto_ad (this)
 @ brief Advance the package More...
 
subroutine sto_fc (this, kiter, hold, hnew, matrix_sln, idxglo, rhs)
 @ brief Fill A and right-hand side for the package More...
 
subroutine sto_fn (this, kiter, hold, hnew, matrix_sln, idxglo, rhs)
 @ brief Fill Newton-Raphson terms in A and right-hand side for the package More...
 
subroutine sto_cq (this, flowja, hnew, hold)
 @ brief Calculate flows for package More...
 
subroutine sto_bd (this, isuppress_output, model_budget)
 @ brief Model budget calculation for package More...
 
subroutine sto_save_model_flows (this, icbcfl, icbcun)
 @ brief Save model flows for package More...
 
subroutine sto_da (this)
 @ brief Deallocate package memory More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine allocate_arrays (this, nodes)
 @ brief Allocate package arrays More...
 
subroutine source_options (this)
 @ brief Source input options for package More...
 
subroutine log_options (this, found)
 @ brief Log found options for package More...
 
subroutine source_data (this)
 @ brief Source input data for package More...
 
subroutine save_old_ss_sy (this)
 @ brief Save old storage property values More...
 

Variables

character(len=lenbudtxt), dimension(2) budtxt = [' STO-SS', ' STO-SY']
 

Detailed Description

This module contains the methods used to add the effects of storage on the groundwater flow equation. The contribution of specific storage and specific yield can be represented.

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwfstomodule::allocate_arrays ( class(gwfstotype), target  this,
integer(i4b), intent(in)  nodes 
)

Allocate and initialize STO package arrays.

Parameters
thisGwfStoType object
[in]nodesactive model nodes

Definition at line 732 of file gwf-sto.f90.

733  ! -- modules
735  ! -- dummy variables
736  class(GwfStoType), target :: this !< GwfStoType object
737  integer(I4B), intent(in) :: nodes !< active model nodes
738  ! -- local variables
739  integer(I4B) :: n
740  !
741  ! -- Allocate arrays
742  call mem_allocate(this%iconvert, nodes, 'ICONVERT', this%memoryPath)
743  call mem_allocate(this%ss, nodes, 'SS', this%memoryPath)
744  call mem_allocate(this%sy, nodes, 'SY', this%memoryPath)
745  call mem_allocate(this%strgss, nodes, 'STRGSS', this%memoryPath)
746  call mem_allocate(this%strgsy, nodes, 'STRGSY', this%memoryPath)
747  !
748  ! -- set input context pointers
749  if (this%inunit > 0) then
750  call mem_setptr(this%iper, 'IPER', this%input_mempath)
751  call mem_setptr(this%storage, 'STORAGE', this%input_mempath)
752  end if
753  !
754  ! -- Initialize arrays
755  this%iss = 0
756  do n = 1, nodes
757  this%iconvert(n) = 1
758  this%ss(n) = dzero
759  this%sy(n) = dzero
760  this%strgss(n) = dzero
761  this%strgsy(n) = dzero
762  if (this%integratechanges /= 0) then
763  this%oldss(n) = dzero
764  if (this%iusesy /= 0) then
765  this%oldsy(n) = dzero
766  end if
767  end if
768  end do

◆ allocate_scalars()

subroutine gwfstomodule::allocate_scalars ( class(gwfstotype this)

Allocate and initialize scalars for the STO package. The base numerical package allocate scalars method is also called.

Parameters
thisGwfStoType object

Definition at line 698 of file gwf-sto.f90.

699  ! -- modules
701  ! -- dummy variables
702  class(GwfStoType) :: this !< GwfStoType object
703  !
704  ! -- allocate scalars in NumericalPackageType
705  call this%NumericalPackageType%allocate_scalars()
706  !
707  ! -- allocate scalars
708  call mem_allocate(this%istor_coef, 'ISTOR_COEF', this%memoryPath)
709  call mem_allocate(this%iconf_ss, 'ICONF_SS', this%memoryPath)
710  call mem_allocate(this%iorig_ss, 'IORIG_SS', this%memoryPath)
711  call mem_allocate(this%iusesy, 'IUSESY', this%memoryPath)
712  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
713  call mem_allocate(this%integratechanges, 'INTEGRATECHANGES', &
714  this%memoryPath)
715  call mem_allocate(this%intvs, 'INTVS', this%memoryPath)
716  !
717  ! -- initialize scalars
718  this%istor_coef = 0
719  this%iconf_ss = 0
720  this%iorig_ss = 0
721  this%iusesy = 0
722  this%satomega = dzero
723  this%integratechanges = 0
724  this%intvs = 0

◆ log_options()

subroutine gwfstomodule::log_options ( class(gwfstotype this,
type(gwfstoparamfoundtype), intent(in)  found 
)

Log options block for STO package.

Parameters
thisGwfStoType object

Definition at line 838 of file gwf-sto.f90.

839  ! -- modules
842  ! -- dummy variables
843  class(GwfStoType) :: this !< GwfStoType object
844  type(GwfStoParamFoundType), intent(in) :: found
845  ! -- local variables
846  ! -- formats
847  character(len=*), parameter :: fmtisvflow = &
848  "(4x,'CELL-BY-CELL FLOW INFORMATION WILL BE SAVED TO BINARY FILE &
849  &WHENEVER ICBCFL IS NOT ZERO.')"
850  character(len=*), parameter :: fmtflow = &
851  &"(4x, 'FLOWS WILL BE SAVED TO FILE: ', a, /4x, 'OPENED ON UNIT: ', I7)"
852  character(len=*), parameter :: fmtorigss = &
853  "(4X,'ORIGINAL_SPECIFIC_STORAGE OPTION:',/, &
854  &1X,'The original specific storage formulation will be used')"
855  character(len=*), parameter :: fmtstoc = &
856  "(4X,'STORAGECOEFFICIENT OPTION:',/, &
857  &1X,'Read storage coefficient rather than specific storage')"
858  character(len=*), parameter :: fmtconfss = &
859  "(4X,'SS_CONFINED_ONLY OPTION:',/, &
860  &1X,'Specific storage changes only occur under confined conditions')"
861  !
862  write (this%iout, '(1x,a)') 'PROCESSING STORAGE OPTIONS'
863  !
864  if (found%ipakcb) then
865  write (this%iout, fmtisvflow)
866  end if
867  !
868  if (found%istor_coef) then
869  write (this%iout, fmtstoc)
870  end if
871  !
872  if (found%ss_confined_only) then
873  write (this%iout, fmtconfss)
874  end if
875  !
876  if (found%iorig_ss) then
877  write (this%iout, fmtorigss)
878  end if
879  !
880  if (found%iconf_ss) then
881  write (this%iout, fmtconfss)
882  end if
883  !
884  write (this%iout, '(1x,a)') 'END OF STORAGE OPTIONS'

◆ save_old_ss_sy()

subroutine gwfstomodule::save_old_ss_sy ( class(gwfstotype this)

Save ss and sy values from the previous time step for use with storage integration when integratechanges is non-zero.

Parameters
thisGwfStoType object

Definition at line 988 of file gwf-sto.f90.

989  ! -- modules
991  ! -- dummy variables
992  class(GwfStoType) :: this !< GwfStoType object
993  ! -- local variables
994  integer(I4B) :: n
995  !
996  ! -- Allocate TVS arrays if needed
997  if (.not. associated(this%oldss)) then
998  call mem_allocate(this%oldss, this%dis%nodes, 'OLDSS', this%memoryPath)
999  end if
1000  if (this%iusesy == 1 .and. .not. associated(this%oldsy)) then
1001  call mem_allocate(this%oldsy, this%dis%nodes, 'OLDSY', this%memoryPath)
1002  end if
1003  !
1004  ! -- Save current specific storage
1005  do n = 1, this%dis%nodes
1006  this%oldss(n) = this%ss(n)
1007  end do
1008  !
1009  ! -- Save current specific yield, if used
1010  if (this%iusesy == 1) then
1011  do n = 1, this%dis%nodes
1012  this%oldsy(n) = this%sy(n)
1013  end do
1014  end if

◆ source_data()

subroutine gwfstomodule::source_data ( class(gwfstotype this)

Source griddata block parameters for STO package.

Parameters
thisGwfStoType object

Definition at line 892 of file gwf-sto.f90.

893  ! -- modules
896  ! -- dummy variables
897  class(GwfStotype) :: this !< GwfStoType object
898  ! -- local variables
899  character(len=LINELENGTH) :: cellstr
900  logical(LGP) :: isconv
901  character(len=24), dimension(4) :: aname
902  integer(I4B) :: n
903  integer(I4B), dimension(:), pointer, contiguous :: map
904  type(GwfStoParamFoundType) :: found
905  !
906  ! -- initialize data
907  data aname(1)/' ICONVERT'/
908  data aname(2)/' SPECIFIC STORAGE'/
909  data aname(3)/' SPECIFIC YIELD'/
910  !
911  ! -- set map to reduce data input arrays
912  map => null()
913  if (this%dis%nodes < this%dis%nodesuser) map => this%dis%nodeuser
914  !
915  ! -- copy data from input context
916  call mem_set_value(this%iconvert, 'ICONVERT', this%input_mempath, map, &
917  found%iconvert)
918  call mem_set_value(this%ss, 'SS', this%input_mempath, map, found%ss)
919  call mem_set_value(this%sy, 'SY', this%input_mempath, map, found%sy)
920  !
921  ! -- Check for ICONVERT
922  if (found%iconvert) then
923  isconv = .false.
924  do n = 1, this%dis%nodes
925  if (this%iconvert(n) /= 0) then
926  isconv = .true.
927  this%iusesy = 1
928  exit
929  end if
930  end do
931  else
932  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
933  trim(adjustl(aname(1))), ' not found.'
934  call store_error(errmsg)
935  end if
936  !
937  ! -- Check for SS
938  if (found%ss) then
939  ! no-op
940  else
941  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
942  trim(adjustl(aname(2))), ' not found.'
943  call store_error(errmsg)
944  end if
945  !
946  ! -- Check for SY
947  if (found%sy) then
948  ! no-op
949  else
950  if (isconv) then
951  write (errmsg, '(a, a, a)') 'Error in GRIDDATA block: ', &
952  trim(adjustl(aname(3))), ' not found.'
953  call store_error(errmsg)
954  end if
955  end if
956  !
957  if (count_errors() > 0) then
958  call store_error_filename(this%input_fname)
959  end if
960  !
961  ! -- Check SS and SY for negative values
962  do n = 1, this%dis%nodes
963  if (this%ss(n) < dzero) then
964  call this%dis%noder_to_string(n, cellstr)
965  write (errmsg, '(a,2(1x,a),1x,g0,1x,a)') &
966  'Error in SS DATA: SS value in cell', trim(adjustl(cellstr)), &
967  'is less than zero (', this%ss(n), ').'
968  call store_error(errmsg)
969  end if
970  if (found%sy) then
971  if (this%sy(n) < dzero) then
972  call this%dis%noder_to_string(n, cellstr)
973  write (errmsg, '(a,2(1x,a),1x,g0,1x,a)') &
974  'Error in SY DATA: SY value in cell', trim(adjustl(cellstr)), &
975  'is less than zero (', this%sy(n), ').'
976  call store_error(errmsg)
977  end if
978  end if
979  end do
Here is the call graph for this function:

◆ source_options()

subroutine gwfstomodule::source_options ( class(gwfstotype this)

Source options block parameters for STO package.

Parameters
thisGwfStoType object

Definition at line 776 of file gwf-sto.f90.

777  ! -- modules
778  use constantsmodule, only: lenmempath
784  ! -- dummy variables
785  class(GwfStoType) :: this !< GwfStoType object
786  ! -- local variables
787  type(GwfStoParamFoundType) :: found
788  type(CharacterStringType), dimension(:), pointer, contiguous :: tvs6_mempaths
789  character(len=LINELENGTH) :: tvs6_filename
790  character(len=LENMEMPATH) :: tvs6_mempath
791  !
792  ! -- source package input
793  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
794  call mem_set_value(this%istor_coef, 'ISTOR_COEF', this%input_mempath, &
795  found%istor_coef)
796  call mem_set_value(this%iconf_ss, 'SS_CONFINED_ONLY', this%input_mempath, &
797  found%ss_confined_only)
798  call mem_set_value(this%iorig_ss, 'IORIG_SS', this%input_mempath, &
799  found%iorig_ss)
800  call mem_set_value(this%iconf_ss, 'ICONF_SS', this%input_mempath, &
801  found%iconf_ss)
802  !
803  if (found%ipakcb) then
804  this%ipakcb = -1
805  end if
806  !
807  if (found%ss_confined_only) then
808  this%iorig_ss = 0
809  end if
810  !
811  ! -- TVS6 subpackage
812  if (filein_fname(tvs6_filename, 'TVS6_FILENAME', &
813  this%input_mempath, this%input_fname)) then
814  call mem_setptr(tvs6_mempaths, 'TVS6_MEMPATH', this%input_mempath)
815  tvs6_mempath = tvs6_mempaths(1)
816  this%intvs = 1 ! tvs active
817  call tvs_cr(this%tvs, this%name_model, tvs6_mempath, this%intvs, this%iout)
818  end if
819  !
820  if (found%iconf_ss) then
821  this%iorig_ss = 0
822  end if
823  !
824  ! -- log found options
825  call this%log_options(found)
826  !
827  ! -- set omega value used for saturation calculations
828  if (this%inewton > 0) then
829  this%satomega = dem6
830  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ sto_ad()

subroutine gwfstomodule::sto_ad ( class(gwfstotype this)

Advance data in the STO package.

Parameters
thisGwfStoType object

Definition at line 203 of file gwf-sto.f90.

204  ! -- modules
205  use tdismodule, only: kstp
206  ! -- dummy variables
207  class(GwfStoType) :: this !< GwfStoType object
208  !
209  ! -- Store ss and sy values from end of last time step if needed
210  if (this%integratechanges /= 0 .and. kstp > 1) then
211  call this%save_old_ss_sy()
212  end if
213  !
214  ! -- TVS
215  if (this%intvs /= 0) then
216  call this%tvs%ad()
217  end if
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24

◆ sto_ar()

subroutine gwfstomodule::sto_ar ( class(gwfstotype this,
class(disbasetype), intent(in), pointer  dis,
integer(i4b), dimension(:), pointer, contiguous  ibound 
)

Method to allocate and read static data for the STO package.

Parameters
thisGwfStoType object
[in]dismodel discretization object
iboundmodel ibound array

Definition at line 103 of file gwf-sto.f90.

104  ! -- modules
107  ! -- dummy variables
108  class(GwfStoType) :: this !< GwfStoType object
109  class(DisBaseType), pointer, intent(in) :: dis !< model discretization object
110  integer(I4B), dimension(:), pointer, contiguous :: ibound !< model ibound array
111  ! -- local variables
112  ! -- formats
113  character(len=*), parameter :: fmtsto = &
114  "(1x,/1x,'STO -- STORAGE PACKAGE, VERSION 1, 5/19/2014', &
115  &' INPUT READ FROM MEMPATH: ', A, //)"
116  !
117  ! --print a message identifying the storage package.
118  write (this%iout, fmtsto) this%input_mempath
119  !
120  ! -- store pointers to arguments that were passed in
121  this%dis => dis
122  this%ibound => ibound
123  !
124  ! -- set pointer to gwf iss
125  call mem_setptr(this%iss, 'ISS', create_mem_path(this%name_model))
126  !
127  ! -- Allocate arrays
128  call this%allocate_arrays(dis%nodes)
129  !!
130  !! -- Register side effect handlers
131  !call this%register_handlers()
132  !
133  ! -- Read storage options
134  call this%source_options()
135  !
136  ! -- read the data block
137  call this%source_data()
138  !
139  ! -- TVS
140  if (this%intvs /= 0) then
141  call this%tvs%ar(this%dis)
142  end if
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:

◆ sto_bd()

subroutine gwfstomodule::sto_bd ( class(gwfstotype this,
integer(i4b), intent(in)  isuppress_output,
type(budgettype), intent(inout)  model_budget 
)

Budget calculation for the STO package components. Components include specific storage and specific yield storage.

Parameters
thisGwfStoType object
[in]isuppress_outputflag to suppress model output
[in,out]model_budgetmodel budget object

Definition at line 571 of file gwf-sto.f90.

572  ! -- modules
573  use tdismodule, only: delt
575  ! -- dummy variables
576  class(GwfStoType) :: this !< GwfStoType object
577  integer(I4B), intent(in) :: isuppress_output !< flag to suppress model output
578  type(BudgetType), intent(inout) :: model_budget !< model budget object
579  ! -- local variables
580  real(DP) :: rin
581  real(DP) :: rout
582  !
583  ! -- Add confined storage rates to model budget
584  call rate_accumulator(this%strgss, rin, rout)
585  call model_budget%addentry(rin, rout, delt, budtxt(1), &
586  isuppress_output, ' STORAGE')
587  !
588  ! -- Add unconfined storage rates to model budget
589  if (this%iusesy == 1) then
590  call rate_accumulator(this%strgsy, rin, rout)
591  call model_budget%addentry(rin, rout, delt, budtxt(2), &
592  isuppress_output, ' STORAGE')
593  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:

◆ sto_cq()

subroutine gwfstomodule::sto_cq ( class(gwfstotype this,
real(dp), dimension(:), intent(inout), contiguous  flowja,
real(dp), dimension(:), intent(in), contiguous  hnew,
real(dp), dimension(:), intent(in), contiguous  hold 
)

Flow calculation for the STO package components. Components include specific storage and specific yield storage.

Parameters
thisGwfStoType object
[in,out]flowjaconnection flows
[in]hnewcurrent head
[in]holdprevious head

Definition at line 446 of file gwf-sto.f90.

447  ! -- modules
448  use tdismodule, only: delt
449  ! -- dummy variables
450  class(GwfStoType) :: this !< GwfStoType object
451  real(DP), dimension(:), contiguous, intent(inout) :: flowja !< connection flows
452  real(DP), dimension(:), contiguous, intent(in) :: hnew !< current head
453  real(DP), dimension(:), contiguous, intent(in) :: hold !< previous head
454  ! -- local variables
455  integer(I4B) :: n
456  integer(I4B) :: idiag
457  real(DP) :: rate
458  real(DP) :: tled
459  real(DP) :: sc1
460  real(DP) :: sc2
461  real(DP) :: rho1
462  real(DP) :: rho2
463  real(DP) :: sc1old
464  real(DP) :: sc2old
465  real(DP) :: rho1old
466  real(DP) :: rho2old
467  real(DP) :: tp
468  real(DP) :: bt
469  real(DP) :: snold
470  real(DP) :: snnew
471  real(DP) :: aterm
472  real(DP) :: rhsterm
473  !
474  ! -- initialize strg arrays
475  do n = 1, this%dis%nodes
476  this%strgss(n) = dzero
477  this%strgsy(n) = dzero
478  end do
479  !
480  ! -- Set strt to zero or calculate terms if not steady-state stress period
481  if (this%iss == 0) then
482  !
483  ! -- set variables
484  tled = done / delt
485  !
486  ! -- Calculate storage change
487  do n = 1, this%dis%nodes
488  if (this%ibound(n) <= 0) cycle
489  ! -- aquifer elevations and thickness
490  tp = this%dis%top(n)
491  bt = this%dis%bot(n)
492  !
493  ! -- aquifer saturation
494  if (this%iconvert(n) == 0) then
495  snold = done
496  snnew = done
497  else
498  snold = squadraticsaturation(tp, bt, hold(n), this%satomega)
499  snnew = squadraticsaturation(tp, bt, hnew(n), this%satomega)
500  end if
501  !
502  ! -- primary storage coefficient
503  sc1 = sscapacity(this%istor_coef, tp, bt, this%dis%area(n), this%ss(n))
504  rho1 = sc1 * tled
505  !
506  if (this%integratechanges /= 0) then
507  ! -- Integration of storage changes (e.g. when using TVS):
508  ! separate the old (start of time step) and new (end of time step)
509  ! primary storage capacities
510  sc1old = sscapacity(this%istor_coef, tp, bt, this%dis%area(n), &
511  this%oldss(n))
512  rho1old = sc1old * tled
513  else
514  ! -- No integration of storage changes: old and new values are
515  ! identical => normal MF6 storage formulation
516  rho1old = rho1
517  end if
518  !
519  ! -- calculate specific storage terms and rate
520  call ssterms(this%iconvert(n), this%iorig_ss, this%iconf_ss, tp, bt, &
521  rho1, rho1old, snnew, snold, hnew(n), hold(n), &
522  aterm, rhsterm, rate)
523  !
524  ! -- save rate
525  this%strgss(n) = rate
526  !
527  ! -- add storage term to flowja
528  idiag = this%dis%con%ia(n)
529  flowja(idiag) = flowja(idiag) + rate
530  !
531  ! -- specific yield
532  rate = dzero
533  if (this%iconvert(n) /= 0) then
534  !
535  ! -- secondary storage coefficient
536  sc2 = sycapacity(this%dis%area(n), this%sy(n))
537  rho2 = sc2 * tled
538  !
539  if (this%integratechanges /= 0) then
540  ! -- Integration of storage changes (e.g. when using TVS):
541  ! separate the old (start of time step) and new (end of time
542  ! step) secondary storage capacities
543  sc2old = sycapacity(this%dis%area(n), this%oldsy(n))
544  rho2old = sc2old * tled
545  else
546  ! -- No integration of storage changes: old and new values are
547  ! identical => normal MF6 storage formulation
548  rho2old = rho2
549  end if
550  !
551  ! -- calculate specific yield storage terms and rate
552  call syterms(tp, bt, rho2, rho2old, snnew, snold, &
553  aterm, rhsterm, rate)
554 
555  end if
556  this%strgsy(n) = rate
557  !
558  ! -- add storage term to flowja
559  idiag = this%dis%con%ia(n)
560  flowja(idiag) = flowja(idiag) + rate
561  end do
562  end if
Here is the call graph for this function:

◆ sto_cr()

subroutine, public gwfstomodule::sto_cr ( type(gwfstotype), pointer  stoobj,
character(len=*), intent(in)  name_model,
character(len=*), intent(in)  mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

Create a new storage (STO) object

Parameters
stoobjGwfStoType object
[in]name_modelname of model
[in]mempathinput context mem path
[in]inunitpackage input file unit
[in]ioutmodel listing file unit

Definition at line 76 of file gwf-sto.f90.

77  ! -- dummy variables
78  type(GwfStoType), pointer :: stoobj !< GwfStoType object
79  character(len=*), intent(in) :: name_model !< name of model
80  character(len=*), intent(in) :: mempath !< input context mem path
81  integer(I4B), intent(in) :: inunit !< package input file unit
82  integer(I4B), intent(in) :: iout !< model listing file unit
83  !
84  ! -- Create the object
85  allocate (stoobj)
86  !
87  ! -- create name and memory path
88  call stoobj%set_names(1, name_model, 'STO', 'STO', mempath)
89  !
90  ! -- Allocate scalars
91  call stoobj%allocate_scalars()
92  !
93  ! -- Set variables
94  stoobj%inunit = inunit
95  stoobj%iout = iout
Here is the caller graph for this function:

◆ sto_da()

subroutine gwfstomodule::sto_da ( class(gwfstotype this)

Deallocate STO package scalars and arrays.

Parameters
thisGwfStoType object

Definition at line 646 of file gwf-sto.f90.

647  ! -- modules
649  ! -- dummy variables
650  class(GwfStoType) :: this !< GwfStoType object
651  !
652  ! -- TVS
653  if (this%intvs /= 0) then
654  call this%tvs%da()
655  deallocate (this%tvs)
656  end if
657  !
658  ! -- Deallocate arrays if package is active
659  if (this%inunit > 0) then
660  call mem_deallocate(this%iconvert)
661  call mem_deallocate(this%ss)
662  call mem_deallocate(this%sy)
663  call mem_deallocate(this%strgss)
664  call mem_deallocate(this%strgsy)
665  !
666  ! -- deallocate TVS arrays
667  if (associated(this%oldss)) then
668  call mem_deallocate(this%oldss)
669  end if
670  if (associated(this%oldsy)) then
671  call mem_deallocate(this%oldsy)
672  end if
673  !
674  ! -- nullify input context pointers
675  nullify (this%iper)
676  nullify (this%storage)
677  end if
678  !
679  ! -- Deallocate scalars
680  call mem_deallocate(this%istor_coef)
681  call mem_deallocate(this%iconf_ss)
682  call mem_deallocate(this%iorig_ss)
683  call mem_deallocate(this%iusesy)
684  call mem_deallocate(this%satomega)
685  call mem_deallocate(this%integratechanges)
686  call mem_deallocate(this%intvs)
687  !
688  ! -- deallocate parent
689  call this%NumericalPackageType%da()

◆ sto_fc()

subroutine gwfstomodule::sto_fc ( class(gwfstotype this,
integer(i4b), intent(in)  kiter,
real(dp), dimension(:), intent(in)  hold,
real(dp), dimension(:), intent(in)  hnew,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs 
)

Fill the coefficient matrix and right-hand side with the STO package terms.

Parameters
thisGwfStoType object
[in]kiterouter iteration number
[in]holdprevious heads
[in]hnewcurrent heads
matrix_slnA matrix
[in]idxgloglobal index model to solution
[in,out]rhsright-hand side

Definition at line 225 of file gwf-sto.f90.

226  ! -- modules
227  use tdismodule, only: delt
228  ! -- dummy variables
229  class(GwfStoType) :: this !< GwfStoType object
230  integer(I4B), intent(in) :: kiter !< outer iteration number
231  real(DP), intent(in), dimension(:) :: hold !< previous heads
232  real(DP), intent(in), dimension(:) :: hnew !< current heads
233  class(MatrixBaseType), pointer :: matrix_sln !< A matrix
234  integer(I4B), intent(in), dimension(:) :: idxglo !< global index model to solution
235  real(DP), intent(inout), dimension(:) :: rhs !< right-hand side
236  ! -- local variables
237  integer(I4B) :: n
238  integer(I4B) :: idiag
239  real(DP) :: tled
240  real(DP) :: sc1
241  real(DP) :: sc2
242  real(DP) :: rho1
243  real(DP) :: rho2
244  real(DP) :: sc1old
245  real(DP) :: sc2old
246  real(DP) :: rho1old
247  real(DP) :: rho2old
248  real(DP) :: tp
249  real(DP) :: bt
250  real(DP) :: snold
251  real(DP) :: snnew
252  real(DP) :: aterm
253  real(DP) :: rhsterm
254  ! -- formats
255  character(len=*), parameter :: fmtsperror = &
256  &"('DETECTED TIME STEP LENGTH OF ZERO. GWF STORAGE PACKAGE CANNOT BE ', &
257  &'USED UNLESS DELT IS NON-ZERO.')"
258  !
259  ! -- test if steady-state stress period
260  if (this%iss /= 0) return
261  !
262  ! -- Ensure time step length is not zero
263  if (delt == dzero) then
264  write (errmsg, fmtsperror)
265  call store_error(errmsg, terminate=.true.)
266  end if
267  !
268  ! -- set variables
269  tled = done / delt
270  !
271  ! -- loop through and calculate storage contribution to hcof and rhs
272  do n = 1, this%dis%nodes
273  idiag = this%dis%con%ia(n)
274  if (this%ibound(n) < 1) cycle
275  !
276  ! -- aquifer elevations and thickness
277  tp = this%dis%top(n)
278  bt = this%dis%bot(n)
279  !
280  ! -- aquifer saturation
281  if (this%iconvert(n) == 0) then
282  snold = done
283  snnew = done
284  else
285  snold = squadraticsaturation(tp, bt, hold(n), this%satomega)
286  snnew = squadraticsaturation(tp, bt, hnew(n), this%satomega)
287  end if
288  !
289  ! -- storage coefficients
290  sc1 = sscapacity(this%istor_coef, tp, bt, this%dis%area(n), this%ss(n))
291  rho1 = sc1 * tled
292  !
293  if (this%integratechanges /= 0) then
294  ! -- Integration of storage changes (e.g. when using TVS):
295  ! separate the old (start of time step) and new (end of time step)
296  ! primary storage capacities
297  sc1old = sscapacity(this%istor_coef, tp, bt, this%dis%area(n), &
298  this%oldss(n))
299  rho1old = sc1old * tled
300  else
301  ! -- No integration of storage changes: old and new values are
302  ! identical => normal MF6 storage formulation
303  rho1old = rho1
304  end if
305  !
306  ! -- calculate specific storage terms
307  call ssterms(this%iconvert(n), this%iorig_ss, this%iconf_ss, tp, bt, &
308  rho1, rho1old, snnew, snold, hnew(n), hold(n), &
309  aterm, rhsterm)
310  !
311  ! -- add specific storage terms to amat and rhs
312  call matrix_sln%add_value_pos(idxglo(idiag), aterm)
313  rhs(n) = rhs(n) + rhsterm
314  !
315  ! -- specific yield
316  if (this%iconvert(n) /= 0) then
317  rhsterm = dzero
318  !
319  ! -- secondary storage coefficient
320  sc2 = sycapacity(this%dis%area(n), this%sy(n))
321  rho2 = sc2 * tled
322  !
323  if (this%integratechanges /= 0) then
324  ! -- Integration of storage changes (e.g. when using TVS):
325  ! separate the old (start of time step) and new (end of time step)
326  ! secondary storage capacities
327  sc2old = sycapacity(this%dis%area(n), this%oldsy(n))
328  rho2old = sc2old * tled
329  else
330  ! -- No integration of storage changes: old and new values are
331  ! identical => normal MF6 storage formulation
332  rho2old = rho2
333  end if
334  !
335  ! -- calculate specific storage terms
336  call syterms(tp, bt, rho2, rho2old, snnew, snold, &
337  aterm, rhsterm)
338 !
339  ! -- add specific yield terms to amat and rhs
340  call matrix_sln%add_value_pos(idxglo(idiag), aterm)
341  rhs(n) = rhs(n) + rhsterm
342  end if
343  end do
Here is the call graph for this function:

◆ sto_fn()

subroutine gwfstomodule::sto_fn ( class(gwfstotype this,
integer(i4b), intent(in)  kiter,
real(dp), dimension(:), intent(in)  hold,
real(dp), dimension(:), intent(in)  hnew,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), dimension(:), intent(in)  idxglo,
real(dp), dimension(:), intent(inout)  rhs 
)

Fill the coefficient matrix and right-hand side with STO package with Newton-Raphson terms.

Parameters
thisGwfStoType object
[in]kiterouter iteration number
[in]holdprevious heads
[in]hnewcurrent heads
matrix_slnA matrix
[in]idxgloglobal index model to solution
[in,out]rhsright-hand side

Definition at line 352 of file gwf-sto.f90.

353  ! -- modules
354  use tdismodule, only: delt
355  ! -- dummy variables
356  class(GwfStoType) :: this !< GwfStoType object
357  integer(I4B), intent(in) :: kiter !< outer iteration number
358  real(DP), intent(in), dimension(:) :: hold !< previous heads
359  real(DP), intent(in), dimension(:) :: hnew !< current heads
360  class(MatrixBaseType), pointer :: matrix_sln !< A matrix
361  integer(I4B), intent(in), dimension(:) :: idxglo !< global index model to solution
362  real(DP), intent(inout), dimension(:) :: rhs !< right-hand side
363  ! -- local variables
364  integer(I4B) :: n
365  integer(I4B) :: idiag
366  real(DP) :: tled
367  real(DP) :: sc1
368  real(DP) :: sc2
369  real(DP) :: rho1
370  real(DP) :: rho2
371  real(DP) :: tp
372  real(DP) :: bt
373  real(DP) :: tthk
374  real(DP) :: h
375  real(DP) :: snnew
376  real(DP) :: derv
377  real(DP) :: rterm
378  real(DP) :: drterm
379  !
380  ! -- test if steady-state stress period
381  if (this%iss /= 0) return
382  !
383  ! -- set variables
384  tled = done / delt
385  !
386  ! -- loop through and calculate storage contribution to hcof and rhs
387  do n = 1, this%dis%nodes
388  idiag = this%dis%con%ia(n)
389  if (this%ibound(n) <= 0) cycle
390  !
391  ! -- aquifer elevations and thickness
392  tp = this%dis%top(n)
393  bt = this%dis%bot(n)
394  tthk = tp - bt
395  h = hnew(n)
396  !
397  ! -- aquifer saturation
398  snnew = squadraticsaturation(tp, bt, h)
399  !
400  ! -- storage coefficients
401  sc1 = sscapacity(this%istor_coef, tp, bt, this%dis%area(n), this%ss(n))
402  sc2 = sycapacity(this%dis%area(n), this%sy(n))
403  rho1 = sc1 * tled
404  rho2 = sc2 * tled
405  !
406  ! -- calculate newton terms for specific storage
407  ! and specific yield
408  if (this%iconvert(n) /= 0) then
409  !
410  ! -- calculate saturation derivative
411  derv = squadraticsaturationderivative(tp, bt, h)
412  !
413  ! -- newton terms for specific storage
414  if (this%iconf_ss == 0) then
415  if (this%iorig_ss == 0) then
416  drterm = -rho1 * derv * (h - bt) + rho1 * tthk * snnew * derv
417  else
418  drterm = -(rho1 * derv * h)
419  end if
420  call matrix_sln%add_value_pos(idxglo(idiag), drterm)
421  rhs(n) = rhs(n) + drterm * h
422  end if
423  !
424  ! -- newton terms for specific yield
425  ! only calculated if the current saturation
426  ! is less than one
427  if (snnew < done) then
428  ! -- calculate newton terms for specific yield
429  if (snnew > dzero) then
430  rterm = -rho2 * tthk * snnew
431  drterm = -rho2 * tthk * derv
432  call matrix_sln%add_value_pos(idxglo(idiag), drterm + rho2)
433  rhs(n) = rhs(n) - rterm + drterm * h + rho2 * bt
434  end if
435  end if
436  end if
437  end do
Here is the call graph for this function:

◆ sto_rp()

subroutine gwfstomodule::sto_rp ( class(gwfstotype this)

Method to read and prepare stress period data for the STO package.

Parameters
thisGwfStoType object

Definition at line 150 of file gwf-sto.f90.

151  ! -- modules
152  use tdismodule, only: kper
153  implicit none
154  ! -- dummy variables
155  class(GwfStoType) :: this !< GwfStoType object
156  ! -- local variables
157  character(len=16) :: css(0:1)
158  ! -- data
159  data css(0)/' TRANSIENT'/
160  data css(1)/' STEADY-STATE'/
161  !
162  ! -- Store ss and sy values from end of last stress period if needed
163  if (this%integratechanges /= 0) then
164  call this%save_old_ss_sy()
165  end if
166  !
167  ! -- confirm package is active
168  if (this%inunit <= 0) return
169  !
170  ! -- confirm loaded iper
171  if (this%iper /= kper) return
172  !
173  write (this%iout, '(//,1x,a)') 'PROCESSING STORAGE PERIOD DATA'
174  !
175  ! -- set period iss
176  if (this%storage == 'STEADY-STATE') then
177  this%iss = 1
178  else if (this%storage == 'TRANSIENT') then
179  this%iss = 0
180  else
181  write (errmsg, '(a,a)') 'Unknown STORAGE data tag: ', &
182  trim(this%storage)
183  call store_error(errmsg)
184  call store_error_filename(this%input_fname)
185  end if
186  !
187  write (this%iout, '(1x,a)') 'END PROCESSING STORAGE PERIOD DATA'
188  !
189  write (this%iout, '(//1X,A,I0,A,A,/)') &
190  'STRESS PERIOD ', kper, ' IS ', trim(adjustl(css(this%iss)))
191  !
192  ! -- TVS
193  if (this%intvs /= 0) then
194  call this%tvs%rp()
195  end if
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ sto_save_model_flows()

subroutine gwfstomodule::sto_save_model_flows ( class(gwfstotype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  icbcun 
)

Save cell-by-cell budget terms for the STO package.

Parameters
thisGwfStoType object
[in]icbcflflag to output budget data
[in]icbcuncell-by-cell file unit number

Definition at line 601 of file gwf-sto.f90.

602  ! -- dummy variables
603  class(GwfStoType) :: this !< GwfStoType object
604  integer(I4B), intent(in) :: icbcfl !< flag to output budget data
605  integer(I4B), intent(in) :: icbcun !< cell-by-cell file unit number
606  ! -- local variables
607  integer(I4B) :: ibinun
608  integer(I4B) :: iprint, nvaluesp, nwidthp
609  character(len=1) :: cdatafmp = ' ', editdesc = ' '
610  real(DP) :: dinact
611  !
612  ! -- Set unit number for binary output
613  if (this%ipakcb < 0) then
614  ibinun = icbcun
615  elseif (this%ipakcb == 0) then
616  ibinun = 0
617  else
618  ibinun = this%ipakcb
619  end if
620  if (icbcfl == 0) ibinun = 0
621  !
622  ! -- Record the storage rates if requested
623  if (ibinun /= 0) then
624  iprint = 0
625  dinact = dzero
626  !
627  ! -- storage(ss)
628  call this%dis%record_array(this%strgss, this%iout, iprint, -ibinun, &
629  budtxt(1), cdatafmp, nvaluesp, &
630  nwidthp, editdesc, dinact)
631  !
632  ! -- storage(sy)
633  if (this%iusesy == 1) then
634  call this%dis%record_array(this%strgsy, this%iout, iprint, -ibinun, &
635  budtxt(2), cdatafmp, nvaluesp, &
636  nwidthp, editdesc, dinact)
637  end if
638  end if

Variable Documentation

◆ budtxt

character(len=lenbudtxt), dimension(2) gwfstomodule::budtxt = [' STO-SS', ' STO-SY']

Definition at line 27 of file gwf-sto.f90.

27  character(len=LENBUDTXT), dimension(2) :: budtxt = & !< text labels for budget terms
28  &[' STO-SS', ' STO-SY']