MODFLOW 6  version 6.8.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 725 of file gwt.f90.

726  ! -- modules
728  ! -- dummy
729  class(GwtModelType) :: this
730  character(len=*), intent(in) :: modelname
731  !
732  ! -- allocate parent class scalars
733  call this%allocate_tsp_scalars(modelname)
734  !
735  ! -- allocate additional members specific to GWT model type
736  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
737  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
738  !
739  this%inmst = 0
740  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 825 of file gwt.f90.

826  class(*), pointer :: model !< The object to be cast
827  class(GwtModelType), pointer :: gwtmodel !< The GWT model
828  !
829  gwtmodel => null()
830  if (.not. associated(model)) return
831  select type (model)
832  type is (gwtmodeltype)
833  gwtmodel => model
834  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 839 of file gwt.f90.

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

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

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

545  use constantsmodule, only: dzero
546  ! -- dummy
547  class(GwtModelType) :: 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%inmst > 0) call this%mst%mst_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
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 676 of file gwt.f90.

677  ! -- modules
678  use constantsmodule, only: lenbudtxt
679  use tdismodule, only: delt
680  ! -- dummy
681  class(GwtModelType) :: 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)
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 475 of file gwt.f90.

476  ! -- dummy
477  class(GwtModelType) :: this
478  integer(I4B), intent(in) :: innertot
479  integer(I4B), intent(in) :: kiter
480  integer(I4B), intent(in) :: iend
481  integer(I4B), intent(in) :: icnvgmod
482  character(len=LENPAKLOC), intent(inout) :: cpak
483  integer(I4B), intent(inout) :: ipak
484  real(DP), intent(inout) :: dpak
485  ! -- local
486  ! -- formats
487  !
488  ! -- If mover is on, then at least 2 outers required
489  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 406 of file gwt.f90.

407  ! -- modules
408  ! -- dummy
409  class(GwtModelType) :: this
410  integer(I4B), intent(in) :: kiter
411  ! -- local
412  class(BndType), pointer :: packobj
413  integer(I4B) :: ip
414  !
415  ! -- Call package cf routines
416  do ip = 1, this%bndlist%Count()
417  packobj => getbndfromlist(this%bndlist, ip)
418  call packobj%bnd_cf()
419  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 496 of file gwt.f90.

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

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

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

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

693  class(GwtModelType) :: this
694  ! -- local
695  integer(I4B) :: iasym
696  integer(I4B) :: ip
697  class(BndType), pointer :: packobj
698  !
699  ! -- Start by setting iasym to zero
700  iasym = 0
701  !
702  ! -- ADV
703  if (this%inadv > 0) then
704  if (this%adv%iasym /= 0) iasym = 1
705  end if
706  !
707  ! -- DSP
708  if (this%indsp > 0) then
709  if (this%dsp%ixt3d /= 0) iasym = 1
710  end if
711  !
712  ! -- Check for any packages that introduce matrix asymmetry
713  do ip = 1, this%bndlist%Count()
714  packobj => getbndfromlist(this%bndlist, ip)
715  if (packobj%iasym /= 0) iasym = 1
716  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 594 of file gwt.f90.

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

◆ 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 576 of file gwt.f90.

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

◆ 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_log()
322  call packobj%bnd_rp_obs()
323  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 747 of file gwt.f90.

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