MODFLOW 6  version 6.6.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_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 672 of file gwe.f90.

673  ! -- modules
675  ! -- dummy
676  class(GweModelType) :: this
677  character(len=*), intent(in) :: modelname
678  !
679  ! -- Allocate parent class scalars
680  call this%allocate_tsp_scalars(modelname)
681  !
682  ! -- Allocate members that are part of model class
683  call mem_allocate(this%inest, 'INEST', this%memoryPath)
684  call mem_allocate(this%incnd, 'INCND', this%memoryPath)
685  !
686  this%inest = 0
687  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 770 of file gwe.f90.

771  ! -- dummy
772  class(*), pointer :: model !< The object to be cast
773  ! -- return
774  class(GweModelType), pointer :: gwemodel !< The GWE model
775  !
776  gwemodel => null()
777  if (.not. associated(model)) return
778  select type (model)
779  type is (gwemodeltype)
780  gwemodel => model
781  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 786 of file gwe.f90.

788  ! -- modules
791  ! -- dummy
792  class(GweModelType) :: this
793  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
794  type(CharacterStringType), dimension(:), contiguous, &
795  pointer, intent(inout) :: pkgtypes
796  type(CharacterStringType), dimension(:), contiguous, &
797  pointer, intent(inout) :: pkgnames
798  type(CharacterStringType), dimension(:), contiguous, &
799  pointer, intent(inout) :: mempaths
800  integer(I4B), dimension(:), contiguous, &
801  pointer, intent(inout) :: inunits
802  ! -- local
803  integer(I4B) :: ipakid, ipaknum
804  character(len=LENFTYPE) :: pkgtype, bndptype
805  character(len=LENPACKAGENAME) :: pkgname
806  character(len=LENMEMPATH) :: mempath
807  integer(I4B), pointer :: inunit
808  integer(I4B) :: n
809  !
810  if (allocated(bndpkgs)) then
811  !
812  ! -- Create stress packages
813  ipakid = 1
814  bndptype = ''
815  do n = 1, size(bndpkgs)
816  !
817  pkgtype = pkgtypes(bndpkgs(n))
818  pkgname = pkgnames(bndpkgs(n))
819  mempath = mempaths(bndpkgs(n))
820  inunit => inunits(bndpkgs(n))
821  !
822  if (bndptype /= pkgtype) then
823  ipaknum = 1
824  bndptype = pkgtype
825  end if
826  !
827  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
828  inunit, this%iout)
829  ipakid = ipakid + 1
830  ipaknum = ipaknum + 1
831  end do
832  !
833  ! -- Cleanup
834  deallocate (bndpkgs)
835  end if
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
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 840 of file gwe.f90.

841  ! -- modules
848  use gweestmodule, only: est_cr
849  use gwecndmodule, only: cnd_cr
850  ! -- dummy
851  class(GweModelType) :: this
852  integer(I4B), intent(in) :: indis
853  ! -- local
854  type(CharacterStringType), dimension(:), contiguous, &
855  pointer :: pkgtypes => null()
856  type(CharacterStringType), dimension(:), contiguous, &
857  pointer :: pkgnames => null()
858  type(CharacterStringType), dimension(:), contiguous, &
859  pointer :: mempaths => null()
860  integer(I4B), dimension(:), contiguous, &
861  pointer :: inunits => null()
862  character(len=LENMEMPATH) :: model_mempath
863  character(len=LENFTYPE) :: pkgtype
864  character(len=LENPACKAGENAME) :: pkgname
865  character(len=LENMEMPATH) :: mempath
866  integer(I4B), pointer :: inunit
867  integer(I4B), dimension(:), allocatable :: bndpkgs
868  integer(I4B) :: n
869  character(len=LENMEMPATH) :: mempathcnd = ''
870  !
871  ! -- Set input memory paths, input/model and input/model/namfile
872  model_mempath = create_mem_path(component=this%name, context=idm_context)
873  !
874  ! -- Set pointers to model path package info
875  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
876  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
877  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
878  call mem_setptr(inunits, 'INUNITS', model_mempath)
879  !
880  do n = 1, size(pkgtypes)
881  !
882  ! -- Attributes for this input package
883  pkgtype = pkgtypes(n)
884  pkgname = pkgnames(n)
885  mempath = mempaths(n)
886  inunit => inunits(n)
887  !
888  ! -- Create dis package as it is a prerequisite for other packages
889  select case (pkgtype)
890  case ('EST6')
891  this%inest = inunit
892  case ('CND6')
893  this%incnd = 1
894  mempathcnd = mempath
895  case ('CTP6', 'ESL6', 'LKE6', 'SFE6', &
896  'MWE6', 'UZE6', 'API6')
897  call expandarray(bndpkgs)
898  bndpkgs(size(bndpkgs)) = n
899  case default
900  ! TODO
901  end select
902  end do
903  !
904  ! -- Create packages that are tied directly to model
905  call est_cr(this%est, this%name, this%inest, this%iout, this%fmi, &
906  this%eqnsclfac, this%gwecommon)
907  call cnd_cr(this%cnd, this%name, mempathcnd, this%incnd, this%iout, &
908  this%fmi, this%eqnsclfac, this%gwecommon)
909  !
910  ! -- Check to make sure that required ftype's have been specified
911  call this%ftype_check(indis, this%inest)
912  !
913  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, inunit, iout, fmi, eqnsclfac, gwecommon)
@ brief Create a new EST package object
Definition: gwe-est.f90:87
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 197 of file gwe.f90.

198  ! -- modules
199  use sparsemodule, only: sparsematrix
200  ! -- dummy
201  class(GweModelType) :: this
202  type(sparsematrix), intent(inout) :: sparse
203  ! -- local
204  class(BndType), pointer :: packobj
205  integer(I4B) :: ip
206  !
207  ! -- Add the internal connections of this model to sparse
208  call this%dis%dis_ac(this%moffset, sparse)
209  if (this%incnd > 0) &
210  call this%cnd%cnd_ac(this%moffset, sparse)
211  !
212  ! -- Add any package connections
213  do ip = 1, this%bndlist%Count()
214  packobj => getbndfromlist(this%bndlist, ip)
215  call packobj%bnd_ac(this%moffset, sparse)
216  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 322 of file gwe.f90.

323  ! -- modules
325  ! -- dummy
326  class(GweModelType) :: this
327  class(BndType), pointer :: packobj
328  ! -- local
329  integer(I4B) :: irestore
330  integer(I4B) :: ip, n
331  !
332  ! -- Reset state variable
333  irestore = 0
334  if (ifailedstepretry > 0) irestore = 1
335  if (irestore == 0) then
336  !
337  ! -- Copy x into xold
338  do n = 1, this%dis%nodes
339  if (this%ibound(n) == 0) then
340  this%xold(n) = dzero
341  else
342  this%xold(n) = this%x(n)
343  end if
344  end do
345  else
346  !
347  ! -- Copy xold into x if this time step is a redo
348  do n = 1, this%dis%nodes
349  this%x(n) = this%xold(n)
350  end do
351  end if
352  !
353  ! -- Advance fmi
354  call this%fmi%fmi_ad(this%x)
355  !
356  ! -- Advance
357  if (this%incnd > 0) call this%cnd%cnd_ad()
358  if (this%inssm > 0) call this%ssm%ssm_ad()
359  do ip = 1, this%bndlist%Count()
360  packobj => getbndfromlist(this%bndlist, ip)
361  call packobj%bnd_ad()
362  if (isimcheck > 0) then
363  call packobj%bnd_ck()
364  end if
365  end do
366  !
367  ! -- Push simulated values to preceding time/subtime step
368  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 249 of file gwe.f90.

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

508  use constantsmodule, only: dzero
509  ! -- dummy
510  class(GweModelType) :: this
511  integer(I4B), intent(in) :: icnvg
512  integer(I4B), intent(in) :: isuppress_output
513  ! -- local
514  integer(I4B) :: ip
515  class(BndType), pointer :: packobj
516  !
517  ! -- Save the solution convergence flag
518  this%icnvg = icnvg
519  !
520  ! -- Budget routines (start by resetting). Sole purpose of this section
521  ! is to add in and outs to model budget. All ins and out for a model
522  ! should be added here to this%budget. In a subsequent exchange call,
523  ! exchange flows might also be added.
524  call this%budget%reset()
525  if (this%inest > 0) call this%est%est_bd(isuppress_output, this%budget)
526  if (this%inssm > 0) call this%ssm%ssm_bd(isuppress_output, this%budget)
527  if (this%infmi > 0) call this%fmi%fmi_bd(isuppress_output, this%budget)
528  if (this%inmvt > 0) call this%mvt%mvt_bd(this%x, this%x)
529  do ip = 1, this%bndlist%Count()
530  packobj => getbndfromlist(this%bndlist, ip)
531  call packobj%bnd_bd(this%budget)
532  end do
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
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 623 of file gwe.f90.

624  ! -- modules
625  use constantsmodule, only: lenbudtxt
626  use tdismodule, only: delt
627  ! -- dummy
628  class(GweModelType) :: this
629  real(DP), dimension(:, :), intent(in) :: budterm
630  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
631  character(len=*), intent(in) :: rowlabel
632  !
633  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 440 of file gwe.f90.

441  ! -- dummy
442  class(GweModelType) :: this
443  integer(I4B), intent(in) :: innertot
444  integer(I4B), intent(in) :: kiter
445  integer(I4B), intent(in) :: iend
446  integer(I4B), intent(in) :: icnvgmod
447  character(len=LENPAKLOC), intent(inout) :: cpak
448  integer(I4B), intent(inout) :: ipak
449  real(DP), intent(inout) :: dpak
450  !
451  ! -- If mover is on, then at least 2 outers required
452  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 376 of file gwe.f90.

377  ! -- dummy
378  class(GweModelType) :: this
379  integer(I4B), intent(in) :: kiter
380  ! -- local
381  class(BndType), pointer :: packobj
382  integer(I4B) :: ip
383  !
384  ! -- Call package cf routines
385  do ip = 1, this%bndlist%Count()
386  packobj => getbndfromlist(this%bndlist, ip)
387  call packobj%bnd_cf()
388  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 459 of file gwe.f90.

460  ! -- modules
461  use sparsemodule, only: csr_diagsum
462  ! -- dummy
463  class(GweModelType) :: this
464  integer(I4B), intent(in) :: icnvg
465  integer(I4B), intent(in) :: isuppress_output
466  ! -- local
467  integer(I4B) :: i
468  integer(I4B) :: ip
469  class(BndType), pointer :: packobj
470  !
471  ! -- Construct the flowja array. Flowja is calculated each time, even if
472  ! output is suppressed. (flowja is positive into a cell.) The diagonal
473  ! position of the flowja array will contain the flow residual after
474  ! these routines are called, so each package is responsible for adding
475  ! its flow to this diagonal position.
476  do i = 1, this%nja
477  this%flowja(i) = dzero
478  end do
479  if (this%inadv > 0) call this%adv%adv_cq(this%x, this%flowja)
480  if (this%incnd > 0) call this%cnd%cnd_cq(this%x, this%flowja)
481  if (this%inest > 0) call this%est%est_cq(this%dis%nodes, this%x, this%xold, &
482  this%flowja)
483  if (this%inssm > 0) call this%ssm%ssm_cq(this%flowja)
484  if (this%infmi > 0) call this%fmi%fmi_cq(this%x, this%flowja)
485  !
486  ! -- Go through packages and call cq routines. cf() routines are called
487  ! first to regenerate non-linear terms to be consistent with the final
488  ! conc solution.
489  do ip = 1, this%bndlist%Count()
490  packobj => getbndfromlist(this%bndlist, ip)
491  call packobj%bnd_cf()
492  call packobj%bnd_cq(this%x, this%flowja)
493  end do
494  !
495  ! -- Finalize calculation of flowja by adding face flows to the diagonal.
496  ! This results in the flow residual being stored in the diagonal
497  ! position for each cell.
498  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 95 of file gwe.f90.

96  ! -- modules
97  use listsmodule, only: basemodellist
103  use budgetmodule, only: budget_cr
105  ! -- dummy
106  character(len=*), intent(in) :: filename !< input file
107  integer(I4B), intent(in) :: id !< consecutive model number listed in mfsim.nam
108  character(len=*), intent(in) :: modelname !< name of the model
109  ! -- local
110  integer(I4B) :: indis
111  type(GweModelType), pointer :: this
112  class(BaseModelType), pointer :: model
113  !
114  ! -- Allocate a new GWE Model (this) and add it to basemodellist
115  allocate (this)
116  !
117  ! -- Set memory path before allocation in memory manager can be done
118  this%memoryPath = create_mem_path(modelname)
119  !
120  ! -- Allocate scalars and add model to basemodellist
121  call this%allocate_scalars(modelname)
122  !
123  ! -- Set labels for transport model - needed by create_packages() below
124  call this%set_tsp_labels(this%macronym, dvt, dvu, dvua)
125  !
126  model => this
127  call addbasemodeltolist(basemodellist, model)
128  !
129  ! -- Instantiate shared data container
130  call gweshared_dat_cr(this%gwecommon)
131  !
132  ! -- Call parent class routine
133  call this%tsp_cr(filename, id, modelname, 'GWE', indis)
134  !
135  ! -- Create model packages
136  call this%create_packages(indis)
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:161
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 556 of file gwe.f90.

557  ! -- modules
561  ! -- dummy
562  class(GweModelType) :: this
563  ! -- local
564  integer(I4B) :: ip
565  class(BndType), pointer :: packobj
566  !
567  ! -- Deallocate idm memory
568  call memorystore_remove(this%name, 'NAM', idm_context)
569  call memorystore_remove(component=this%name, context=idm_context)
570  !
571  ! -- Internal flow packages deallocate
572  call this%dis%dis_da()
573  call this%ic%ic_da()
574  call this%fmi%fmi_da()
575  call this%adv%adv_da()
576  call this%cnd%cnd_da()
577  call this%ssm%ssm_da()
578  call this%est%est_da()
579  call this%mvt%mvt_da()
580  call this%budget%budget_da()
581  call this%oc%oc_da()
582  call this%obs%obs_da()
583  call this%gwecommon%gweshared_dat_da()
584  !
585  ! -- Internal package objects
586  deallocate (this%dis)
587  deallocate (this%ic)
588  deallocate (this%fmi)
589  deallocate (this%adv)
590  deallocate (this%cnd)
591  deallocate (this%ssm)
592  deallocate (this%est)
593  deallocate (this%mvt)
594  deallocate (this%budget)
595  deallocate (this%oc)
596  deallocate (this%obs)
597  nullify (this%gwecommon)
598  !
599  ! -- Boundary packages
600  do ip = 1, this%bndlist%Count()
601  packobj => getbndfromlist(this%bndlist, ip)
602  call packobj%bnd_da()
603  deallocate (packobj)
604  end do
605  !
606  ! -- Scalars
607  call mem_deallocate(this%inest)
608  call mem_deallocate(this%incnd)
609  !
610  ! -- Parent class members
611  call this%TransportModelType%tsp_da()
612  !
613  ! -- NumericalModelType
614  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 145 of file gwe.f90.

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

397  ! -- dummy
398  class(GweModelType) :: this
399  integer(I4B), intent(in) :: kiter
400  class(MatrixBaseType), pointer :: matrix_sln
401  integer(I4B), intent(in) :: inwtflag
402  ! -- local
403  class(BndType), pointer :: packobj
404  integer(I4B) :: ip
405  !
406  ! -- Call fc routines
407  call this%fmi%fmi_fc(this%dis%nodes, this%xold, this%nja, matrix_sln, &
408  this%idxglo, this%rhs)
409  if (this%inmvt > 0) then
410  call this%mvt%mvt_fc(this%x, this%x)
411  end if
412  if (this%inest > 0) then
413  call this%est%est_fc(this%dis%nodes, this%xold, this%nja, matrix_sln, &
414  this%idxglo, this%x, this%rhs, kiter)
415  end if
416  if (this%inadv > 0) then
417  call this%adv%adv_fc(this%dis%nodes, matrix_sln, this%idxglo, this%x, &
418  this%rhs)
419  end if
420  if (this%incnd > 0) then
421  call this%cnd%cnd_fc(kiter, this%dis%nodes, this%nja, matrix_sln, &
422  this%idxglo, this%rhs, this%x)
423  end if
424  if (this%inssm > 0) then
425  call this%ssm%ssm_fc(matrix_sln, this%idxglo, this%rhs)
426  end if
427  !
428  ! -- Packages
429  do ip = 1, this%bndlist%Count()
430  packobj => getbndfromlist(this%bndlist, ip)
431  call packobj%bnd_fc(this%rhs, this%ia, this%idxglo, matrix_sln)
432  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 639 of file gwe.f90.

640  class(GweModelType) :: this
641  ! -- local
642  integer(I4B) :: iasym
643  integer(I4B) :: ip
644  class(BndType), pointer :: packobj
645  !
646  ! -- Start by setting iasym to zero
647  iasym = 0
648  !
649  ! -- ADV
650  if (this%inadv > 0) then
651  if (this%adv%iasym /= 0) iasym = 1
652  end if
653  !
654  ! -- CND
655  if (this%incnd > 0) then
656  if (this%cnd%ixt3d /= 0) iasym = 1
657  end if
658  !
659  ! -- Check for any packages that introduce matrix asymmetry
660  do ip = 1, this%bndlist%Count()
661  packobj => getbndfromlist(this%bndlist, ip)
662  if (packobj%iasym /= 0) iasym = 1
663  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 222 of file gwe.f90.

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

Save and print flows

Definition at line 539 of file gwe.f90.

540  ! dummy
541  class(GweModelType) :: this
542  integer(I4B), intent(in) :: icbcfl
543  integer(I4B), intent(in) :: ibudfl
544  integer(I4B), intent(in) :: icbcun
545  ! local
546 
547  if (this%inest > 0) call this%est%est_ot_flow(icbcfl, icbcun)
548  call this%TransportModelType%tsp_ot_flow(icbcfl, ibudfl, icbcun)
549 

◆ gwe_rp()

subroutine gwemodule::gwe_rp ( class(gwemodeltype this)

This subroutine calls the attached packages' read and prepare routines

Definition at line 292 of file gwe.f90.

293  ! -- modules
294  use tdismodule, only: readnewdata
295  ! -- dummy
296  class(GweModelType) :: this
297  ! -- local
298  class(BndType), pointer :: packobj
299  integer(I4B) :: ip
300  !
301  ! -- In fmi, check for mvt and mvrbudobj consistency
302  call this%fmi%fmi_rp(this%inmvt)
303  if (this%inmvt > 0) call this%mvt%mvt_rp()
304  !
305  ! -- Check with TDIS on whether or not it is time to RP
306  if (.not. readnewdata) return
307  !
308  ! -- Read and prepare
309  if (this%inoc > 0) call this%oc%oc_rp()
310  if (this%inssm > 0) call this%ssm%ssm_rp()
311  do ip = 1, this%bndlist%Count()
312  packobj => getbndfromlist(this%bndlist, ip)
313  call packobj%bnd_rp()
314  call packobj%bnd_rp_obs()
315  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 695 of file gwe.f90.

697  ! -- modules
698  use constantsmodule, only: linelength
699  use simmodule, only: store_error
700  use gwectpmodule, only: ctp_create
701  use gweeslmodule, only: esl_create
702  use gwelkemodule, only: lke_create
703  use gwesfemodule, only: sfe_create
704  use gwemwemodule, only: mwe_create
705  use gweuzemodule, only: uze_create
706  use apimodule, only: api_create
707  ! -- dummy
708  class(GweModelType) :: this
709  character(len=*), intent(in) :: filtyp
710  character(len=LINELENGTH) :: errmsg
711  integer(I4B), intent(in) :: ipakid
712  integer(I4B), intent(in) :: ipaknum
713  character(len=*), intent(in) :: pakname
714  character(len=*), intent(in) :: mempath
715  integer(I4B), intent(in) :: inunit
716  integer(I4B), intent(in) :: iout
717  ! -- local
718  class(BndType), pointer :: packobj
719  class(BndType), pointer :: packobj2
720  integer(I4B) :: ip
721  !
722  ! -- This part creates the package object
723  select case (filtyp)
724  case ('CTP6')
725  call ctp_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
726  pakname, this%depvartype, mempath)
727  case ('ESL6')
728  call esl_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
729  pakname, this%gwecommon)
730  case ('LKE6')
731  call lke_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
732  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
733  this%depvartype, this%depvarunit, this%depvarunitabbrev)
734  case ('SFE6')
735  call sfe_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
736  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
737  this%depvartype, this%depvarunit, this%depvarunitabbrev)
738  case ('MWE6')
739  call mwe_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
740  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
741  this%depvartype, this%depvarunit, this%depvarunitabbrev)
742  case ('UZE6')
743  call uze_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
744  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
745  this%depvartype, this%depvarunit, this%depvarunitabbrev)
746  !case('API6')
747  ! call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
748  ! pakname)
749  case default
750  write (errmsg, *) 'Invalid package type: ', filtyp
751  call store_error(errmsg, terminate=.true.)
752  end select
753  !
754  ! -- Packages is the bndlist that is associated with the parent model
755  ! -- The following statement puts a pointer to this package in the ipakid
756  ! -- position of packages.
757  do ip = 1, this%bndlist%Count()
758  packobj2 => getbndfromlist(this%bndlist, ip)
759  if (packobj2%packName == pakname) then
760  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
761  'already exists: ', trim(pakname)
762  call store_error(errmsg, terminate=.true.)
763  end if
764  end do
765  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)
@ brief Create a new package object
Definition: gwf-api.f90:49
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)
Create an energy source loading package.
Definition: gwe-esl.f90:48
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:103
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:96
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:102
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:94
Here is the call graph for this function:

Variable Documentation

◆ dvt

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

Definition at line 28 of file gwe.f90.

28  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 29 of file gwe.f90.

29  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 30 of file gwe.f90.

30  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 70 of file gwe.f90.

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

◆ gwe_multipkg

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

Definition at line 83 of file gwe.f90.

83  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 69 of file gwe.f90.

69  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 82 of file gwe.f90.

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

◆ niunit_gwe

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

Definition at line 89 of file gwe.f90.

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