MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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_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 (this)
 Print and/or save model 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 734 of file gwt.f90.

735  ! -- modules
737  ! -- dummy
738  class(GwtModelType) :: this
739  character(len=*), intent(in) :: modelname
740  !
741  ! -- allocate parent class scalars
742  call this%allocate_tsp_scalars(modelname)
743  !
744  ! -- allocate additional members specific to GWT model type
745  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
746  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
747  !
748  this%inmst = 0
749  this%indsp = 0
750  !
751  ! -- Return
752  return

◆ castasgwtmodel()

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

Definition at line 840 of file gwt.f90.

841  class(*), pointer :: model !< The object to be cast
842  class(GwtModelType), pointer :: gwtmodel !< The GWT model
843  !
844  gwtmodel => null()
845  if (.not. associated(model)) return
846  select type (model)
847  type is (gwtmodeltype)
848  gwtmodel => model
849  end select
850  !
851  ! -- Return
852  return
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 857 of file gwt.f90.

859  ! -- modules
862  ! -- dummy
863  class(GwtModelType) :: this
864  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
865  type(CharacterStringType), dimension(:), contiguous, &
866  pointer, intent(inout) :: pkgtypes
867  type(CharacterStringType), dimension(:), contiguous, &
868  pointer, intent(inout) :: pkgnames
869  type(CharacterStringType), dimension(:), contiguous, &
870  pointer, intent(inout) :: mempaths
871  integer(I4B), dimension(:), contiguous, &
872  pointer, intent(inout) :: inunits
873  ! -- local
874  integer(I4B) :: ipakid, ipaknum
875  character(len=LENFTYPE) :: pkgtype, bndptype
876  character(len=LENPACKAGENAME) :: pkgname
877  character(len=LENMEMPATH) :: mempath
878  integer(I4B), pointer :: inunit
879  integer(I4B) :: n
880  !
881  if (allocated(bndpkgs)) then
882  !
883  ! -- create stress packages
884  ipakid = 1
885  bndptype = ''
886  do n = 1, size(bndpkgs)
887  !
888  pkgtype = pkgtypes(bndpkgs(n))
889  pkgname = pkgnames(bndpkgs(n))
890  mempath = mempaths(bndpkgs(n))
891  inunit => inunits(bndpkgs(n))
892  !
893  if (bndptype /= pkgtype) then
894  ipaknum = 1
895  bndptype = pkgtype
896  end if
897  !
898  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
899  inunit, this%iout)
900  ipakid = ipakid + 1
901  ipaknum = ipaknum + 1
902  end do
903  !
904  ! -- cleanup
905  deallocate (bndpkgs)
906  end if
907  !
908  ! -- Return
909  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_gwt_packages()

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

Definition at line 914 of file gwt.f90.

915  ! -- modules
922  use gwtmstmodule, only: mst_cr
923  use gwtdspmodule, only: dsp_cr
924  ! -- dummy
925  class(GwtModelType) :: this
926  integer(I4B), intent(in) :: indis
927  ! -- local
928  type(CharacterStringType), dimension(:), contiguous, &
929  pointer :: pkgtypes => null()
930  type(CharacterStringType), dimension(:), contiguous, &
931  pointer :: pkgnames => null()
932  type(CharacterStringType), dimension(:), contiguous, &
933  pointer :: mempaths => null()
934  integer(I4B), dimension(:), contiguous, &
935  pointer :: inunits => null()
936  character(len=LENMEMPATH) :: model_mempath
937  character(len=LENFTYPE) :: pkgtype
938  character(len=LENPACKAGENAME) :: pkgname
939  character(len=LENMEMPATH) :: mempath
940  integer(I4B), pointer :: inunit
941  integer(I4B), dimension(:), allocatable :: bndpkgs
942  integer(I4B) :: n
943  character(len=LENMEMPATH) :: mempathdsp = ''
944  !
945  ! -- set input memory paths, input/model and input/model/namfile
946  model_mempath = create_mem_path(component=this%name, context=idm_context)
947  !
948  ! -- set pointers to model path package info
949  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
950  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
951  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
952  call mem_setptr(inunits, 'INUNITS', model_mempath)
953  !
954  do n = 1, size(pkgtypes)
955  !
956  ! attributes for this input package
957  pkgtype = pkgtypes(n)
958  pkgname = pkgnames(n)
959  mempath = mempaths(n)
960  inunit => inunits(n)
961  !
962  ! -- create dis package as it is a prerequisite for other packages
963  select case (pkgtype)
964  case ('MST6')
965  this%inmst = inunit
966  case ('DSP6')
967  this%indsp = 1
968  mempathdsp = mempath
969  case ('CNC6', 'SRC6', 'LKT6', 'SFT6', &
970  'MWT6', 'UZT6', 'IST6', 'API6')
971  call expandarray(bndpkgs)
972  bndpkgs(size(bndpkgs)) = n
973  case default
974  ! TODO
975  end select
976  end do
977  !
978  ! -- Create packages that are tied directly to model
979  call mst_cr(this%mst, this%name, this%inmst, this%iout, this%fmi)
980  call dsp_cr(this%dsp, this%name, mempathdsp, this%indsp, this%iout, &
981  this%fmi)
982  !
983  ! -- Check to make sure that required ftype's have been specified
984  call this%ftype_check(indis, this%inmst)
985  !
986  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
987  !
988  ! -- Return
989  return
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:98
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 201 of file gwt.f90.

202  ! -- modules
203  use sparsemodule, only: sparsematrix
204  ! -- dummy
205  class(GwtModelType) :: this
206  type(sparsematrix), intent(inout) :: sparse
207  ! -- local
208  class(BndType), pointer :: packobj
209  integer(I4B) :: ip
210  !
211  ! -- Add the internal connections of this model to sparse
212  call this%dis%dis_ac(this%moffset, sparse)
213  if (this%indsp > 0) &
214  call this%dsp%dsp_ac(this%moffset, sparse)
215  !
216  ! -- Add any package connections
217  do ip = 1, this%bndlist%Count()
218  packobj => getbndfromlist(this%bndlist, ip)
219  call packobj%bnd_ac(this%moffset, sparse)
220  end do
221  !
222  ! -- Return
223  return
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 344 of file gwt.f90.

345  ! -- modules
347  ! -- dummy
348  class(GwtModelType) :: 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%indsp > 0) call this%dsp%dsp_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:

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

260  ! -- modules
261  use constantsmodule, only: dhnoflo
262  ! -- dummy
263  class(GwtModelType) :: this
264  ! -- locals
265  integer(I4B) :: ip
266  class(BndType), pointer :: packobj
267  !
268  ! -- Allocate and read modules attached to model
269  call this%fmi%fmi_ar(this%ibound)
270  if (this%inmvt > 0) call this%mvt%mvt_ar()
271  if (this%inic > 0) call this%ic%ic_ar(this%x)
272  if (this%inmst > 0) call this%mst%mst_ar(this%dis, this%ibound)
273  if (this%inadv > 0) call this%adv%adv_ar(this%dis, this%ibound)
274  if (this%indsp > 0) call this%dsp%dsp_ar(this%ibound, this%mst%thetam)
275  if (this%inssm > 0) call this%ssm%ssm_ar(this%dis, this%ibound, this%x)
276  if (this%inobs > 0) call this%obs%tsp_obs_ar(this%ic, this%x, this%flowja)
277  !
278  ! -- Set governing equation scale factor. Note that this scale factor
279  ! -- cannot be set arbitrarily. For solute transport, it must be set
280  ! -- to 1. Setting it to a different value will NOT automatically
281  ! -- scale all the terms of the governing equation correctly by that
282  ! -- value. This is because much of the coding in the associated
283  ! -- packages implicitly assumes the governing equation for solute
284  ! -- transport is scaled by 1. (effectively unscaled).
285  this%eqnsclfac = done
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:

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

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

680  ! -- modules
681  use constantsmodule, only: lenbudtxt
682  use tdismodule, only: delt
683  ! -- dummy
684  class(GwtModelType) :: this
685  real(DP), dimension(:, :), intent(in) :: budterm
686  character(len=LENBUDTXT), dimension(:), intent(in) :: budtxt
687  character(len=*), intent(in) :: rowlabel
688  !
689  call this%budget%addentry(budterm, delt, budtxt, rowlabel=rowlabel)
690  !
691  ! -- Return
692  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

◆ 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)
490  !
491  ! -- Return
492  return

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

401  ! -- modules
402  ! -- dummy
403  class(GwtModelType) :: 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:

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

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

98  ! -- modules
99  use listsmodule, only: basemodellist
105  use budgetmodule, only: budget_cr
106  ! -- dummy
107  character(len=*), intent(in) :: filename !< input file
108  integer(I4B), intent(in) :: id !< consecutive model number listed in mfsim.nam
109  character(len=*), intent(in) :: modelname !< name of the model
110  ! -- local
111  integer(I4B) :: indis
112  type(GwtModelType), pointer :: this
113  class(BaseModelType), pointer :: model
114  !
115  ! -- Allocate a new GWT Model (this)
116  allocate (this)
117  !
118  ! -- Set memory path before allocation in memory manager can be done
119  this%memoryPath = create_mem_path(modelname)
120  !
121  ! -- Allocate scalars and add model to basemodellist
122  call this%allocate_scalars(modelname)
123  !
124  ! -- set labels for transport model - needed by create_packages() below
125  call this%set_tsp_labels(this%macronym, dvt, dvu, dvua)
126  !
127  model => this
128  call addbasemodeltolist(basemodellist, model)
129  !
130  ! -- Call parent class routine
131  call this%tsp_cr(filename, id, modelname, 'GWT', indis)
132  !
133  ! -- create model packages
134  call this%create_packages(indis)
135  !
136  ! -- Return
137  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
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 memorylist_remove(this%name, 'NAM', idm_context)
625  call memorylist_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()
668  !
669  ! -- Return
670  return
subroutine, public memorylist_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)
194  !
195  ! -- Return
196  return
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_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 423 of file gwt.f90.

424  ! -- modules
425  ! -- dummy
426  class(GwtModelType) :: 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%inmst > 0) then
441  call this%mst%mst_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%indsp > 0) then
449  call this%dsp%dsp_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:

◆ gwt_get_iasym()

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

Definition at line 698 of file gwt.f90.

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

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

◆ gwt_ot()

subroutine gwtmodule::gwt_ot ( class(gwtmodeltype this)

Call the parent class output routine

Definition at line 585 of file gwt.f90.

586  ! -- dummy
587  class(GwtModelType) :: this
588  ! -- local
589  integer(I4B) :: icbcfl
590  integer(I4B) :: icbcun
591  !
592  !
593  ! -- Initialize
594  icbcfl = 0
595  !
596  ! -- Because mst belongs to gwt, call mst_ot_flow directly (and not from parent)
597  if (this%oc%oc_save('BUDGET')) icbcfl = 1
598  icbcun = this%oc%oc_save_unit('BUDGET')
599  if (this%inmst > 0) call this%mst%mst_ot_flow(icbcfl, icbcun)
600  !
601  ! -- Call parent class _ot routines.
602  call this%tsp_ot(this%inmst)
603  !
604  ! -- Return
605  return

◆ gwt_rp()

subroutine gwtmodule::gwt_rp ( class(gwtmodeltype this)

Call the read and prepare routines of the attached packages

Definition at line 311 of file gwt.f90.

312  ! -- modules
313  use tdismodule, only: readnewdata
314  ! -- dummy
315  class(GwtModelType) :: 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 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 759 of file gwt.f90.

761  ! -- modules
762  use constantsmodule, only: linelength
763  use simmodule, only: store_error
764  use gwtcncmodule, only: cnc_create
765  use gwtsrcmodule, only: src_create
766  use gwtistmodule, only: ist_create
767  use gwtlktmodule, only: lkt_create
768  use gwtsftmodule, only: sft_create
769  use gwtmwtmodule, only: mwt_create
770  use gwtuztmodule, only: uzt_create
771  use apimodule, only: api_create
772  ! -- dummy
773  class(GwtModelType) :: this
774  character(len=*), intent(in) :: filtyp
775  character(len=LINELENGTH) :: errmsg
776  integer(I4B), intent(in) :: ipakid
777  integer(I4B), intent(in) :: ipaknum
778  character(len=*), intent(in) :: pakname
779  character(len=*), intent(in) :: mempath
780  integer(I4B), intent(in) :: inunit
781  integer(I4B), intent(in) :: iout
782  ! -- local
783  class(BndType), pointer :: packobj
784  class(BndType), pointer :: packobj2
785  integer(I4B) :: ip
786  !
787  ! -- This part creates the package object
788  select case (filtyp)
789  case ('CNC6')
790  call cnc_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
791  pakname, dvt, mempath)
792  case ('SRC6')
793  call src_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
794  this%depvartype, pakname)
795  case ('LKT6')
796  call lkt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
797  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
798  this%depvarunit, this%depvarunitabbrev)
799  case ('SFT6')
800  call sft_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
801  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
802  this%depvarunit, this%depvarunitabbrev)
803  case ('MWT6')
804  call mwt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
805  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
806  this%depvarunit, this%depvarunitabbrev)
807  case ('UZT6')
808  call uzt_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
809  pakname, this%fmi, this%eqnsclfac, this%depvartype, &
810  this%depvarunit, this%depvarunitabbrev)
811  case ('IST6')
812  call ist_create(packobj, ipakid, ipaknum, inunit, iout, this%name, &
813  pakname, this%fmi, this%mst)
814  case ('API6')
815  call api_create(packobj, ipakid, ipaknum, inunit, iout, this%name, pakname)
816  case default
817  write (errmsg, *) 'Invalid package type: ', filtyp
818  call store_error(errmsg, terminate=.true.)
819  end select
820  !
821  ! -- Packages is the bndlist that is associated with the parent model
822  ! -- The following statement puts a pointer to this package in the ipakid
823  ! -- position of packages.
824  do ip = 1, this%bndlist%Count()
825  packobj2 => getbndfromlist(this%bndlist, ip)
826  if (packobj2%packName == pakname) then
827  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
828  'already exists: ', trim(pakname)
829  call store_error(errmsg, terminate=.true.)
830  end if
831  end do
832  call addbndtolist(this%bndlist, packobj)
833  !
834  ! -- Return
835  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 cnc_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, depvartype, mempath)
Create a new constant concentration or temperature package.
Definition: gwt-cnc.f90:58
– @ 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:109
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 31 of file gwt.f90.

31  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 32 of file gwt.f90.

32  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 33 of file gwt.f90.

33  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 72 of file gwt.f90.

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

◆ gwt_multipkg

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

Definition at line 85 of file gwt.f90.

85  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 71 of file gwt.f90.

71  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 84 of file gwt.f90.

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

◆ niunit_gwt

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

Definition at line 91 of file gwt.f90.

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