MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 (this)
 GWE Model Output. 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 731 of file gwe.f90.

732  ! -- modules
734  ! -- dummy
735  class(GweModelType) :: this
736  character(len=*), intent(in) :: modelname
737  !
738  ! -- Allocate parent class scalars
739  call this%allocate_tsp_scalars(modelname)
740  !
741  ! -- Allocate members that are part of model class
742  call mem_allocate(this%inest, 'INEST', this%memoryPath)
743  call mem_allocate(this%incnd, 'INCND', this%memoryPath)
744  !
745  this%inest = 0
746  this%incnd = 0
747  !
748  ! -- Return
749  return

◆ castasgwemodel()

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

Definition at line 835 of file gwe.f90.

836  ! -- dummy
837  class(*), pointer :: model !< The object to be cast
838  ! -- return
839  class(GweModelType), pointer :: gwemodel !< The GWE model
840  !
841  gwemodel => null()
842  if (.not. associated(model)) return
843  select type (model)
844  type is (gwemodeltype)
845  gwemodel => model
846  end select
847  !
848  ! -- Return
849  return
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 854 of file gwe.f90.

856  ! -- modules
859  ! -- dummy
860  class(GweModelType) :: this
861  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
862  type(CharacterStringType), dimension(:), contiguous, &
863  pointer, intent(inout) :: pkgtypes
864  type(CharacterStringType), dimension(:), contiguous, &
865  pointer, intent(inout) :: pkgnames
866  type(CharacterStringType), dimension(:), contiguous, &
867  pointer, intent(inout) :: mempaths
868  integer(I4B), dimension(:), contiguous, &
869  pointer, intent(inout) :: inunits
870  ! -- local
871  integer(I4B) :: ipakid, ipaknum
872  character(len=LENFTYPE) :: pkgtype, bndptype
873  character(len=LENPACKAGENAME) :: pkgname
874  character(len=LENMEMPATH) :: mempath
875  integer(I4B), pointer :: inunit
876  integer(I4B) :: n
877  !
878  if (allocated(bndpkgs)) then
879  !
880  ! -- Create stress packages
881  ipakid = 1
882  bndptype = ''
883  do n = 1, size(bndpkgs)
884  !
885  pkgtype = pkgtypes(bndpkgs(n))
886  pkgname = pkgnames(bndpkgs(n))
887  mempath = mempaths(bndpkgs(n))
888  inunit => inunits(bndpkgs(n))
889  !
890  if (bndptype /= pkgtype) then
891  ipaknum = 1
892  bndptype = pkgtype
893  end if
894  !
895  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
896  inunit, this%iout)
897  ipakid = ipakid + 1
898  ipaknum = ipaknum + 1
899  end do
900  !
901  ! -- Cleanup
902  deallocate (bndpkgs)
903  end if
904  !
905  ! -- Return
906  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
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 911 of file gwe.f90.

912  ! -- modules
919  use gweestmodule, only: est_cr
920  use gwecndmodule, only: cnd_cr
921  ! -- dummy
922  class(GweModelType) :: this
923  integer(I4B), intent(in) :: indis
924  ! -- local
925  type(CharacterStringType), dimension(:), contiguous, &
926  pointer :: pkgtypes => null()
927  type(CharacterStringType), dimension(:), contiguous, &
928  pointer :: pkgnames => null()
929  type(CharacterStringType), dimension(:), contiguous, &
930  pointer :: mempaths => null()
931  integer(I4B), dimension(:), contiguous, &
932  pointer :: inunits => null()
933  character(len=LENMEMPATH) :: model_mempath
934  character(len=LENFTYPE) :: pkgtype
935  character(len=LENPACKAGENAME) :: pkgname
936  character(len=LENMEMPATH) :: mempath
937  integer(I4B), pointer :: inunit
938  integer(I4B), dimension(:), allocatable :: bndpkgs
939  integer(I4B) :: n
940  character(len=LENMEMPATH) :: mempathcnd = ''
941  !
942  ! -- Set input memory paths, input/model and input/model/namfile
943  model_mempath = create_mem_path(component=this%name, context=idm_context)
944  !
945  ! -- Set pointers to model path package info
946  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
947  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
948  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
949  call mem_setptr(inunits, 'INUNITS', model_mempath)
950  !
951  do n = 1, size(pkgtypes)
952  !
953  ! -- Attributes for this input package
954  pkgtype = pkgtypes(n)
955  pkgname = pkgnames(n)
956  mempath = mempaths(n)
957  inunit => inunits(n)
958  !
959  ! -- Create dis package as it is a prerequisite for other packages
960  select case (pkgtype)
961  case ('EST6')
962  this%inest = inunit
963  case ('CND6')
964  this%incnd = 1
965  mempathcnd = mempath
966  case ('CTP6', 'ESL6', 'LKE6', 'SFE6', &
967  'MWE6', 'UZE6', 'API6')
968  call expandarray(bndpkgs)
969  bndpkgs(size(bndpkgs)) = n
970  case default
971  ! TODO
972  end select
973  end do
974  !
975  ! -- Create packages that are tied directly to model
976  call est_cr(this%est, this%name, this%inest, this%iout, this%fmi, &
977  this%eqnsclfac, this%gwecommon)
978  call cnd_cr(this%cnd, this%name, mempathcnd, this%incnd, this%iout, &
979  this%fmi, this%eqnsclfac, this%gwecommon)
980  !
981  ! -- Check to make sure that required ftype's have been specified
982  call this%ftype_check(indis, this%inest)
983  !
984  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
985  !
986  ! -- Return
987  return
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:89
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 207 of file gwe.f90.

208  ! -- modules
209  use sparsemodule, only: sparsematrix
210  ! -- dummy
211  class(GweModelType) :: this
212  type(sparsematrix), intent(inout) :: sparse
213  ! -- local
214  class(BndType), pointer :: packobj
215  integer(I4B) :: ip
216  !
217  ! -- Add the internal connections of this model to sparse
218  call this%dis%dis_ac(this%moffset, sparse)
219  if (this%incnd > 0) &
220  call this%cnd%cnd_ac(this%moffset, sparse)
221  !
222  ! -- Add any package connections
223  do ip = 1, this%bndlist%Count()
224  packobj => getbndfromlist(this%bndlist, ip)
225  call packobj%bnd_ac(this%moffset, sparse)
226  end do
227  !
228  ! -- Return
229  return
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 344 of file gwe.f90.

345  ! -- modules
347  ! -- dummy
348  class(GweModelType) :: this
349  class(BndType), pointer :: packobj
350  ! -- local
351  integer(I4B) :: irestore
352  integer(I4B) :: ip, n
353  !
354  ! -- Reset state variable
355  irestore = 0
356  if (ifailedstepretry > 0) irestore = 1
357  if (irestore == 0) then
358  !
359  ! -- Copy x into xold
360  do n = 1, this%dis%nodes
361  if (this%ibound(n) == 0) then
362  this%xold(n) = dzero
363  else
364  this%xold(n) = this%x(n)
365  end if
366  end do
367  else
368  !
369  ! -- Copy xold into x if this time step is a redo
370  do n = 1, this%dis%nodes
371  this%x(n) = this%xold(n)
372  end do
373  end if
374  !
375  ! -- Advance fmi
376  call this%fmi%fmi_ad(this%x)
377  !
378  ! -- Advance
379  if (this%incnd > 0) call this%cnd%cnd_ad()
380  if (this%inssm > 0) call this%ssm%ssm_ad()
381  do ip = 1, this%bndlist%Count()
382  packobj => getbndfromlist(this%bndlist, ip)
383  call packobj%bnd_ad()
384  if (isimcheck > 0) then
385  call packobj%bnd_ck()
386  end if
387  end do
388  !
389  ! -- Push simulated values to preceding time/subtime step
390  call this%obs%obs_ad()
391  !
392  ! -- Return
393  return
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 265 of file gwe.f90.

266  ! -- modules
267  use constantsmodule, only: dhnoflo
268  ! -- dummy
269  class(GweModelType) :: this
270  ! -- locals
271  integer(I4B) :: ip
272  class(BndType), pointer :: packobj
273  !
274  ! -- Allocate and read modules attached to model
275  call this%fmi%fmi_ar(this%ibound)
276  if (this%inmvt > 0) call this%mvt%mvt_ar()
277  if (this%inic > 0) call this%ic%ic_ar(this%x)
278  if (this%inest > 0) call this%est%est_ar(this%dis, this%ibound)
279  if (this%inadv > 0) call this%adv%adv_ar(this%dis, this%ibound)
280  if (this%incnd > 0) call this%cnd%cnd_ar(this%ibound, this%est%porosity)
281  if (this%inssm > 0) call this%ssm%ssm_ar(this%dis, this%ibound, this%x)
282  if (this%inobs > 0) call this%obs%tsp_obs_ar(this%ic, this%x, this%flowja)
283  !
284  ! -- Set governing equation scale factor
285  this%eqnsclfac = this%gwecommon%gwerhow * this%gwecommon%gwecpw
286  !
287  ! -- Call dis_ar to write binary grid file
288  !call this%dis%dis_ar(this%npf%icelltype)
289  !
290  ! -- Set up output control
291  call this%oc%oc_ar(this%x, this%dis, dhnoflo, this%depvartype)
292  call this%budget%set_ibudcsv(this%oc%ibudcsv)
293  !
294  ! -- Package input files now open, so allocate and read
295  do ip = 1, this%bndlist%Count()
296  packobj => getbndfromlist(this%bndlist, ip)
297  call packobj%set_pointers(this%dis%nodes, this%ibound, this%x, &
298  this%xold, this%flowja)
299  ! -- Read and allocate package
300  call packobj%bnd_ar()
301  end do
302  !
303  ! -- Return
304  return
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:92
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 544 of file gwe.f90.

545  use constantsmodule, only: dzero
546  ! -- dummy
547  class(GweModelType) :: this
548  integer(I4B), intent(in) :: icnvg
549  integer(I4B), intent(in) :: isuppress_output
550  ! -- local
551  integer(I4B) :: ip
552  class(BndType), pointer :: packobj
553  !
554  ! -- Save the solution convergence flag
555  this%icnvg = icnvg
556  !
557  ! -- Budget routines (start by resetting). Sole purpose of this section
558  ! is to add in and outs to model budget. All ins and out for a model
559  ! should be added here to this%budget. In a subsequent exchange call,
560  ! exchange flows might also be added.
561  call this%budget%reset()
562  if (this%inest > 0) call this%est%est_bd(isuppress_output, this%budget)
563  if (this%inssm > 0) call this%ssm%ssm_bd(isuppress_output, this%budget)
564  if (this%infmi > 0) call this%fmi%fmi_bd(isuppress_output, this%budget)
565  if (this%inmvt > 0) call this%mvt%mvt_bd(this%x, this%x)
566  do ip = 1, this%bndlist%Count()
567  packobj => getbndfromlist(this%bndlist, ip)
568  call packobj%bnd_bd(this%budget)
569  end do
570  !
571  ! -- Return
572  return
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
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 676 of file gwe.f90.

677  ! -- modules
678  use constantsmodule, only: lenbudtxt
679  use tdismodule, only: delt
680  ! -- dummy
681  class(GweModelType) :: this
682  real(DP), dimension(:, :), intent(in) :: budterm
683  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
684  character(len=*), intent(in) :: rowlabel
685  !
686  call this%budget%addentry(budterm, delt, budtxt, rowlabel=rowlabel)
687  !
688  ! -- Return
689  return
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:36
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 471 of file gwe.f90.

472  ! -- dummy
473  class(GweModelType) :: this
474  integer(I4B), intent(in) :: innertot
475  integer(I4B), intent(in) :: kiter
476  integer(I4B), intent(in) :: iend
477  integer(I4B), intent(in) :: icnvgmod
478  character(len=LENPAKLOC), intent(inout) :: cpak
479  integer(I4B), intent(inout) :: ipak
480  real(DP), intent(inout) :: dpak
481  !
482  ! -- If mover is on, then at least 2 outers required
483  if (this%inmvt > 0) call this%mvt%mvt_cc(kiter, iend, icnvgmod, cpak, dpak)
484  !
485  ! -- Return
486  return

◆ 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 401 of file gwe.f90.

402  ! -- dummy
403  class(GweModelType) :: this
404  integer(I4B), intent(in) :: kiter
405  ! -- local
406  class(BndType), pointer :: packobj
407  integer(I4B) :: ip
408  !
409  ! -- Call package cf routines
410  do ip = 1, this%bndlist%Count()
411  packobj => getbndfromlist(this%bndlist, ip)
412  call packobj%bnd_cf()
413  end do
414  !
415  ! -- Return
416  return
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 493 of file gwe.f90.

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

607  ! -- modules
611  ! -- dummy
612  class(GweModelType) :: this
613  ! -- local
614  integer(I4B) :: ip
615  class(BndType), pointer :: packobj
616  !
617  ! -- Deallocate idm memory
618  call memorylist_remove(this%name, 'NAM', idm_context)
619  call memorylist_remove(component=this%name, context=idm_context)
620  !
621  ! -- Internal flow packages deallocate
622  call this%dis%dis_da()
623  call this%ic%ic_da()
624  call this%fmi%fmi_da()
625  call this%adv%adv_da()
626  call this%cnd%cnd_da()
627  call this%ssm%ssm_da()
628  call this%est%est_da()
629  call this%mvt%mvt_da()
630  call this%budget%budget_da()
631  call this%oc%oc_da()
632  call this%obs%obs_da()
633  call this%gwecommon%gweshared_dat_da()
634  !
635  ! -- Internal package objects
636  deallocate (this%dis)
637  deallocate (this%ic)
638  deallocate (this%fmi)
639  deallocate (this%adv)
640  deallocate (this%cnd)
641  deallocate (this%ssm)
642  deallocate (this%est)
643  deallocate (this%mvt)
644  deallocate (this%budget)
645  deallocate (this%oc)
646  deallocate (this%obs)
647  nullify (this%gwecommon)
648  !
649  ! -- Boundary packages
650  do ip = 1, this%bndlist%Count()
651  packobj => getbndfromlist(this%bndlist, ip)
652  call packobj%bnd_da()
653  deallocate (packobj)
654  end do
655  !
656  ! -- Scalars
657  call mem_deallocate(this%inest)
658  call mem_deallocate(this%incnd)
659  !
660  ! -- Parent class members
661  call this%TransportModelType%tsp_da()
662  !
663  ! -- NumericalModelType
664  call this%NumericalModelType%model_da()
665  !
666  ! -- Return
667  return
subroutine, public memorylist_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 148 of file gwe.f90.

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

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

◆ gwe_get_iasym()

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

Definition at line 695 of file gwe.f90.

696  class(GweModelType) :: this
697  ! -- local
698  integer(I4B) :: iasym
699  integer(I4B) :: ip
700  class(BndType), pointer :: packobj
701  !
702  ! -- Start by setting iasym to zero
703  iasym = 0
704  !
705  ! -- ADV
706  if (this%inadv > 0) then
707  if (this%adv%iasym /= 0) iasym = 1
708  end if
709  !
710  ! -- CND
711  if (this%incnd > 0) then
712  if (this%cnd%ixt3d /= 0) iasym = 1
713  end if
714  !
715  ! -- Check for any packages that introduce matrix asymmetry
716  do ip = 1, this%bndlist%Count()
717  packobj => getbndfromlist(this%bndlist, ip)
718  if (packobj%iasym /= 0) iasym = 1
719  end do
720  !
721  ! -- Return
722  return
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 235 of file gwe.f90.

236  ! -- dummy
237  class(GweModelType) :: this
238  class(MatrixBaseType), pointer :: matrix_sln !< global system matrix
239  ! -- local
240  class(BndType), pointer :: packobj
241  integer(I4B) :: ip
242  !
243  ! -- Find the position of each connection in the global ia, ja structure
244  ! and store them in idxglo.
245  call this%dis%dis_mc(this%moffset, this%idxglo, matrix_sln)
246  !
247  if (this%incnd > 0) call this%cnd%cnd_mc(this%moffset, matrix_sln)
248  !
249  ! -- Map any package connections
250  do ip = 1, this%bndlist%Count()
251  packobj => getbndfromlist(this%bndlist, ip)
252  call packobj%bnd_mc(this%moffset, matrix_sln)
253  end do
254  !
255  ! -- Return
256  return
Here is the call graph for this function:

◆ gwe_ot()

subroutine gwemodule::gwe_ot ( class(gwemodeltype this)

This subroutine calls the parent class output routine.

Definition at line 579 of file gwe.f90.

580  ! -- dummy
581  class(GweModelType) :: this
582  ! -- local
583  integer(I4B) :: icbcfl
584  integer(I4B) :: icbcun
585  ! -- formats
586  !
587  ! -- Initialize
588  icbcfl = 0
589  !
590  ! -- Because est belongs to gwe, call est_ot_flow directly (and not from parent)
591  if (this%oc%oc_save('BUDGET')) icbcfl = 1
592  icbcun = this%oc%oc_save_unit('BUDGET')
593  if (this%inest > 0) call this%est%est_ot_flow(icbcfl, icbcun)
594  !
595  ! -- Call parent class _ot routines.
596  call this%tsp_ot(this%inest)
597  !
598  ! -- Return
599  return

◆ gwe_rp()

subroutine gwemodule::gwe_rp ( class(gwemodeltype this)

This subroutine calls the attached packages' read and prepare routines

Definition at line 311 of file gwe.f90.

312  ! -- modules
313  use tdismodule, only: readnewdata
314  ! -- dummy
315  class(GweModelType) :: this
316  ! -- local
317  class(BndType), pointer :: packobj
318  integer(I4B) :: ip
319  !
320  ! -- In fmi, check for mvt and mvrbudobj consistency
321  call this%fmi%fmi_rp(this%inmvt)
322  if (this%inmvt > 0) call this%mvt%mvt_rp()
323  !
324  ! -- Check with TDIS on whether or not it is time to RP
325  if (.not. readnewdata) return
326  !
327  ! -- Read and prepare
328  if (this%inoc > 0) call this%oc%oc_rp()
329  if (this%inssm > 0) call this%ssm%ssm_rp()
330  do ip = 1, this%bndlist%Count()
331  packobj => getbndfromlist(this%bndlist, ip)
332  call packobj%bnd_rp()
333  call packobj%bnd_rp_obs()
334  end do
335  !
336  ! -- Return
337  return
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 757 of file gwe.f90.

759  ! -- modules
760  use constantsmodule, only: linelength
761  use simmodule, only: store_error
762  use gwectpmodule, only: ctp_create
763  use gweeslmodule, only: esl_create
764  use gwelkemodule, only: lke_create
765  use gwesfemodule, only: sfe_create
766  use gwemwemodule, only: mwe_create
767  use gweuzemodule, only: uze_create
768  use apimodule, only: api_create
769  ! -- dummy
770  class(GweModelType) :: this
771  character(len=*), intent(in) :: filtyp
772  character(len=LINELENGTH) :: errmsg
773  integer(I4B), intent(in) :: ipakid
774  integer(I4B), intent(in) :: ipaknum
775  character(len=*), intent(in) :: pakname
776  character(len=*), intent(in) :: mempath
777  integer(I4B), intent(in) :: inunit
778  integer(I4B), intent(in) :: iout
779  ! -- local
780  class(BndType), pointer :: packobj
781  class(BndType), pointer :: packobj2
782  integer(I4B) :: ip
783  !
784  ! -- This part creates the package object
785  select case (filtyp)
786  case ('CTP6')
787  call ctp_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
788  pakname, this%depvartype, mempath)
789  case ('ESL6')
790  call esl_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
791  pakname, this%gwecommon)
792  case ('LKE6')
793  call lke_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
794  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
795  this%depvartype, this%depvarunit, this%depvarunitabbrev)
796  case ('SFE6')
797  call sfe_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
798  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
799  this%depvartype, this%depvarunit, this%depvarunitabbrev)
800  case ('MWE6')
801  call mwe_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
802  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
803  this%depvartype, this%depvarunit, this%depvarunitabbrev)
804  case ('UZE6')
805  call uze_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
806  pakname, this%fmi, this%eqnsclfac, this%gwecommon, &
807  this%depvartype, this%depvarunit, this%depvarunitabbrev)
808  !case('API6')
809  ! call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
810  ! pakname)
811  case default
812  write (errmsg, *) 'Invalid package type: ', filtyp
813  call store_error(errmsg, terminate=.true.)
814  end select
815  !
816  ! -- Packages is the bndlist that is associated with the parent model
817  ! -- The following statement puts a pointer to this package in the ipakid
818  ! -- position of packages.
819  do ip = 1, this%bndlist%Count()
820  packobj2 => getbndfromlist(this%bndlist, ip)
821  if (packobj2%packName == pakname) then
822  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
823  'already exists: ', trim(pakname)
824  call store_error(errmsg, terminate=.true.)
825  end if
826  end do
827  call addbndtolist(this%bndlist, packobj)
828  !
829  ! -- Return
830  return
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:58
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