MODFLOW 6  version 6.6.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 733 of file gwf-sto.f90.

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

700  ! -- modules
702  ! -- dummy variables
703  class(GwfStoType) :: this !< GwfStoType object
704  !
705  ! -- allocate scalars in NumericalPackageType
706  call this%NumericalPackageType%allocate_scalars()
707  !
708  ! -- allocate scalars
709  call mem_allocate(this%istor_coef, 'ISTOR_COEF', this%memoryPath)
710  call mem_allocate(this%iconf_ss, 'ICONF_SS', this%memoryPath)
711  call mem_allocate(this%iorig_ss, 'IORIG_SS', this%memoryPath)
712  call mem_allocate(this%iusesy, 'IUSESY', this%memoryPath)
713  call mem_allocate(this%satomega, 'SATOMEGA', this%memoryPath)
714  call mem_allocate(this%integratechanges, 'INTEGRATECHANGES', &
715  this%memoryPath)
716  call mem_allocate(this%intvs, 'INTVS', this%memoryPath)
717  !
718  ! -- initialize scalars
719  this%istor_coef = 0
720  this%iconf_ss = 0
721  this%iorig_ss = 0
722  this%iusesy = 0
723  this%satomega = dzero
724  this%integratechanges = 0
725  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 837 of file gwf-sto.f90.

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

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

◆ source_data()

subroutine gwfstomodule::source_data ( class(gwfstotype this)

Source griddata block parameters for STO package.

Parameters
thisGwfStoType object

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

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

778  ! -- modules
779  use constantsmodule, only: lenmempath
783  ! -- dummy variables
784  class(GwfStoType) :: this !< GwfStoType object
785  ! -- local variables
786  type(GwfStoParamFoundType) :: found
787  character(len=LENMEMPATH) :: tvs6_mempath !< mempath of loaded subpackage
788  character(len=LINELENGTH) :: fname
789  !
790  ! -- source package input
791  call mem_set_value(this%ipakcb, 'IPAKCB', this%input_mempath, found%ipakcb)
792  call mem_set_value(this%istor_coef, 'ISTOR_COEF', this%input_mempath, &
793  found%istor_coef)
794  call mem_set_value(this%iconf_ss, 'SS_CONFINED_ONLY', this%input_mempath, &
795  found%ss_confined_only)
796  call mem_set_value(tvs6_mempath, 'TVS6_MEMPATH', this%input_mempath, &
797  found%tvs6_filename)
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  ! -- enforce 0 or 1 TVS6_FILENAME entries in option block
812  if (filein_fname(fname, 'TVS6_FILENAME', this%input_mempath, &
813  this%input_fname)) then
814  this%intvs = getunit()
815  call openfile(this%intvs, this%iout, fname, 'TVS')
816  call tvs_cr(this%tvs, this%name_model, this%intvs, this%iout)
817  end if
818  !
819  if (found%iconf_ss) then
820  this%iorig_ss = 0
821  end if
822  !
823  ! -- log found options
824  call this%log_options(found)
825  !
826  ! -- set omega value used for saturation calculations
827  if (this%inewton > 0) then
828  this%satomega = dem6
829  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
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
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 204 of file gwf-sto.f90.

205  ! -- modules
206  use tdismodule, only: kstp
207  ! -- dummy variables
208  class(GwfStoType) :: this !< GwfStoType object
209  !
210  ! -- Store ss and sy values from end of last time step if needed
211  if (this%integratechanges /= 0 .and. kstp > 1) then
212  call this%save_old_ss_sy()
213  end if
214  !
215  ! -- TVS
216  if (this%intvs /= 0) then
217  call this%tvs%ad()
218  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 104 of file gwf-sto.f90.

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

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

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

78  ! -- dummy variables
79  type(GwfStoType), pointer :: stoobj !< GwfStoType object
80  character(len=*), intent(in) :: name_model !< name of model
81  character(len=*), intent(in) :: mempath !< input context mem path
82  integer(I4B), intent(in) :: inunit !< package input file unit
83  integer(I4B), intent(in) :: iout !< model listing file unit
84  !
85  ! -- Create the object
86  allocate (stoobj)
87  !
88  ! -- create name and memory path
89  call stoobj%set_names(1, name_model, 'STO', 'STO', mempath)
90  !
91  ! -- Allocate scalars
92  call stoobj%allocate_scalars()
93  !
94  ! -- Set variables
95  stoobj%inunit = inunit
96  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 647 of file gwf-sto.f90.

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

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

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

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

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

Variable Documentation

◆ budtxt

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

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

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