MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
gwtmodule Module Reference

Data Types

type  gwtmodeltype
 

Functions/Subroutines

subroutine, public gwt_cr (filename, id, modelname)
 Create a new groundwater transport model object. More...
 
subroutine gwt_df (this)
 Define packages of the GWT model. More...
 
subroutine gwt_ac (this, sparse)
 Add the internal connections of this model to the sparse matrix. More...
 
subroutine gwt_mc (this, matrix_sln)
 Map the positions of the GWT model connections in the numerical solution coefficient matrix. More...
 
subroutine gwt_ar (this)
 GWT Model Allocate and Read. More...
 
subroutine gwt_rp (this)
 GWT Model Read and Prepare. More...
 
subroutine gwt_dt (this)
 GWT Model time step size. More...
 
subroutine gwt_ad (this)
 GWT Model Time Step Advance. More...
 
subroutine gwt_cf (this, kiter)
 GWT Model calculate coefficients. More...
 
subroutine gwt_fc (this, kiter, matrix_sln, inwtflag)
 GWT Model fill coefficients. More...
 
subroutine gwt_cc (this, innertot, kiter, iend, icnvgmod, cpak, ipak, dpak)
 GWT Model Final Convergence Check. More...
 
subroutine gwt_cq (this, icnvg, isuppress_output)
 GWT Model calculate flow. More...
 
subroutine gwt_bd (this, icnvg, isuppress_output)
 GWT Model Budget. More...
 
subroutine gwt_ot_flow (this, icbcfl, ibudfl, icbcun)
 GWT model output routine. More...
 
subroutine gwt_ot_dv (this, idvsave, idvprint, ipflag)
 GWT model dependent variable output. More...
 
subroutine gwt_da (this)
 Deallocate. More...
 
subroutine gwt_bdentry (this, budterm, budtxt, rowlabel)
 GroundWater Transport Model Budget Entry. More...
 
integer(i4b) function gwt_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(gwtmodeltype) function, pointer, public castasgwtmodel (model)
 Cast to GwtModelType. More...
 
subroutine create_bndpkgs (this, bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
 Source package info and begin to process. More...
 
subroutine create_gwt_packages (this, indis)
 Source package info and begin to process. More...
 

Variables

character(len=lenvarname), parameter dvt = 'CONCENTRATION '
 dependent variable type, varies based on model type More...
 
character(len=lenvarname), parameter dvu = 'MASS '
 dependent variable unit of measure, either "mass" or "energy" More...
 
character(len=lenvarname), parameter dvua = 'M '
 abbreviation of the dependent variable unit of measure, either "M" or "E" More...
 
integer(i4b), parameter, public gwt_nbasepkg = 50
 GWT base package array descriptors. More...
 
character(len=lenpackagetype), dimension(gwt_nbasepkg), public gwt_basepkg
 
integer(i4b), parameter, public gwt_nmultipkg = 50
 GWT multi package array descriptors. More...
 
character(len=lenpackagetype), dimension(gwt_nmultipkg), public gwt_multipkg
 
integer(i4b), parameter niunit_gwt = GWT_NBASEPKG + GWT_NMULTIPKG
 

Function/Subroutine Documentation

◆ allocate_scalars()

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

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

Definition at line 724 of file gwt.f90.

725  ! -- modules
727  ! -- dummy
728  class(GwtModelType) :: this
729  character(len=*), intent(in) :: modelname
730  !
731  ! -- allocate parent class scalars
732  call this%allocate_tsp_scalars(modelname)
733  !
734  ! -- allocate additional members specific to GWT model type
735  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
736  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
737  !
738  this%inmst = 0
739  this%indsp = 0

◆ castasgwtmodel()

class(gwtmodeltype) function, pointer, public gwtmodule::castasgwtmodel ( class(*), pointer  model)
Parameters
modelThe object to be cast
Returns
The GWT model

Definition at line 823 of file gwt.f90.

824  class(*), pointer :: model !< The object to be cast
825  class(GwtModelType), pointer :: gwtmodel !< The GWT model
826  !
827  gwtmodel => null()
828  if (.not. associated(model)) return
829  select type (model)
830  type is (gwtmodeltype)
831  gwtmodel => model
832  end select
Here is the caller graph for this function:

◆ create_bndpkgs()

subroutine gwtmodule::create_bndpkgs ( class(gwtmodeltype 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 837 of file gwt.f90.

839  ! -- modules
842  ! -- dummy
843  class(GwtModelType) :: this
844  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
845  type(CharacterStringType), dimension(:), contiguous, &
846  pointer, intent(inout) :: pkgtypes
847  type(CharacterStringType), dimension(:), contiguous, &
848  pointer, intent(inout) :: pkgnames
849  type(CharacterStringType), dimension(:), contiguous, &
850  pointer, intent(inout) :: mempaths
851  integer(I4B), dimension(:), contiguous, &
852  pointer, intent(inout) :: inunits
853  ! -- local
854  integer(I4B) :: ipakid, ipaknum
855  character(len=LENFTYPE) :: pkgtype, bndptype
856  character(len=LENPACKAGENAME) :: pkgname
857  character(len=LENMEMPATH) :: mempath
858  integer(I4B), pointer :: inunit
859  integer(I4B) :: n
860  !
861  if (allocated(bndpkgs)) then
862  !
863  ! -- create stress packages
864  ipakid = 1
865  bndptype = ''
866  do n = 1, size(bndpkgs)
867  !
868  pkgtype = pkgtypes(bndpkgs(n))
869  pkgname = pkgnames(bndpkgs(n))
870  mempath = mempaths(bndpkgs(n))
871  inunit => inunits(bndpkgs(n))
872  !
873  if (bndptype /= pkgtype) then
874  ipaknum = 1
875  bndptype = pkgtype
876  end if
877  !
878  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
879  inunit, this%iout)
880  ipakid = ipakid + 1
881  ipaknum = ipaknum + 1
882  end do
883  !
884  ! -- cleanup
885  deallocate (bndpkgs)
886  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_gwt_packages()

subroutine gwtmodule::create_gwt_packages ( class(gwtmodeltype this,
integer(i4b), intent(in)  indis 
)

Definition at line 891 of file gwt.f90.

892  ! -- modules
899  use gwtmstmodule, only: mst_cr
900  use gwtdspmodule, only: dsp_cr
901  ! -- dummy
902  class(GwtModelType) :: this
903  integer(I4B), intent(in) :: indis
904  ! -- local
905  type(CharacterStringType), dimension(:), contiguous, &
906  pointer :: pkgtypes => null()
907  type(CharacterStringType), dimension(:), contiguous, &
908  pointer :: pkgnames => null()
909  type(CharacterStringType), dimension(:), contiguous, &
910  pointer :: mempaths => null()
911  integer(I4B), dimension(:), contiguous, &
912  pointer :: inunits => null()
913  character(len=LENMEMPATH) :: model_mempath
914  character(len=LENFTYPE) :: pkgtype
915  character(len=LENPACKAGENAME) :: pkgname
916  character(len=LENMEMPATH) :: mempath
917  integer(I4B), pointer :: inunit
918  integer(I4B), dimension(:), allocatable :: bndpkgs
919  integer(I4B) :: n
920  character(len=LENMEMPATH) :: mempathdsp = ''
921  !
922  ! -- set input memory paths, input/model and input/model/namfile
923  model_mempath = create_mem_path(component=this%name, context=idm_context)
924  !
925  ! -- set pointers to model path package info
926  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
927  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
928  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
929  call mem_setptr(inunits, 'INUNITS', model_mempath)
930  !
931  do n = 1, size(pkgtypes)
932  !
933  ! attributes for this input package
934  pkgtype = pkgtypes(n)
935  pkgname = pkgnames(n)
936  mempath = mempaths(n)
937  inunit => inunits(n)
938  !
939  ! -- create dis package as it is a prerequisite for other packages
940  select case (pkgtype)
941  case ('MST6')
942  this%inmst = inunit
943  case ('DSP6')
944  this%indsp = 1
945  mempathdsp = mempath
946  case ('CNC6', 'SRC6', 'LKT6', 'SFT6', &
947  'MWT6', 'UZT6', 'IST6', 'API6')
948  call expandarray(bndpkgs)
949  bndpkgs(size(bndpkgs)) = n
950  case default
951  ! TODO
952  end select
953  end do
954  !
955  ! -- Create packages that are tied directly to model
956  call mst_cr(this%mst, this%name, this%inmst, this%iout, this%fmi)
957  call dsp_cr(this%dsp, this%name, mempathdsp, this%indsp, this%iout, &
958  this%fmi)
959  !
960  ! -- Check to make sure that required ftype's have been specified
961  call this%ftype_check(indis, this%inmst)
962  !
963  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
subroutine, public dsp_cr(dspobj, name_model, input_mempath, inunit, iout, fmi)
Create a DSP object.
Definition: gwt-dsp.f90:78
– @ brief Mobile Storage and Transfer (MST) Module
Definition: gwt-mst.f90:10
subroutine, public mst_cr(mstobj, name_model, inunit, iout, fmi)
@ brief Create a new package object
Definition: gwt-mst.f90:103
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:

◆ gwt_ac()

subroutine gwtmodule::gwt_ac ( class(gwtmodeltype this,
type(sparsematrix), intent(inout)  sparse 
)

Definition at line 198 of file gwt.f90.

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

◆ gwt_ad()

subroutine gwtmodule::gwt_ad ( class(gwtmodeltype this)

Call the advance subroutines of the attached packages

Definition at line 352 of file gwt.f90.

353  ! -- modules
355  ! -- dummy
356  class(GwtModelType) :: this
357  class(BndType), pointer :: packobj
358  ! -- local
359  integer(I4B) :: irestore
360  integer(I4B) :: ip, n
361  !
362  ! -- Reset state variable
363  irestore = 0
364  if (ifailedstepretry > 0) irestore = 1
365  if (irestore == 0) then
366  !
367  ! -- copy x into xold
368  do n = 1, this%dis%nodes
369  if (this%ibound(n) == 0) then
370  this%xold(n) = dzero
371  else
372  this%xold(n) = this%x(n)
373  end if
374  end do
375  else
376  !
377  ! -- copy xold into x if this time step is a redo
378  do n = 1, this%dis%nodes
379  this%x(n) = this%xold(n)
380  end do
381  end if
382  !
383  ! -- Advance fmi
384  call this%fmi%fmi_ad(this%x)
385  !
386  ! -- Advance
387  if (this%indsp > 0) call this%dsp%dsp_ad()
388  if (this%inssm > 0) call this%ssm%ssm_ad()
389  do ip = 1, this%bndlist%Count()
390  packobj => getbndfromlist(this%bndlist, ip)
391  call packobj%bnd_ad()
392  if (isimcheck > 0) then
393  call packobj%bnd_ck()
394  end if
395  end do
396  !
397  ! -- Push simulated values to preceding time/subtime step
398  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:

◆ gwt_ar()

subroutine gwtmodule::gwt_ar ( class(gwtmodeltype 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 250 of file gwt.f90.

251  ! -- modules
252  use constantsmodule, only: dhnoflo
253  ! -- dummy
254  class(GwtModelType) :: this
255  ! -- locals
256  integer(I4B) :: ip
257  class(BndType), pointer :: packobj
258  !
259  ! -- Allocate and read modules attached to model
260  call this%fmi%fmi_ar(this%ibound)
261  if (this%inmvt > 0) call this%mvt%mvt_ar()
262  if (this%inic > 0) call this%ic%ic_ar(this%x)
263  if (this%inmst > 0) call this%mst%mst_ar(this%dis, this%ibound)
264  if (this%inadv > 0) call this%adv%adv_ar(this%dis, this%ibound)
265  if (this%indsp > 0) call this%dsp%dsp_ar(this%ibound, this%mst%thetam)
266  if (this%inssm > 0) call this%ssm%ssm_ar(this%dis, this%ibound, this%x)
267  if (this%inobs > 0) call this%obs%tsp_obs_ar(this%ic, this%x, this%flowja)
268  !
269  ! -- Set governing equation scale factor. Note that this scale factor
270  ! -- cannot be set arbitrarily. For solute transport, it must be set
271  ! -- to 1. Setting it to a different value will NOT automatically
272  ! -- scale all the terms of the governing equation correctly by that
273  ! -- value. This is because much of the coding in the associated
274  ! -- packages implicitly assumes the governing equation for solute
275  ! -- transport is scaled by 1. (effectively unscaled).
276  this%eqnsclfac = done
277  !
278  ! -- Call dis_ar to write binary grid file
279  !call this%dis%dis_ar(this%npf%icelltype)
280  !
281  ! -- set up output control
282  call this%oc%oc_ar(this%x, this%dis, dhnoflo, this%depvartype)
283  call this%budget%set_ibudcsv(this%oc%ibudcsv)
284  !
285  ! -- Package input files now open, so allocate and read
286  do ip = 1, this%bndlist%Count()
287  packobj => getbndfromlist(this%bndlist, ip)
288  call packobj%set_pointers(this%dis%nodes, this%ibound, this%x, &
289  this%xold, this%flowja)
290  ! -- Read and allocate package
291  call packobj%bnd_ar()
292  end do
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
Here is the call graph for this function:

◆ gwt_bd()

subroutine gwtmodule::gwt_bd ( class(gwtmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)

This subroutine: (1) calculates intercell flows (flowja) (2) calculates package contributions to the model budget

Definition at line 543 of file gwt.f90.

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

◆ gwt_bdentry()

subroutine gwtmodule::gwt_bdentry ( class(gwtmodeltype 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 gwt model object so that the exchange object could add its contributions.

Definition at line 675 of file gwt.f90.

676  ! -- modules
677  use constantsmodule, only: lenbudtxt
678  use tdismodule, only: delt
679  ! -- dummy
680  class(GwtModelType) :: this
681  real(DP), dimension(:, :), intent(in) :: budterm
682  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
683  character(len=*), intent(in) :: rowlabel
684  !
685  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

◆ gwt_cc()

subroutine gwtmodule::gwt_cc ( class(gwtmodeltype 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, call the MVR convergence check subroutines to force at least 2 outer iterations. The other advanced transport packages are solved in the matrix equations directly which means the solver is completing the necessary checks thereby eliminating need to call package cc routines. That is, no need to loop over active packages and run: call packobjbnd_cc(iend, icnvg, hclose, rclose)

Definition at line 474 of file gwt.f90.

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

◆ gwt_cf()

subroutine gwtmodule::gwt_cf ( class(gwtmodeltype this,
integer(i4b), intent(in)  kiter 
)

Call the calculate coefficients subroutines of the attached packages

Definition at line 405 of file gwt.f90.

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

◆ gwt_cq()

subroutine gwtmodule::gwt_cq ( class(gwtmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)
private

Call the intercell flows (flow ja) subroutine

Definition at line 495 of file gwt.f90.

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

◆ gwt_cr()

subroutine, public gwtmodule::gwt_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 100 of file gwt.f90.

101  ! -- modules
102  use listsmodule, only: basemodellist
108  use budgetmodule, only: budget_cr
109  ! -- dummy
110  character(len=*), intent(in) :: filename !< input file
111  integer(I4B), intent(in) :: id !< consecutive model number listed in mfsim.nam
112  character(len=*), intent(in) :: modelname !< name of the model
113  ! -- local
114  integer(I4B) :: indis
115  type(GwtModelType), pointer :: this
116  class(BaseModelType), pointer :: model
117  !
118  ! -- Allocate a new GWT Model (this)
119  allocate (this)
120  !
121  ! -- Set memory path before allocation in memory manager can be done
122  this%memoryPath = create_mem_path(modelname)
123  !
124  ! -- Allocate scalars and add model to basemodellist
125  call this%allocate_scalars(modelname)
126  !
127  ! -- set labels for transport model - needed by create_packages() below
128  call this%set_tsp_labels(this%macronym, dvt, dvu, dvua)
129  !
130  model => this
131  call addbasemodeltolist(basemodellist, model)
132  !
133  ! -- Call parent class routine
134  call this%tsp_cr(filename, id, modelname, 'GWT', indis)
135  !
136  ! -- create model packages
137  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
type(listtype), public basemodellist
Definition: mf6lists.f90:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gwt_da()

subroutine gwtmodule::gwt_da ( class(gwtmodeltype this)
private

Deallocate memory at conclusion of model run

Definition at line 611 of file gwt.f90.

612  ! -- modules
616  ! -- dummy
617  class(GwtModelType) :: this
618  ! -- local
619  integer(I4B) :: ip
620  class(BndType), pointer :: packobj
621  !
622  ! -- Deallocate idm memory
623  call memorystore_remove(this%name, 'NAM', idm_context)
624  call memorystore_remove(component=this%name, context=idm_context)
625  !
626  ! -- Internal flow packages deallocate
627  call this%dis%dis_da()
628  call this%ic%ic_da()
629  call this%fmi%fmi_da()
630  call this%adv%adv_da()
631  call this%dsp%dsp_da()
632  call this%ssm%ssm_da()
633  call this%mst%mst_da()
634  call this%mvt%mvt_da()
635  call this%budget%budget_da()
636  call this%oc%oc_da()
637  call this%obs%obs_da()
638  !
639  ! -- Internal package objects
640  deallocate (this%dis)
641  deallocate (this%ic)
642  deallocate (this%dsp)
643  deallocate (this%ssm)
644  deallocate (this%mst)
645  deallocate (this%adv)
646  deallocate (this%mvt)
647  deallocate (this%budget)
648  deallocate (this%oc)
649  deallocate (this%obs)
650  !
651  ! -- Boundary packages
652  do ip = 1, this%bndlist%Count()
653  packobj => getbndfromlist(this%bndlist, ip)
654  call packobj%bnd_da()
655  deallocate (packobj)
656  end do
657  !
658  ! -- Scalars
659  call mem_deallocate(this%indsp)
660  call mem_deallocate(this%inmst)
661  !
662  ! -- Parent class members
663  call this%TransportModelType%tsp_da()
664  !
665  ! -- NumericalModelType
666  call this%NumericalModelType%model_da()
subroutine, public memorystore_remove(component, subcomponent, context)
Here is the call graph for this function:

◆ gwt_df()

subroutine gwtmodule::gwt_df ( class(gwtmodeltype this)

This subroutine defines a gwt model type. Steps include: (1) call df routines for each package (2) set variables and pointers

Definition at line 146 of file gwt.f90.

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

◆ gwt_dt()

subroutine gwtmodule::gwt_dt ( class(gwtmodeltype 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 331 of file gwt.f90.

332  use tdismodule, only: kstp, kper
334  ! dummy
335  class(GwtModelType) :: this
336  ! local
337  real(DP) :: dtmax
338  character(len=LINELENGTH) :: msg
339  dtmax = dnodata
340 
341  ! advection package courant stability
342  call this%adv%adv_dt(dtmax, msg, this%mst%thetam)
343  if (msg /= '') then
344  call ats_submit_delt(kstp, kper, dtmax, msg)
345  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:

◆ gwt_fc()

subroutine gwtmodule::gwt_fc ( class(gwtmodeltype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln,
integer(i4b), intent(in)  inwtflag 
)
private

Call the fill coefficients subroutines attached packages

Definition at line 425 of file gwt.f90.

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

◆ gwt_get_iasym()

integer(i4b) function gwtmodule::gwt_get_iasym ( class(gwtmodeltype this)

Definition at line 691 of file gwt.f90.

692  class(GwtModelType) :: this
693  ! -- local
694  integer(I4B) :: iasym
695  integer(I4B) :: ip
696  class(BndType), pointer :: packobj
697  !
698  ! -- Start by setting iasym to zero
699  iasym = 0
700  !
701  ! -- ADV
702  if (this%inadv > 0) then
703  if (this%adv%iasym /= 0) iasym = 1
704  end if
705  !
706  ! -- DSP
707  if (this%indsp > 0) then
708  if (this%dsp%ixt3d /= 0) iasym = 1
709  end if
710  !
711  ! -- Check for any packages that introduce matrix asymmetry
712  do ip = 1, this%bndlist%Count()
713  packobj => getbndfromlist(this%bndlist, ip)
714  if (packobj%iasym /= 0) iasym = 1
715  end do
Here is the call graph for this function:

◆ gwt_mc()

subroutine gwtmodule::gwt_mc ( class(gwtmodeltype this,
class(matrixbasetype), pointer  matrix_sln 
)
Parameters
matrix_slnglobal system matrix

Definition at line 223 of file gwt.f90.

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

◆ gwt_ot_dv()

subroutine gwtmodule::gwt_ot_dv ( class(gwtmodeltype this,
integer(i4b), intent(in)  idvsave,
integer(i4b), intent(in)  idvprint,
integer(i4b), intent(inout)  ipflag 
)
private

Definition at line 593 of file gwt.f90.

594  class(GwtModelType) :: this
595  integer(I4B), intent(in) :: idvsave
596  integer(I4B), intent(in) :: idvprint
597  integer(I4B), intent(inout) :: ipflag
598 
599  ! call GWT dependent variable output routines
600  if (this%inmst > 0) call this%mst%mst_ot_dv(idvsave)
601 
602  ! call general transport model dependent variable output routines
603  call this%TransportModelType%tsp_ot_dv(idvsave, idvprint, ipflag)
604 

◆ gwt_ot_flow()

subroutine gwtmodule::gwt_ot_flow ( class(gwtmodeltype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(in)  icbcun 
)

Save and print flows

Definition at line 575 of file gwt.f90.

576  ! dummy
577  class(GwtModelType) :: this
578  integer(I4B), intent(in) :: icbcfl
579  integer(I4B), intent(in) :: ibudfl
580  integer(I4B), intent(in) :: icbcun
581  ! local
582 
583  ! call GWT flow output routines
584  if (this%inmst > 0) call this%mst%mst_ot_flow(icbcfl, icbcun)
585 
586  ! call general transport model flow output routines
587  call this%TransportModelType%tsp_ot_flow(icbcfl, ibudfl, icbcun)
588 

◆ gwt_rp()

subroutine gwtmodule::gwt_rp ( class(gwtmodeltype this)

Call the read and prepare routines of the attached packages

Definition at line 299 of file gwt.f90.

300  ! -- modules
301  use tdismodule, only: readnewdata
302  ! -- dummy
303  class(GwtModelType) :: this
304  ! -- local
305  class(BndType), pointer :: packobj
306  integer(I4B) :: ip
307  !
308  ! -- In fmi, check for mvt and mvrbudobj consistency
309  call this%fmi%fmi_rp(this%inmvt)
310  if (this%inmvt > 0) call this%mvt%mvt_rp()
311  !
312  ! -- Check with TDIS on whether or not it is time to RP
313  if (.not. readnewdata) return
314  !
315  ! -- Read and prepare
316  if (this%inoc > 0) call this%oc%oc_rp()
317  if (this%inssm > 0) call this%ssm%ssm_rp()
318  do ip = 1, this%bndlist%Count()
319  packobj => getbndfromlist(this%bndlist, ip)
320  call packobj%bnd_rp()
321  call packobj%bnd_rp_obs()
322  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 gwtmodule::package_create ( class(gwtmodeltype 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 
)

Call the package create routines for packages activated by the user.

Definition at line 746 of file gwt.f90.

748  ! -- modules
749  use simmodule, only: store_error
750  use gwtcncmodule, only: cnc_create
751  use gwtsrcmodule, only: src_create
752  use gwtistmodule, only: ist_create
753  use gwtlktmodule, only: lkt_create
754  use gwtsftmodule, only: sft_create
755  use gwtmwtmodule, only: mwt_create
756  use gwtuztmodule, only: uzt_create
757  use apimodule, only: api_create
758  ! -- dummy
759  class(GwtModelType) :: this
760  character(len=*), intent(in) :: filtyp
761  character(len=LINELENGTH) :: errmsg
762  integer(I4B), intent(in) :: ipakid
763  integer(I4B), intent(in) :: ipaknum
764  character(len=*), intent(in) :: pakname
765  character(len=*), intent(in) :: mempath
766  integer(I4B), intent(in) :: inunit
767  integer(I4B), intent(in) :: iout
768  ! -- local
769  class(BndType), pointer :: packobj
770  class(BndType), pointer :: packobj2
771  integer(I4B) :: ip
772  !
773  ! -- This part creates the package object
774  select case (filtyp)
775  case ('CNC6')
776  call cnc_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
777  pakname, dvt, mempath)
778  case ('SRC6')
779  call src_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
780  this%depvartype, pakname)
781  case ('LKT6')
782  call lkt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
783  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
784  this%depvarunit, this%depvarunitabbrev)
785  case ('SFT6')
786  call sft_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
787  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
788  this%depvarunit, this%depvarunitabbrev)
789  case ('MWT6')
790  call mwt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
791  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
792  this%depvarunit, this%depvarunitabbrev)
793  case ('UZT6')
794  call uzt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
795  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
796  this%depvarunit, this%depvarunitabbrev)
797  case ('IST6')
798  call ist_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
799  pakname, this%fmi, this%mst)
800  case ('API6')
801  call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
802  case default
803  write (errmsg, *) 'Invalid package type: ', filtyp
804  call store_error(errmsg, terminate=.true.)
805  end select
806  !
807  ! -- Packages is the bndlist that is associated with the parent model
808  ! -- The following statement puts a pointer to this package in the ipakid
809  ! -- position of packages.
810  do ip = 1, this%bndlist%Count()
811  packobj2 => getbndfromlist(this%bndlist, ip)
812  if (packobj2%packName == pakname) then
813  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
814  'already exists: ', trim(pakname)
815  call store_error(errmsg, terminate=.true.)
816  end if
817  end do
818  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 cnc_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype, mempath)
Create a new constant concentration or temperature package.
Definition: gwt-cnc.f90:57
– @ brief Immobile Storage and Transfer (IST) Module
Definition: gwt-ist.f90:14
subroutine, public ist_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, mst)
@ brief Create a new package object
Definition: gwt-ist.f90:117
subroutine, public lkt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, dvt, dvu, dvua)
Create a new lkt package.
Definition: gwt-lkt.f90:99
subroutine, public mwt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, dvt, dvu, dvua)
Create new MWT package.
Definition: gwt-mwt.f90:92
subroutine, public sft_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, dvt, dvu, dvua)
Create a new sft package.
Definition: gwt-sft.f90:96
subroutine, public src_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype)
Create a source loading package.
Definition: gwt-src.f90:47
subroutine, public uzt_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi, eqnsclfac, dvt, dvu, dvua)
Create a new UZT package.
Definition: gwt-uzt.f90:85
Here is the call graph for this function:

Variable Documentation

◆ dvt

character(len=lenvarname), parameter gwtmodule::dvt = 'CONCENTRATION '
private

Definition at line 32 of file gwt.f90.

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

◆ dvu

character(len=lenvarname), parameter gwtmodule::dvu = 'MASS '
private

Definition at line 33 of file gwt.f90.

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

◆ dvua

character(len=lenvarname), parameter gwtmodule::dvua = 'M '
private

Definition at line 34 of file gwt.f90.

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

◆ gwt_basepkg

character(len=lenpackagetype), dimension(gwt_nbasepkg), public gwtmodule::gwt_basepkg

Definition at line 75 of file gwt.f90.

75  character(len=LENPACKAGETYPE), dimension(GWT_NBASEPKG) :: GWT_BASEPKG

◆ gwt_multipkg

character(len=lenpackagetype), dimension(gwt_nmultipkg), public gwtmodule::gwt_multipkg

Definition at line 88 of file gwt.f90.

88  character(len=LENPACKAGETYPE), dimension(GWT_NMULTIPKG) :: GWT_MULTIPKG

◆ gwt_nbasepkg

integer(i4b), parameter, public gwtmodule::gwt_nbasepkg = 50

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

Definition at line 74 of file gwt.f90.

74  integer(I4B), parameter :: GWT_NBASEPKG = 50

◆ gwt_nmultipkg

integer(i4b), parameter, public gwtmodule::gwt_nmultipkg = 50

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

Definition at line 87 of file gwt.f90.

87  integer(I4B), parameter :: GWT_NMULTIPKG = 50

◆ niunit_gwt

integer(i4b), parameter gwtmodule::niunit_gwt = GWT_NBASEPKG + GWT_NMULTIPKG
private

Definition at line 94 of file gwt.f90.

94  integer(I4B), parameter :: NIUNIT_GWT = gwt_nbasepkg + gwt_nmultipkg