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

Data Types

type  gwemodeltype
 

Functions/Subroutines

subroutine, public gwe_cr (filename, id, modelname)
 Create a new groundwater energy transport model object. More...
 
subroutine gwe_df (this)
 Define packages of the GWE model. More...
 
subroutine gwe_ac (this, sparse)
 Add the internal connections of this model to the sparse matrix. More...
 
subroutine gwe_mc (this, matrix_sln)
 Map the positions of the GWE model connections in the numerical solution coefficient matrix. More...
 
subroutine gwe_ar (this)
 GWE Model Allocate and Read. More...
 
subroutine gwe_rp (this)
 GWE Model Read and Prepare. More...
 
subroutine gwe_dt (this)
 GWT Model time step size. More...
 
subroutine gwe_ad (this)
 GWE Model Time Step Advance. More...
 
subroutine gwe_cf (this, kiter)
 GWE Model calculate coefficients. More...
 
subroutine gwe_fc (this, kiter, matrix_sln, inwtflag)
 GWE Model fill coefficients. More...
 
subroutine gwe_cc (this, innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
 GWE Model Final Convergence Check. More...
 
subroutine gwe_cq (this, icnvg, isuppress_output)
 GWE Model calculate flow. More...
 
subroutine gwe_bd (this, icnvg, isuppress_output)
 GWE Model Budget. More...
 
subroutine gwe_ot_flow (this, icbcfl, ibudfl, icbcun)
 GWE model output routine. More...
 
subroutine gwe_da (this)
 Deallocate. More...
 
subroutine gwe_bdentry (this, budterm, budtxt, rowlabel)
 GroundWater Energy Transport Model Budget Entry. More...
 
integer(i4b) function gwe_get_iasym (this)
 return 1 if any package causes the matrix to be asymmetric. Otherwise return 0. More...
 
subroutine allocate_scalars (this, modelname)
 Allocate memory for non-allocatable members. More...
 
subroutine package_create (this, filtyp, ipakid, ipaknum, pakname, mempath, inunit, iout)
 Create boundary condition packages for this model. More...
 
class(gwemodeltype) function, pointer, public castasgwemodel (model)
 Cast to GweModelType. More...
 
subroutine create_bndpkgs (this, bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
 Source package info and begin to process. More...
 
subroutine create_gwe_packages (this, indis)
 Source package info and begin to process. More...
 

Variables

character(len=lenvarname), parameter dvt = 'TEMPERATURE '
 dependent variable type, varies based on model type More...
 
character(len=lenvarname), parameter dvu = 'ENERGY '
 dependent variable unit of measure, either "mass" or "energy" More...
 
character(len=lenvarname), parameter dvua = 'E '
 abbreviation of the dependent variable unit of measure, either "M" or "E" More...
 
integer(i4b), parameter, public gwe_nbasepkg = 50
 GWE base package array descriptors. More...
 
character(len=lenpackagetype), dimension(gwe_nbasepkg), public gwe_basepkg
 
integer(i4b), parameter, public gwe_nmultipkg = 50
 GWE multi package array descriptors. More...
 
character(len=lenpackagetype), dimension(gwe_nmultipkg), public gwe_multipkg
 
integer(i4b), parameter niunit_gwe = GWE_NBASEPKG + GWE_NMULTIPKG
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine gwemodule::allocate_scalars ( class(gwemodeltype this,
character(len=*), intent(in)  modelname 
)
private

A subroutine for allocating the scalars specific to the GWE model type. Additional scalars used by the parent class are allocated by the parent class.

Definition at line 697 of file gwe.f90.

698  ! -- modules
700  ! -- dummy
701  class(GweModelType) :: this
702  character(len=*), intent(in) :: modelname
703  !
704  ! -- Allocate parent class scalars
705  call this%allocate_tsp_scalars(modelname)
706  !
707  ! -- Allocate members that are part of model class
708  call mem_allocate(this%inest, 'INEST', this%memoryPath)
709  call mem_allocate(this%incnd, 'INCND', this%memoryPath)
710  !
711  this%inest = 0
712  this%incnd = 0

◆ castasgwemodel()

class(gwemodeltype) function, pointer, public gwemodule::castasgwemodel ( class(*), pointer  model)
Parameters
modelThe object to be cast
Returns
The GWE model

Definition at line 794 of file gwe.f90.

795  ! -- dummy
796  class(*), pointer :: model !< The object to be cast
797  ! -- return
798  class(GweModelType), pointer :: gwemodel !< The GWE model
799  !
800  gwemodel => null()
801  if (.not. associated(model)) return
802  select type (model)
803  type is (gwemodeltype)
804  gwemodel => model
805  end select
Here is the caller graph for this function:

◆ create_bndpkgs()

subroutine gwemodule::create_bndpkgs ( class(gwemodeltype this,
integer(i4b), dimension(:), intent(inout), allocatable  bndpkgs,
type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  pkgtypes,
type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  pkgnames,
type(characterstringtype), dimension(:), intent(inout), pointer, contiguous  mempaths,
integer(i4b), dimension(:), intent(inout), pointer, contiguous  inunits 
)
private

Definition at line 810 of file gwe.f90.

812  ! -- modules
815  ! -- dummy
816  class(GweModelType) :: this
817  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
818  type(CharacterStringType), dimension(:), contiguous, &
819  pointer, intent(inout) :: pkgtypes
820  type(CharacterStringType), dimension(:), contiguous, &
821  pointer, intent(inout) :: pkgnames
822  type(CharacterStringType), dimension(:), contiguous, &
823  pointer, intent(inout) :: mempaths
824  integer(I4B), dimension(:), contiguous, &
825  pointer, intent(inout) :: inunits
826  ! -- local
827  integer(I4B) :: ipakid, ipaknum
828  character(len=LENFTYPE) :: pkgtype, bndptype
829  character(len=LENPACKAGENAME) :: pkgname
830  character(len=LENMEMPATH) :: mempath
831  integer(I4B), pointer :: inunit
832  integer(I4B) :: n
833  !
834  if (allocated(bndpkgs)) then
835  !
836  ! -- Create stress packages
837  ipakid = 1
838  bndptype = ''
839  do n = 1, size(bndpkgs)
840  !
841  pkgtype = pkgtypes(bndpkgs(n))
842  pkgname = pkgnames(bndpkgs(n))
843  mempath = mempaths(bndpkgs(n))
844  inunit => inunits(bndpkgs(n))
845  !
846  if (bndptype /= pkgtype) then
847  ipaknum = 1
848  bndptype = pkgtype
849  end if
850  !
851  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
852  inunit, this%iout)
853  ipakid = ipakid + 1
854  ipaknum = ipaknum + 1
855  end do
856  !
857  ! -- Cleanup
858  deallocate (bndpkgs)
859  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23

◆ create_gwe_packages()

subroutine gwemodule::create_gwe_packages ( class(gwemodeltype this,
integer(i4b), intent(in)  indis 
)

Definition at line 864 of file gwe.f90.

865  ! -- modules
872  use gweestmodule, only: est_cr
873  use gwecndmodule, only: cnd_cr
874  ! -- dummy
875  class(GweModelType) :: this
876  integer(I4B), intent(in) :: indis
877  ! -- local
878  type(CharacterStringType), dimension(:), contiguous, &
879  pointer :: pkgtypes => null()
880  type(CharacterStringType), dimension(:), contiguous, &
881  pointer :: pkgnames => null()
882  type(CharacterStringType), dimension(:), contiguous, &
883  pointer :: mempaths => null()
884  integer(I4B), dimension(:), contiguous, &
885  pointer :: inunits => null()
886  character(len=LENMEMPATH) :: model_mempath
887  character(len=LENFTYPE) :: pkgtype
888  character(len=LENPACKAGENAME) :: pkgname
889  character(len=LENMEMPATH) :: mempath
890  integer(I4B), pointer :: inunit
891  integer(I4B), dimension(:), allocatable :: bndpkgs
892  integer(I4B) :: n
893  character(len=LENMEMPATH) :: mempathcnd = ''
894  character(len=LENMEMPATH) :: mempathest = ''
895  !
896  ! -- Set input memory paths, input/model and input/model/namfile
897  model_mempath = create_mem_path(component=this%name, context=idm_context)
898  !
899  ! -- Set pointers to model path package info
900  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
901  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
902  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
903  call mem_setptr(inunits, 'INUNITS', model_mempath)
904  !
905  do n = 1, size(pkgtypes)
906  !
907  ! -- Attributes for this input package
908  pkgtype = pkgtypes(n)
909  pkgname = pkgnames(n)
910  mempath = mempaths(n)
911  inunit => inunits(n)
912  !
913  ! -- Create dis package as it is a prerequisite for other packages
914  select case (pkgtype)
915  case ('EST6')
916  this%inest = 1
917  mempathest = mempath
918  case ('CND6')
919  this%incnd = 1
920  mempathcnd = mempath
921  case ('CTP6', 'ESL6', 'LKE6', 'SFE6', &
922  'MWE6', 'UZE6', 'API6')
923  call expandarray(bndpkgs)
924  bndpkgs(size(bndpkgs)) = n
925  case default
926  ! TODO
927  end select
928  end do
929  !
930  ! -- Create packages that are tied directly to model
931  call est_cr(this%est, this%name, mempathest, this%inest, this%iout, &
932  this%fmi, this%eqnsclfac, this%gwecommon)
933  call cnd_cr(this%cnd, this%name, mempathcnd, this%incnd, this%iout, &
934  this%fmi, this%eqnsclfac, this%gwecommon)
935  !
936  ! -- Check to make sure that required ftype's have been specified
937  call this%ftype_check(indis, this%inest)
938  !
939  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
subroutine, public cnd_cr(cndobj, name_model, input_mempath, inunit, iout, fmi, eqnsclfac, gwecommon)
Create a new CND object.
Definition: gwe-cnd.f90:86
@ brief Energy Storage and Transfer (EST) Module
Definition: gwe-est.f90:13
subroutine, public est_cr(estobj, name_model, input_mempath, inunit, iout, fmi, eqnsclfac, gwecommon)
@ brief Create a new EST package object
Definition: gwe-est.f90:105
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ gwe_ac()

subroutine gwemodule::gwe_ac ( class(gwemodeltype this,
type(sparsematrix), intent(inout)  sparse 
)

Definition at line 199 of file gwe.f90.

200  ! -- modules
201  use sparsemodule, only: sparsematrix
202  ! -- dummy
203  class(GweModelType) :: this
204  type(sparsematrix), intent(inout) :: sparse
205  ! -- local
206  class(BndType), pointer :: packobj
207  integer(I4B) :: ip
208  !
209  ! -- Add the internal connections of this model to sparse
210  call this%dis%dis_ac(this%moffset, sparse)
211  if (this%incnd > 0) &
212  call this%cnd%cnd_ac(this%moffset, sparse)
213  !
214  ! -- Add any package connections
215  do ip = 1, this%bndlist%Count()
216  packobj => getbndfromlist(this%bndlist, ip)
217  call packobj%bnd_ac(this%moffset, sparse)
218  end do
Here is the call graph for this function:

◆ gwe_ad()

subroutine gwemodule::gwe_ad ( class(gwemodeltype this)

This subroutine calls the attached packages' advance subroutines

Definition at line 348 of file gwe.f90.

349  ! -- modules
351  ! -- dummy
352  class(GweModelType) :: this
353  class(BndType), pointer :: packobj
354  ! -- local
355  integer(I4B) :: irestore
356  integer(I4B) :: ip, n
357  !
358  ! -- Reset state variable
359  irestore = 0
360  if (ifailedstepretry > 0) irestore = 1
361  if (irestore == 0) then
362  !
363  ! -- Copy x into xold
364  do n = 1, this%dis%nodes
365  if (this%ibound(n) == 0) then
366  this%xold(n) = dzero
367  else
368  this%xold(n) = this%x(n)
369  end if
370  end do
371  else
372  !
373  ! -- Copy xold into x if this time step is a redo
374  do n = 1, this%dis%nodes
375  this%x(n) = this%xold(n)
376  end do
377  end if
378  !
379  ! -- Advance fmi
380  call this%fmi%fmi_ad(this%x)
381  !
382  ! -- Advance
383  if (this%incnd > 0) call this%cnd%cnd_ad()
384  if (this%inssm > 0) call this%ssm%ssm_ad()
385  do ip = 1, this%bndlist%Count()
386  packobj => getbndfromlist(this%bndlist, ip)
387  call packobj%bnd_ad()
388  if (isimcheck > 0) then
389  call packobj%bnd_ck()
390  end if
391  end do
392  !
393  ! -- Push simulated values to preceding time/subtime step
394  call this%obs%obs_ad()
integer(i4b) isimcheck
simulation input check flag (1) to check input, (0) to ignore checks
integer(i4b) ifailedstepretry
current retry for this time step
Here is the call graph for this function:

◆ gwe_ar()

subroutine gwemodule::gwe_ar ( class(gwemodeltype this)
private

This subroutine:

  • allocates and reads packages that are part of this model,
  • allocates memory for arrays used by this model object

Definition at line 251 of file gwe.f90.

252  ! -- modules
253  use constantsmodule, only: dhnoflo
254  ! -- dummy
255  class(GweModelType) :: this
256  ! -- locals
257  integer(I4B) :: ip
258  class(BndType), pointer :: packobj
259  !
260  ! -- Allocate and read modules attached to model
261  call this%fmi%fmi_ar(this%ibound)
262  if (this%inmvt > 0) call this%mvt%mvt_ar()
263  if (this%inic > 0) call this%ic%ic_ar(this%x)
264  if (this%inest > 0) call this%est%est_ar(this%dis, this%ibound)
265  if (this%inadv > 0) call this%adv%adv_ar(this%dis, this%ibound)
266  if (this%incnd > 0) call this%cnd%cnd_ar(this%ibound, this%est%porosity)
267  if (this%inssm > 0) call this%ssm%ssm_ar(this%dis, this%ibound, this%x)
268  if (this%inobs > 0) call this%obs%tsp_obs_ar(this%ic, this%x, this%flowja)
269  !
270  ! -- Set governing equation scale factor
271  this%eqnsclfac = this%gwecommon%gwerhow * this%gwecommon%gwecpw
272  !
273  ! -- Call dis_ar to write binary grid file
274  !call this%dis%dis_ar(this%npf%icelltype)
275  !
276  ! -- Set up output control
277  call this%oc%oc_ar(this%x, this%dis, dhnoflo, this%depvartype)
278  call this%budget%set_ibudcsv(this%oc%ibudcsv)
279  !
280  ! -- Package input files now open, so allocate and read
281  do ip = 1, this%bndlist%Count()
282  packobj => getbndfromlist(this%bndlist, ip)
283  call packobj%set_pointers(this%dis%nodes, this%ibound, this%x, &
284  this%xold, this%flowja)
285  ! -- Read and allocate package
286  call packobj%bnd_ar()
287  end do
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
Here is the call graph for this function:

◆ gwe_bd()

subroutine gwemodule::gwe_bd ( class(gwemodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)

This subroutine:

  • calculates intercell flows (flowja)
  • calculates package contributions to the model budget

Definition at line 533 of file gwe.f90.

534  ! -- dummy
535  class(GweModelType) :: this
536  integer(I4B), intent(in) :: icnvg
537  integer(I4B), intent(in) :: isuppress_output
538  ! -- local
539  integer(I4B) :: ip
540  class(BndType), pointer :: packobj
541  !
542  ! -- Save the solution convergence flag
543  this%icnvg = icnvg
544  !
545  ! -- Budget routines (start by resetting). Sole purpose of this section
546  ! is to add in and outs to model budget. All ins and out for a model
547  ! should be added here to this%budget. In a subsequent exchange call,
548  ! exchange flows might also be added.
549  call this%budget%reset()
550  if (this%inest > 0) call this%est%est_bd(isuppress_output, this%budget)
551  if (this%inssm > 0) call this%ssm%ssm_bd(isuppress_output, this%budget)
552  if (this%infmi > 0) call this%fmi%fmi_bd(isuppress_output, this%budget)
553  if (this%inmvt > 0) call this%mvt%mvt_bd(this%x, this%x)
554  do ip = 1, this%bndlist%Count()
555  packobj => getbndfromlist(this%bndlist, ip)
556  call packobj%bnd_bd(this%budget)
557  end do
Here is the call graph for this function:

◆ gwe_bdentry()

subroutine gwemodule::gwe_bdentry ( class(gwemodeltype this,
real(dp), dimension(:, :), intent(in)  budterm,
character(len=lenbudtxt), dimension(:), intent(in)  budtxt,
character(len=*), intent(in)  rowlabel 
)

This subroutine adds a budget entry to the flow budget. It was added as a method for the gwe model object so that the exchange object could add its contributions.

Definition at line 648 of file gwe.f90.

649  ! -- modules
650  use constantsmodule, only: lenbudtxt
651  use tdismodule, only: delt
652  ! -- dummy
653  class(GweModelType) :: this
654  real(DP), dimension(:, :), intent(in) :: budterm
655  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
656  character(len=*), intent(in) :: rowlabel
657  !
658  call this%budget%addentry(budterm, delt, budtxt, rowlabel=rowlabel)
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29

◆ gwe_cc()

subroutine gwemodule::gwe_cc ( class(gwemodeltype this,
integer(i4b), intent(in)  innertot,
integer(i4b), intent(in)  kiter,
integer(i4b), intent(in)  iend,
integer(i4b), intent(in)  icnvgmod,
character(len=lenpakloc), intent(inout)  cpak,
integer(i4b), intent(inout)  ipak,
real(dp), intent(inout)  dpak 
)
private

If MVR/MVT is active, this subroutine calls the MVR convergence check subroutines.

Definition at line 466 of file gwe.f90.

467  ! -- dummy
468  class(GweModelType) :: this
469  integer(I4B), intent(in) :: innertot
470  integer(I4B), intent(in) :: kiter
471  integer(I4B), intent(in) :: iend
472  integer(I4B), intent(in) :: icnvgmod
473  character(len=LENPAKLOC), intent(inout) :: cpak
474  integer(I4B), intent(inout) :: ipak
475  real(DP), intent(inout) :: dpak
476  !
477  ! -- If mover is on, then at least 2 outers required
478  if (this%inmvt > 0) call this%mvt%mvt_cc(kiter, iend, icnvgmod, cpak, dpak)

◆ gwe_cf()

subroutine gwemodule::gwe_cf ( class(gwemodeltype this,
integer(i4b), intent(in)  kiter 
)

This subroutine calls the attached packages' calculate coefficients subroutines

Definition at line 402 of file gwe.f90.

403  ! -- dummy
404  class(GweModelType) :: this
405  integer(I4B), intent(in) :: kiter
406  ! -- local
407  class(BndType), pointer :: packobj
408  integer(I4B) :: ip
409  !
410  ! -- Call package cf routines
411  do ip = 1, this%bndlist%Count()
412  packobj => getbndfromlist(this%bndlist, ip)
413  call packobj%bnd_cf()
414  end do
Here is the call graph for this function:

◆ gwe_cq()

subroutine gwemodule::gwe_cq ( class(gwemodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)
private

This subroutine calls the attached packages' intercell flows (flow ja)

Definition at line 485 of file gwe.f90.

486  ! -- modules
487  use sparsemodule, only: csr_diagsum
488  ! -- dummy
489  class(GweModelType) :: this
490  integer(I4B), intent(in) :: icnvg
491  integer(I4B), intent(in) :: isuppress_output
492  ! -- local
493  integer(I4B) :: i
494  integer(I4B) :: ip
495  class(BndType), pointer :: packobj
496  !
497  ! -- Construct the flowja array. Flowja is calculated each time, even if
498  ! output is suppressed. (flowja is positive into a cell.) The diagonal
499  ! position of the flowja array will contain the flow residual after
500  ! these routines are called, so each package is responsible for adding
501  ! its flow to this diagonal position.
502  do i = 1, this%nja
503  this%flowja(i) = dzero
504  end do
505  if (this%inadv > 0) call this%adv%adv_cq(this%x, this%flowja)
506  if (this%incnd > 0) call this%cnd%cnd_cq(this%x, this%flowja)
507  if (this%inest > 0) call this%est%est_cq(this%dis%nodes, this%x, this%xold, &
508  this%flowja)
509  if (this%inssm > 0) call this%ssm%ssm_cq(this%flowja)
510  if (this%infmi > 0) call this%fmi%fmi_cq(this%x, this%flowja)
511  !
512  ! -- Go through packages and call cq routines. cf() routines are called
513  ! first to regenerate non-linear terms to be consistent with the final
514  ! conc solution.
515  do ip = 1, this%bndlist%Count()
516  packobj => getbndfromlist(this%bndlist, ip)
517  call packobj%bnd_cf()
518  call packobj%bnd_cq(this%x, this%flowja)
519  end do
520  !
521  ! -- Finalize calculation of flowja by adding face flows to the diagonal.
522  ! This results in the flow residual being stored in the diagonal
523  ! position for each cell.
524  call csr_diagsum(this%dis%con%ia, this%flowja)
subroutine csr_diagsum(ia, flowja)
Definition: Sparse.f90:263
Here is the call graph for this function:

◆ gwe_cr()

subroutine, public gwemodule::gwe_cr ( character(len=*), intent(in)  filename,
integer(i4b), intent(in)  id,
character(len=*), intent(in)  modelname 
)
Parameters
[in]filenameinput file
[in]idconsecutive model number listed in mfsim.nam
[in]modelnamename of the model

Definition at line 97 of file gwe.f90.

98  ! -- modules
99  use listsmodule, only: basemodellist
105  use budgetmodule, only: budget_cr
107  ! -- dummy
108  character(len=*), intent(in) :: filename !< input file
109  integer(I4B), intent(in) :: id !< consecutive model number listed in mfsim.nam
110  character(len=*), intent(in) :: modelname !< name of the model
111  ! -- local
112  integer(I4B) :: indis
113  type(GweModelType), pointer :: this
114  class(BaseModelType), pointer :: model
115  !
116  ! -- Allocate a new GWE Model (this) and add it to basemodellist
117  allocate (this)
118  !
119  ! -- Set memory path before allocation in memory manager can be done
120  this%memoryPath = create_mem_path(modelname)
121  !
122  ! -- Allocate scalars and add model to basemodellist
123  call this%allocate_scalars(modelname)
124  !
125  ! -- Set labels for transport model - needed by create_packages() below
126  call this%set_tsp_labels(this%macronym, dvt, dvu, dvua)
127  !
128  model => this
129  call addbasemodeltolist(basemodellist, model)
130  !
131  ! -- Instantiate shared data container
132  call gweshared_dat_cr(this%gwecommon)
133  !
134  ! -- Call parent class routine
135  call this%tsp_cr(filename, id, modelname, 'GWE', indis)
136  !
137  ! -- Create model packages
138  call this%create_packages(indis)
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:160
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public budget_cr(this, name_model)
@ brief Create a new budget object
Definition: Budget.f90:84
subroutine, public gweshared_dat_cr(gwe_dat)
Allocate the shared data.
type(listtype), public basemodellist
Definition: mf6lists.f90:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gwe_da()

subroutine gwemodule::gwe_da ( class(gwemodeltype this)
private

Deallocate memory at conclusion of model run

Definition at line 581 of file gwe.f90.

582  ! -- modules
586  ! -- dummy
587  class(GweModelType) :: this
588  ! -- local
589  integer(I4B) :: ip
590  class(BndType), pointer :: packobj
591  !
592  ! -- Deallocate idm memory
593  call memorystore_remove(this%name, 'NAM', idm_context)
594  call memorystore_remove(component=this%name, context=idm_context)
595  !
596  ! -- Internal flow packages deallocate
597  call this%dis%dis_da()
598  call this%ic%ic_da()
599  call this%fmi%fmi_da()
600  call this%adv%adv_da()
601  call this%cnd%cnd_da()
602  call this%ssm%ssm_da()
603  call this%est%est_da()
604  call this%mvt%mvt_da()
605  call this%budget%budget_da()
606  call this%oc%oc_da()
607  call this%obs%obs_da()
608  call this%gwecommon%gweshared_dat_da()
609  !
610  ! -- Internal package objects
611  deallocate (this%dis)
612  deallocate (this%ic)
613  deallocate (this%fmi)
614  deallocate (this%adv)
615  deallocate (this%cnd)
616  deallocate (this%ssm)
617  deallocate (this%est)
618  deallocate (this%mvt)
619  deallocate (this%budget)
620  deallocate (this%oc)
621  deallocate (this%obs)
622  nullify (this%gwecommon)
623  !
624  ! -- Boundary packages
625  do ip = 1, this%bndlist%Count()
626  packobj => getbndfromlist(this%bndlist, ip)
627  call packobj%bnd_da()
628  deallocate (packobj)
629  end do
630  !
631  ! -- Scalars
632  call mem_deallocate(this%inest)
633  call mem_deallocate(this%incnd)
634  !
635  ! -- Parent class members
636  call this%TransportModelType%tsp_da()
637  !
638  ! -- NumericalModelType
639  call this%NumericalModelType%model_da()
subroutine, public memorystore_remove(component, subcomponent, context)
Here is the call graph for this function:

◆ gwe_df()

subroutine gwemodule::gwe_df ( class(gwemodeltype this)

This subroutine defines a gwe model type. Steps include:

  • call df routines for each package
  • set variables and pointers

Definition at line 147 of file gwe.f90.

148  ! -- modules
149  use simmodule, only: store_error
150  ! -- dummy
151  class(GweModelType) :: this
152  ! -- local
153  integer(I4B) :: ip
154  class(BndType), pointer :: packobj
155  !
156  ! -- Define packages and utility objects
157  call this%dis%dis_df()
158  call this%fmi%fmi_df(this%dis, 0)
159  if (this%inmvt > 0) call this%mvt%mvt_df(this%dis)
160  if (this%inadv > 0) call this%adv%adv_df()
161  if (this%incnd > 0) call this%cnd%cnd_df(this%dis)
162  if (this%inssm > 0) call this%ssm%ssm_df()
163  call this%oc%oc_df()
164  call this%budget%budget_df(niunit_gwe, this%depvarunit, &
165  this%depvarunitabbrev)
166  !
167  ! -- Check for SSM package
168  if (this%inssm == 0) then
169  if (this%fmi%nflowpack > 0) then
170  call store_error('Flow model has boundary packages, but there &
171  &is no SSM package. The SSM package must be activated.', &
172  terminate=.true.)
173  end if
174  end if
175  !
176  ! -- Assign or point model members to dis members
177  this%neq = this%dis%nodes
178  this%nja = this%dis%nja
179  this%ia => this%dis%con%ia
180  this%ja => this%dis%con%ja
181  !
182  ! -- Allocate model arrays, now that neq and nja are assigned
183  call this%allocate_arrays()
184  !
185  ! -- Define packages and assign iout for time series managers
186  do ip = 1, this%bndlist%Count()
187  packobj => getbndfromlist(this%bndlist, ip)
188  call packobj%bnd_df(this%neq, this%dis)
189  packobj%TsManager%iout = this%iout
190  packobj%TasManager%iout = this%iout
191  end do
192  !
193  ! -- Store information needed for observations
194  call this%obs%obs_df(this%iout, this%name, 'GWE', this%dis)
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:

◆ gwe_dt()

subroutine gwemodule::gwe_dt ( class(gwemodeltype this)

Calculate the maximum allowable time step size subject to time-step constraints. If adaptive time steps are used, then the time step used will be no larger than dtmax calculated here.

Definition at line 327 of file gwe.f90.

328  use tdismodule, only: kstp, kper
330  ! dummy
331  class(GweModelType) :: this
332  ! local
333  real(DP) :: dtmax
334  character(len=LINELENGTH) :: msg
335  dtmax = dnodata
336 
337  ! advection package courant stability
338  call this%adv%adv_dt(dtmax, msg, this%est%porosity)
339  if (msg /= '') then
340  call ats_submit_delt(kstp, kper, dtmax, msg)
341  end if
subroutine, public ats_submit_delt(kstp, kper, dt, sloc, idir)
@ brief Allow and external caller to submit preferred time step
Definition: ats.f90:493
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ gwe_fc()

subroutine gwemodule::gwe_fc ( class(gwemodeltype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), intent(in)  inwtflag 
)
private

This subroutine calls the attached packages' fill coefficients subroutines

Definition at line 422 of file gwe.f90.

423  ! -- dummy
424  class(GweModelType) :: this
425  integer(I4B), intent(in) :: kiter
426  class(MatrixBaseType), pointer :: matrix_sln
427  integer(I4B), intent(in) :: inwtflag
428  ! -- local
429  class(BndType), pointer :: packobj
430  integer(I4B) :: ip
431  !
432  ! -- Call fc routines
433  call this%fmi%fmi_fc(this%dis%nodes, this%xold, this%nja, matrix_sln, &
434  this%idxglo, this%rhs)
435  if (this%inmvt > 0) then
436  call this%mvt%mvt_fc(this%x, this%x)
437  end if
438  if (this%inest > 0) then
439  call this%est%est_fc(this%dis%nodes, this%xold, this%nja, matrix_sln, &
440  this%idxglo, this%x, this%rhs, kiter)
441  end if
442  if (this%inadv > 0) then
443  call this%adv%adv_fc(this%dis%nodes, matrix_sln, this%idxglo, this%x, &
444  this%rhs)
445  end if
446  if (this%incnd > 0) then
447  call this%cnd%cnd_fc(kiter, this%dis%nodes, this%nja, matrix_sln, &
448  this%idxglo, this%rhs, this%x)
449  end if
450  if (this%inssm > 0) then
451  call this%ssm%ssm_fc(matrix_sln, this%idxglo, this%rhs)
452  end if
453  !
454  ! -- Packages
455  do ip = 1, this%bndlist%Count()
456  packobj => getbndfromlist(this%bndlist, ip)
457  call packobj%bnd_fc(this%rhs, this%ia, this%idxglo, matrix_sln)
458  end do
Here is the call graph for this function:

◆ gwe_get_iasym()

integer(i4b) function gwemodule::gwe_get_iasym ( class(gwemodeltype this)

Definition at line 664 of file gwe.f90.

665  class(GweModelType) :: this
666  ! -- local
667  integer(I4B) :: iasym
668  integer(I4B) :: ip
669  class(BndType), pointer :: packobj
670  !
671  ! -- Start by setting iasym to zero
672  iasym = 0
673  !
674  ! -- ADV
675  if (this%inadv > 0) then
676  if (this%adv%iasym /= 0) iasym = 1
677  end if
678  !
679  ! -- CND
680  if (this%incnd > 0) then
681  if (this%cnd%ixt3d /= 0) iasym = 1
682  end if
683  !
684  ! -- Check for any packages that introduce matrix asymmetry
685  do ip = 1, this%bndlist%Count()
686  packobj => getbndfromlist(this%bndlist, ip)
687  if (packobj%iasym /= 0) iasym = 1
688  end do
Here is the call graph for this function:

◆ gwe_mc()

subroutine gwemodule::gwe_mc ( class(gwemodeltype this,
class(matrixbasetype), pointer  matrix_sln 
)
Parameters
matrix_slnglobal system matrix

Definition at line 224 of file gwe.f90.

225  ! -- dummy
226  class(GweModelType) :: this
227  class(MatrixBaseType), pointer :: matrix_sln !< global system matrix
228  ! -- local
229  class(BndType), pointer :: packobj
230  integer(I4B) :: ip
231  !
232  ! -- Find the position of each connection in the global ia, ja structure
233  ! and store them in idxglo.
234  call this%dis%dis_mc(this%moffset, this%idxglo, matrix_sln)
235  !
236  if (this%incnd > 0) call this%cnd%cnd_mc(this%moffset, matrix_sln)
237  !
238  ! -- Map any package connections
239  do ip = 1, this%bndlist%Count()
240  packobj => getbndfromlist(this%bndlist, ip)
241  call packobj%bnd_mc(this%moffset, matrix_sln)
242  end do
Here is the call graph for this function:

◆ gwe_ot_flow()

subroutine gwemodule::gwe_ot_flow ( class(gwemodeltype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(in)  icbcun 
)
private

Save and print flows

Definition at line 564 of file gwe.f90.

565  ! dummy
566  class(GweModelType) :: this
567  integer(I4B), intent(in) :: icbcfl
568  integer(I4B), intent(in) :: ibudfl
569  integer(I4B), intent(in) :: icbcun
570  ! local
571 
572  if (this%inest > 0) call this%est%est_ot_flow(icbcfl, icbcun)
573  call this%TransportModelType%tsp_ot_flow(icbcfl, ibudfl, icbcun)
574 

◆ gwe_rp()

subroutine gwemodule::gwe_rp ( class(gwemodeltype this)

This subroutine calls the attached packages' read and prepare routines

Definition at line 294 of file gwe.f90.

295  ! -- modules
296  use tdismodule, only: readnewdata
297  ! -- dummy
298  class(GweModelType) :: this
299  ! -- local
300  class(BndType), pointer :: packobj
301  integer(I4B) :: ip
302  !
303  ! -- In fmi, check for mvt and mvrbudobj consistency
304  call this%fmi%fmi_rp(this%inmvt)
305  if (this%inmvt > 0) call this%mvt%mvt_rp()
306  !
307  ! -- Check with TDIS on whether or not it is time to RP
308  if (.not. readnewdata) return
309  !
310  ! -- Read and prepare
311  if (this%inoc > 0) call this%oc%oc_rp()
312  if (this%inssm > 0) call this%ssm%ssm_rp()
313  do ip = 1, this%bndlist%Count()
314  packobj => getbndfromlist(this%bndlist, ip)
315  call packobj%bnd_rp()
316  call packobj%bnd_rp_log()
317  call packobj%bnd_rp_obs()
318  end do
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:26
Here is the call graph for this function:

◆ package_create()

subroutine gwemodule::package_create ( class(gwemodeltype this,
character(len=*), intent(in)  filtyp,
integer(i4b), intent(in)  ipakid,
integer(i4b), intent(in)  ipaknum,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  mempath,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout 
)

This subroutine calls the package create routines for packages activated by the user.

Definition at line 720 of file gwe.f90.

722  ! -- modules
723  use simmodule, only: store_error
724  use gwectpmodule, only: ctp_create
725  use gweeslmodule, only: esl_create
726  use gwelkemodule, only: lke_create
727  use gwesfemodule, only: sfe_create
728  use gwemwemodule, only: mwe_create
729  use gweuzemodule, only: uze_create
730  use apimodule, only: api_create
731  ! -- dummy
732  class(GweModelType) :: this
733  character(len=*), intent(in) :: filtyp
734  character(len=LINELENGTH) :: errmsg
735  integer(I4B), intent(in) :: ipakid
736  integer(I4B), intent(in) :: ipaknum
737  character(len=*), intent(in) :: pakname
738  character(len=*), intent(in) :: mempath
739  integer(I4B), intent(in) :: inunit
740  integer(I4B), intent(in) :: iout
741  ! -- local
742  class(BndType), pointer :: packobj
743  class(BndType), pointer :: packobj2
744  integer(I4B) :: ip
745  !
746  ! -- This part creates the package object
747  select case (filtyp)
748  case ('CTP6')
749  call ctp_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
750  pakname, this%depvartype, mempath)
751  case ('ESL6')
752  call esl_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
753  pakname, this%gwecommon, mempath)
754  case ('LKE6')
755  call lke_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
756  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
757  this%depvartype, this%depvarunit, this%depvarunitabbrev)
758  case ('SFE6')
759  call sfe_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
760  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
761  this%depvartype, this%depvarunit, this%depvarunitabbrev)
762  case ('MWE6')
763  call mwe_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
764  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
765  this%depvartype, this%depvarunit, this%depvarunitabbrev)
766  case ('UZE6')
767  call uze_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
768  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
769  this%depvartype, this%depvarunit, this%depvarunitabbrev)
770  !case('API6')
771  ! call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
772  ! pakname)
773  case default
774  write (errmsg, *) 'Invalid package type: ', filtyp
775  call store_error(errmsg, terminate=.true.)
776  end select
777  !
778  ! -- Packages is the bndlist that is associated with the parent model
779  ! -- The following statement puts a pointer to this package in the ipakid
780  ! -- position of packages.
781  do ip = 1, this%bndlist%Count()
782  packobj2 => getbndfromlist(this%bndlist, ip)
783  if (packobj2%packName == pakname) then
784  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
785  'already exists: ', trim(pakname)
786  call store_error(errmsg, terminate=.true.)
787  end if
788  end do
789  call addbndtolist(this%bndlist, packobj)
This module contains the API package methods.
Definition: gwf-api.f90:12
subroutine, public api_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, mempath)
@ brief Create a new package object
Definition: gwf-api.f90:51
subroutine, public ctp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype, mempath)
Create a new constant temperature package.
Definition: gwe-ctp.f90:57
subroutine, public esl_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, gwecommon, input_mempath)
Create an energy source loading package.
Definition: gwe-esl.f90:50
subroutine, public lke_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
Create a new lke package.
Definition: gwe-lke.f90:107
subroutine, public mwe_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
Create new MWE package.
Definition: gwe-mwe.f90:100
subroutine, public sfe_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
Create a new sfe package.
Definition: gwe-sfe.f90:112
subroutine, public uze_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, gwecommon, dvt, dvu, dvua)
Create a new UZE package.
Definition: gwe-uze.f90:96
Here is the call graph for this function:

Variable Documentation

◆ dvt

character(len=lenvarname), parameter gwemodule::dvt = 'TEMPERATURE '
private

Definition at line 29 of file gwe.f90.

29  character(len=LENVARNAME), parameter :: dvt = 'TEMPERATURE ' !< dependent variable type, varies based on model type

◆ dvu

character(len=lenvarname), parameter gwemodule::dvu = 'ENERGY '
private

Definition at line 30 of file gwe.f90.

30  character(len=LENVARNAME), parameter :: dvu = 'ENERGY ' !< dependent variable unit of measure, either "mass" or "energy"

◆ dvua

character(len=lenvarname), parameter gwemodule::dvua = 'E '
private

Definition at line 31 of file gwe.f90.

31  character(len=LENVARNAME), parameter :: dvua = 'E ' !< abbreviation of the dependent variable unit of measure, either "M" or "E"

◆ gwe_basepkg

character(len=lenpackagetype), dimension(gwe_nbasepkg), public gwemodule::gwe_basepkg

Definition at line 72 of file gwe.f90.

72  character(len=LENPACKAGETYPE), dimension(GWE_NBASEPKG) :: GWE_BASEPKG

◆ gwe_multipkg

character(len=lenpackagetype), dimension(gwe_nmultipkg), public gwemodule::gwe_multipkg

Definition at line 85 of file gwe.f90.

85  character(len=LENPACKAGETYPE), dimension(GWE_NMULTIPKG) :: GWE_MULTIPKG

◆ gwe_nbasepkg

integer(i4b), parameter, public gwemodule::gwe_nbasepkg = 50

GWE6 model base package types. Only listed packages are candidates for input and these will be loaded in the order specified.

Definition at line 71 of file gwe.f90.

71  integer(I4B), parameter :: GWE_NBASEPKG = 50

◆ gwe_nmultipkg

integer(i4b), parameter, public gwemodule::gwe_nmultipkg = 50

GWE6 model multi-instance package types. Only listed packages are candidates for input and these will be loaded in the order specified.

Definition at line 84 of file gwe.f90.

84  integer(I4B), parameter :: GWE_NMULTIPKG = 50

◆ niunit_gwe

integer(i4b), parameter gwemodule::niunit_gwe = GWE_NBASEPKG + GWE_NMULTIPKG
private

Definition at line 91 of file gwe.f90.

91  integer(I4B), parameter :: NIUNIT_GWE = gwe_nbasepkg + gwe_nmultipkg