MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
prtmodule Module Reference

Data Types

type  prtmodeltype
 Particle tracking (PRT) model. More...
 

Functions/Subroutines

subroutine, public prt_cr (filename, id, modelname)
 Create a new particle tracking model object. More...
 
subroutine prt_df (this)
 Define packages. More...
 
subroutine prt_ar (this)
 Allocate and read. More...
 
subroutine prt_rp (this)
 Read and prepare (calls package read and prepare routines) More...
 
subroutine prt_ad (this)
 Time step advance (calls package advance subroutines) More...
 
subroutine prt_cq (this, icnvg, isuppress_output)
 Calculate intercell flow (flowja) More...
 
subroutine prt_cq_sto (this)
 Calculate particle mass storage. More...
 
subroutine prt_bd (this, icnvg, isuppress_output)
 Calculate flows and budget. More...
 
subroutine prt_ot (this)
 Print and/or save model output. More...
 
subroutine prt_ot_flow (this, icbcfl, ibudfl, icbcun)
 Save flows. More...
 
subroutine prt_ot_saveflow (this, nja, flowja, icbcfl, icbcun)
 Save intercell flows. More...
 
subroutine prt_ot_printflow (this, ibudfl, flowja)
 Print intercell flows. More...
 
subroutine prt_ot_dv (this, idvsave, idvprint, ipflag)
 Print dependent variables. More...
 
subroutine prt_ot_bdsummary (this, ibudfl, ipflag)
 Print budget summary. More...
 
subroutine prt_da (this)
 Deallocate. More...
 
subroutine allocate_scalars (this, modelname)
 Allocate memory for scalars. More...
 
subroutine allocate_arrays (this)
 Allocate arrays. More...
 
subroutine package_create (this, filtyp, ipakid, ipaknum, pakname, mempath, inunit, iout)
 Create boundary condition packages for this model. More...
 
subroutine ftype_check (this, indis)
 Check to make sure required input files have been specified. More...
 
subroutine prt_solve (this)
 Solve the model. More...
 
subroutine create_bndpkgs (this, bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
 Source package info and begin to process. More...
 
subroutine create_packages (this)
 Source package info and begin to process. More...
 
subroutine log_namfile_options (this, found)
 Write model namfile options to list file. More...
 

Variables

integer(i4b), parameter nbditems = 1
 
character(len=lenbudtxt), dimension(nbditemsbudtxt
 
integer(i4b), parameter, public prt_nbasepkg = 50
 PRT base package array descriptors. More...
 
character(len=lenpackagetype), dimension(prt_nbasepkg), public prt_basepkg
 
integer(i4b), parameter, public prt_nmultipkg = 50
 PRT multi package array descriptors. More...
 
character(len=lenpackagetype), dimension(prt_nmultipkg), public prt_multipkg
 
integer(i4b), parameter niunit_prt = PRT_NBASEPKG + PRT_NMULTIPKG
 

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine prtmodule::allocate_arrays ( class(prtmodeltype this)
private

Definition at line 810 of file prt.f90.

812  class(PrtModelType) :: this
813  integer(I4B) :: n
814 
815  ! Allocate arrays in parent type
816  this%nja = this%dis%nja
817  call this%NumericalModelType%allocate_arrays()
818 
819  ! Allocate and initialize arrays
820  call mem_allocate(this%masssto, this%dis%nodes, &
821  'MASSSTO', this%memoryPath)
822  call mem_allocate(this%massstoold, this%dis%nodes, &
823  'MASSSTOOLD', this%memoryPath)
824  call mem_allocate(this%ratesto, this%dis%nodes, &
825  'RATESTO', this%memoryPath)
826  ! explicit model, so these must be manually allocated
827  call mem_allocate(this%x, this%dis%nodes, 'X', this%memoryPath)
828  call mem_allocate(this%rhs, this%dis%nodes, 'RHS', this%memoryPath)
829  call mem_allocate(this%ibound, this%dis%nodes, 'IBOUND', this%memoryPath)
830  do n = 1, this%dis%nodes
831  this%masssto(n) = dzero
832  this%massstoold(n) = dzero
833  this%ratesto(n) = dzero
834  this%x(n) = dzero
835  this%rhs(n) = dzero
836  this%ibound(n) = 1
837  end do

◆ allocate_scalars()

subroutine prtmodule::allocate_scalars ( class(prtmodeltype this,
character(len=*), intent(in)  modelname 
)

Definition at line 781 of file prt.f90.

782  ! dummy
783  class(PrtModelType) :: this
784  character(len=*), intent(in) :: modelname
785 
786  ! allocate members from parent class
787  call this%NumericalModelType%allocate_scalars(modelname)
788 
789  ! allocate members that are part of model class
790  call mem_allocate(this%infmi, 'INFMI', this%memoryPath)
791  call mem_allocate(this%inmip, 'INMIP', this%memoryPath)
792  call mem_allocate(this%inmvt, 'INMVT', this%memoryPath)
793  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
794  call mem_allocate(this%inadv, 'INADV', this%memoryPath)
795  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
796  call mem_allocate(this%inssm, 'INSSM', this%memoryPath)
797  call mem_allocate(this%inoc, 'INOC ', this%memoryPath)
798 
799  this%infmi = 0
800  this%inmip = 0
801  this%inmvt = 0
802  this%inmst = 0
803  this%inadv = 0
804  this%indsp = 0
805  this%inssm = 0
806  this%inoc = 0

◆ create_bndpkgs()

subroutine prtmodule::create_bndpkgs ( class(prtmodeltype 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 
)

Definition at line 997 of file prt.f90.

999  ! modules
1002  ! dummy
1003  class(PrtModelType) :: this
1004  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1005  type(CharacterStringType), dimension(:), contiguous, &
1006  pointer, intent(inout) :: pkgtypes
1007  type(CharacterStringType), dimension(:), contiguous, &
1008  pointer, intent(inout) :: pkgnames
1009  type(CharacterStringType), dimension(:), contiguous, &
1010  pointer, intent(inout) :: mempaths
1011  integer(I4B), dimension(:), contiguous, &
1012  pointer, intent(inout) :: inunits
1013  ! local
1014  integer(I4B) :: ipakid, ipaknum
1015  character(len=LENFTYPE) :: pkgtype, bndptype
1016  character(len=LENPACKAGENAME) :: pkgname
1017  character(len=LENMEMPATH) :: mempath
1018  integer(I4B), pointer :: inunit
1019  integer(I4B) :: n
1020 
1021  if (allocated(bndpkgs)) then
1022  ! create stress packages
1023  ipakid = 1
1024  bndptype = ''
1025  do n = 1, size(bndpkgs)
1026  pkgtype = pkgtypes(bndpkgs(n))
1027  pkgname = pkgnames(bndpkgs(n))
1028  mempath = mempaths(bndpkgs(n))
1029  inunit => inunits(bndpkgs(n))
1030 
1031  if (bndptype /= pkgtype) then
1032  ipaknum = 1
1033  bndptype = pkgtype
1034  end if
1035 
1036  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1037  inunit, this%iout)
1038  ipakid = ipakid + 1
1039  ipaknum = ipaknum + 1
1040  end do
1041 
1042  ! cleanup
1043  deallocate (bndpkgs)
1044  end if
1045 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23

◆ create_packages()

subroutine prtmodule::create_packages ( class(prtmodeltype this)

Definition at line 1049 of file prt.f90.

1050  ! modules
1053  use arrayhandlersmodule, only: expandarray
1054  use memorymanagermodule, only: mem_setptr
1056  use simvariablesmodule, only: idm_context
1057  use budgetmodule, only: budget_cr
1061  use prtmipmodule, only: mip_cr
1062  use prtfmimodule, only: fmi_cr
1063  use prtocmodule, only: oc_cr
1064  ! dummy
1065  class(PrtModelType) :: this
1066  ! local
1067  type(CharacterStringType), dimension(:), contiguous, &
1068  pointer :: pkgtypes => null()
1069  type(CharacterStringType), dimension(:), contiguous, &
1070  pointer :: pkgnames => null()
1071  type(CharacterStringType), dimension(:), contiguous, &
1072  pointer :: mempaths => null()
1073  integer(I4B), dimension(:), contiguous, &
1074  pointer :: inunits => null()
1075  character(len=LENMEMPATH) :: model_mempath
1076  character(len=LENFTYPE) :: pkgtype
1077  character(len=LENPACKAGENAME) :: pkgname
1078  character(len=LENMEMPATH) :: mempath
1079  integer(I4B), pointer :: inunit
1080  integer(I4B), dimension(:), allocatable :: bndpkgs
1081  integer(I4B) :: n
1082  integer(I4B) :: indis = 0 ! DIS enabled flag
1083  character(len=LENMEMPATH) :: mempathmip = ''
1084  character(len=LENMEMPATH) :: mempathfmi = ''
1085 
1086  ! set input memory paths, input/model and input/model/namfile
1087  model_mempath = create_mem_path(component=this%name, context=idm_context)
1088 
1089  ! set pointers to model path package info
1090  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1091  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1092  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1093  call mem_setptr(inunits, 'INUNITS', model_mempath)
1094 
1095  do n = 1, size(pkgtypes)
1096  ! attributes for this input package
1097  pkgtype = pkgtypes(n)
1098  pkgname = pkgnames(n)
1099  mempath = mempaths(n)
1100  inunit => inunits(n)
1101 
1102  ! create dis package first as it is a prerequisite for other packages
1103  select case (pkgtype)
1104  case ('DIS6')
1105  indis = 1
1106  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1107  case ('DISV6')
1108  indis = 1
1109  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1110  case ('DISU6')
1111  indis = 1
1112  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1113  case ('MIP6')
1114  this%inmip = 1
1115  mempathmip = mempath
1116  case ('FMI6')
1117  this%infmi = 1
1118  mempathfmi = mempath
1119  case ('OC6')
1120  this%inoc = inunit
1121  case ('PRP6')
1122  call expandarray(bndpkgs)
1123  bndpkgs(size(bndpkgs)) = n
1124  case default
1125  call pstop(1, "Unrecognized package type: "//pkgtype)
1126  end select
1127  end do
1128 
1129  ! Create budget manager
1130  call budget_cr(this%budget, this%name)
1131 
1132  ! Create tracking method pools
1133  call create_method_pool()
1136 
1137  ! Create packages that are tied directly to model
1138  call mip_cr(this%mip, this%name, mempathmip, this%inmip, this%iout, this%dis)
1139  call fmi_cr(this%fmi, this%name, mempathfmi, this%infmi, this%iout)
1140  call oc_cr(this%oc, this%name, this%inoc, this%iout)
1141 
1142  ! Check to make sure that required ftype's have been specified
1143  call this%ftype_check(indis)
1144 
1145  ! Create boundary packages
1146  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
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
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Cell-level tracking methods.
subroutine, public create_method_cell_pool()
Create the cell method pool.
Model-level tracking methods.
Definition: MethodPool.f90:2
subroutine, public create_method_pool()
Create the method pool.
Definition: MethodPool.f90:18
Subcell-level tracking methods.
subroutine, public create_method_subcell_pool()
Create the subcell method pool.
subroutine, public fmi_cr(fmiobj, name_model, input_mempath, inunit, iout)
Create a new PrtFmi object.
Definition: prt-fmi.f90:38
subroutine, public mip_cr(mip, name_model, input_mempath, inunit, iout, dis)
Create a model input object.
Definition: prt-mip.f90:34
subroutine, public oc_cr(ocobj, name_model, inunit, iout)
@ brief Create an output control object
Definition: prt-oc.f90:52
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
Here is the call graph for this function:

◆ ftype_check()

subroutine prtmodule::ftype_check ( class(prtmodeltype this,
integer(i4b), intent(in)  indis 
)

Definition at line 890 of file prt.f90.

891  ! dummy
892  class(PrtModelType) :: this
893  integer(I4B), intent(in) :: indis
894  ! local
895  character(len=LINELENGTH) :: errmsg
896 
897  ! Check for DIS(u) and MIP. Stop if not present.
898  if (indis == 0) then
899  write (errmsg, '(1x,a)') &
900  'Discretization (DIS6, DISV6, or DISU6) package not specified.'
901  call store_error(errmsg)
902  end if
903  if (this%inmip == 0) then
904  write (errmsg, '(1x,a)') &
905  'Model input (MIP6) package not specified.'
906  call store_error(errmsg)
907  end if
908 
909  if (count_errors() > 0) then
910  write (errmsg, '(1x,a)') 'One or more required package(s) not specified.'
911  call store_error(errmsg)
912  call store_error_filename(this%filename)
913  end if
Here is the call graph for this function:

◆ log_namfile_options()

subroutine prtmodule::log_namfile_options ( class(prtmodeltype this,
type(prtnamparamfoundtype), intent(in)  found 
)

Definition at line 1150 of file prt.f90.

1152  class(PrtModelType) :: this
1153  type(PrtNamParamFoundType), intent(in) :: found
1154 
1155  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1156 
1157  if (found%print_input) then
1158  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1159  'FOR ALL MODEL STRESS PACKAGES'
1160  end if
1161 
1162  if (found%print_flows) then
1163  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1164  'FOR ALL MODEL PACKAGES'
1165  end if
1166 
1167  if (found%save_flows) then
1168  write (this%iout, '(4x,a)') &
1169  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1170  end if
1171 
1172  write (this%iout, '(1x,a)') 'END NAMEFILE OPTIONS:'

◆ package_create()

subroutine prtmodule::package_create ( class(prtmodeltype 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 
)

Definition at line 841 of file prt.f90.

843  ! modules
844  use constantsmodule, only: linelength
845  use prtprpmodule, only: prp_create
846  use apimodule, only: api_create
847  ! dummy
848  class(PrtModelType) :: this
849  character(len=*), intent(in) :: filtyp
850  character(len=LINELENGTH) :: errmsg
851  integer(I4B), intent(in) :: ipakid
852  integer(I4B), intent(in) :: ipaknum
853  character(len=*), intent(in) :: pakname
854  character(len=*), intent(in) :: mempath
855  integer(I4B), intent(in) :: inunit
856  integer(I4B), intent(in) :: iout
857  ! local
858  class(BndType), pointer :: packobj
859  class(BndType), pointer :: packobj2
860  integer(I4B) :: ip
861 
862  ! This part creates the package object
863  select case (filtyp)
864  case ('PRP6')
865  call prp_create(packobj, ipakid, ipaknum, inunit, iout, &
866  this%name, pakname, mempath, this%fmi)
867  case ('API6')
868  call api_create(packobj, ipakid, ipaknum, inunit, iout, &
869  this%name, pakname, mempath)
870  case default
871  write (errmsg, *) 'Invalid package type: ', filtyp
872  call store_error(errmsg, terminate=.true.)
873  end select
874 
875  ! Packages is the bndlist that is associated with the parent model
876  ! The following statement puts a pointer to this package in the ipakid
877  ! position of packages.
878  do ip = 1, this%bndlist%Count()
879  packobj2 => getbndfromlist(this%bndlist, ip)
880  if (packobj2%packName == pakname) then
881  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
882  'already exists: ', trim(pakname)
883  call store_error(errmsg, terminate=.true.)
884  end if
885  end do
886  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 prp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, input_mempath, fmi)
Create a new particle release point package.
Definition: prt-prp.f90:101
Here is the call graph for this function:

◆ prt_ad()

subroutine prtmodule::prt_ad ( class(prtmodeltype this)

Definition at line 346 of file prt.f90.

347  ! modules
349  ! dummy
350  class(PrtModelType) :: this
351  class(BndType), pointer :: packobj
352  ! local
353  integer(I4B) :: irestore
354  integer(I4B) :: ip, n, i
355 
356  ! Reset state variable
357  irestore = 0
358  if (ifailedstepretry > 0) irestore = 1
359 
360  ! Copy masssto into massstoold
361  do n = 1, this%dis%nodes
362  this%massstoold(n) = this%masssto(n)
363  end do
364 
365  ! Advance fmi
366  call this%fmi%fmi_ad()
367 
368  ! Advance
369  do ip = 1, this%bndlist%Count()
370  packobj => getbndfromlist(this%bndlist, ip)
371  call packobj%bnd_ad()
372  if (isimcheck > 0) then
373  call packobj%bnd_ck()
374  end if
375  end do
376  !
377  ! Initialize the flowja array. Flowja is calculated each time,
378  ! even if output is suppressed. (Flowja represents flow of particle
379  ! mass and is positive into a cell. Currently, each particle is assigned
380  ! unit mass.) Flowja is updated continually as particles are tracked
381  ! over the time step and at the end of the time step. The diagonal
382  ! position of the flowja array will contain the flow residual.
383  do i = 1, this%dis%nja
384  this%flowja(i) = dzero
385  end do
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:

◆ prt_ar()

subroutine prtmodule::prt_ar ( class(prtmodeltype this)

(1) allocates and reads packages part of this model, (2) allocates memory for arrays part of this model object

Definition at line 232 of file prt.f90.

233  ! modules
234  use constantsmodule, only: dhnoflo
235  use prtprpmodule, only: prtprptype
236  use prtmipmodule, only: prtmiptype
238  ! dummy
239  class(PrtModelType) :: this
240  ! locals
241  integer(I4B) :: ip, nprp
242  class(BndType), pointer :: packobj
243 
244  ! Set up basic packages
245  call this%fmi%fmi_ar(this%ibound)
246  if (this%inmip > 0) call this%mip%mip_ar()
247 
248  ! Set up output control and budget
249  call this%oc%oc_ar(this%dis, dhnoflo)
250  call this%budget%set_ibudcsv(this%oc%ibudcsv)
251 
252  ! Select tracking events
253  call this%tracks%select_events( &
254  this%oc%trackrelease, &
255  this%oc%trackfeatexit, &
256  this%oc%tracktimestep, &
257  this%oc%trackterminate, &
258  this%oc%trackweaksink, &
259  this%oc%trackusertime, &
260  this%oc%tracksubfexit)
261 
262  ! Set up boundary pkgs and pkg-scoped track files
263  nprp = 0
264  do ip = 1, this%bndlist%Count()
265  packobj => getbndfromlist(this%bndlist, ip)
266  select type (packobj)
267  type is (prtprptype)
268  nprp = nprp + 1
269  call packobj%prp_set_pointers(this%ibound, this%mip%izone)
270  call packobj%bnd_ar()
271  call packobj%bnd_ar()
272  if (packobj%itrkout > 0) then
273  call this%tracks%init_file( &
274  packobj%itrkout, &
275  iprp=nprp)
276  end if
277  if (packobj%itrkcsv > 0) then
278  call this%tracks%init_file( &
279  packobj%itrkcsv, &
280  csv=.true., &
281  iprp=nprp)
282  end if
283  class default
284  call packobj%bnd_ar()
285  end select
286  end do
287 
288  ! Set up model-scoped track files
289  if (this%oc%itrkout > 0) &
290  call this%tracks%init_file(this%oc%itrkout)
291  if (this%oc%itrkcsv > 0) &
292  call this%tracks%init_file(this%oc%itrkcsv, csv=.true.)
293 
294  ! Set up the tracking method
295  select type (dis => this%dis)
296  type is (distype)
297  call method_dis%init( &
298  fmi=this%fmi, &
299  events=this%events, &
300  izone=this%mip%izone, &
301  flowja=this%flowja, &
302  porosity=this%mip%porosity, &
303  retfactor=this%mip%retfactor, &
304  tracktimes=this%oc%tracktimes)
305  this%method => method_dis
306  type is (disvtype)
307  call method_disv%init( &
308  fmi=this%fmi, &
309  events=this%events, &
310  izone=this%mip%izone, &
311  flowja=this%flowja, &
312  porosity=this%mip%porosity, &
313  retfactor=this%mip%retfactor, &
314  tracktimes=this%oc%tracktimes)
315  this%method => method_disv
316  end select
317 
318  ! Subscribe track output manager to events
319  call this%events%subscribe(this%tracks)
320 
321  ! Set verbose tracing if requested
322  if (this%oc%dump_event_trace) this%tracks%iout = 0
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
type(methoddisvtype), pointer, public method_disv
Definition: MethodPool.f90:12
type(methoddistype), pointer, public method_dis
Definition: MethodPool.f90:11
Particle release point (PRP) package.
Definition: prt-prp.f90:38
Here is the call graph for this function:

◆ prt_bd()

subroutine prtmodule::prt_bd ( class(prtmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)

(1) Calculate intercell flows (flowja) (2) Calculate package contributions to model budget

Definition at line 488 of file prt.f90.

489  ! modules
490  use tdismodule, only: delt
491  use budgetmodule, only: rate_accumulator
492  ! dummy
493  class(PrtModelType) :: this
494  integer(I4B), intent(in) :: icnvg
495  integer(I4B), intent(in) :: isuppress_output
496  ! local
497  integer(I4B) :: ip
498  class(BndType), pointer :: packobj
499  real(DP) :: rin
500  real(DP) :: rout
501 
502  ! Budget routines (start by resetting). Sole purpose of this section
503  ! is to add in and outs to model budget. All ins and out for a model
504  ! should be added here to this%budget. In a subsequent exchange call,
505  ! exchange flows might also be added.
506  call this%budget%reset()
507  call rate_accumulator(this%ratesto, rin, rout)
508  call this%budget%addentry(rin, rout, delt, budtxt(1), &
509  isuppress_output, ' PRT')
510  do ip = 1, this%bndlist%Count()
511  packobj => getbndfromlist(this%bndlist, ip)
512  call packobj%bnd_bd(this%budget)
513  end do
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Here is the call graph for this function:

◆ prt_cq()

subroutine prtmodule::prt_cq ( class(prtmodeltype this,
integer(i4b), intent(in)  icnvg,
integer(i4b), intent(in)  isuppress_output 
)

Definition at line 389 of file prt.f90.

390  ! modules
391  use sparsemodule, only: csr_diagsum
392  use tdismodule, only: delt
393  use prtprpmodule, only: prtprptype
394  ! dummy
395  class(PrtModelType) :: this
396  integer(I4B), intent(in) :: icnvg
397  integer(I4B), intent(in) :: isuppress_output
398  ! local
399  integer(I4B) :: i
400  integer(I4B) :: ip
401  class(BndType), pointer :: packobj
402  real(DP) :: tled
403 
404  ! Flowja is calculated each time, even if output is suppressed.
405  ! Flowja represents flow of particle mass and is positive into a cell.
406  ! Currently, each particle is assigned unit mass.
407  !
408  ! Reciprocal of time step size.
409  tled = done / delt
410  !
411  ! Flowja was updated continually as particles were tracked over the
412  ! time step. At this point, flowja contains the net particle mass
413  ! exchanged between cells during the time step. To convert these to
414  ! flow rates (particle mass per time), divide by the time step size.
415  do i = 1, this%dis%nja
416  this%flowja(i) = this%flowja(i) * tled
417  end do
418 
419  ! Particle mass storage
420  call this%prt_cq_sto()
421 
422  ! Go through packages and call cq routines. Just a formality.
423  do ip = 1, this%bndlist%Count()
424  packobj => getbndfromlist(this%bndlist, ip)
425  call packobj%bnd_cq(this%masssto, this%flowja)
426  end do
427 
428  ! Finalize calculation of flowja by adding face flows to the diagonal.
429  ! This results in the flow residual being stored in the diagonal
430  ! position for each cell.
431  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:

◆ prt_cq_sto()

subroutine prtmodule::prt_cq_sto ( class(prtmodeltype this)

Definition at line 435 of file prt.f90.

436  ! modules
437  use tdismodule, only: delt
438  use prtprpmodule, only: prtprptype
439  ! dummy
440  class(PrtModelType) :: this
441  ! local
442  integer(I4B) :: ip
443  class(BndType), pointer :: packobj
444  integer(I4B) :: n
445  integer(I4B) :: np
446  integer(I4B) :: idiag
447  integer(I4B) :: istatus
448  real(DP) :: tled
449  real(DP) :: rate
450 
451  ! Reciprocal of time step size.
452  tled = done / delt
453 
454  ! Particle mass storage rate
455  do n = 1, this%dis%nodes
456  this%masssto(n) = dzero
457  this%ratesto(n) = dzero
458  end do
459  do ip = 1, this%bndlist%Count()
460  packobj => getbndfromlist(this%bndlist, ip)
461  select type (packobj)
462  type is (prtprptype)
463  do np = 1, packobj%nparticles
464  istatus = packobj%particles%istatus(np)
465  ! this may need to change if istatus flags change
466  if ((istatus > 0) .and. (istatus /= term_unreleased)) then
467  n = packobj%particles%itrdomain(np, level_feature)
468  ! Each particle currently assigned unit mass
469  this%masssto(n) = this%masssto(n) + done
470  end if
471  end do
472  end select
473  end do
474  do n = 1, this%dis%nodes
475  rate = -(this%masssto(n) - this%massstoold(n)) * tled
476  this%ratesto(n) = rate
477  idiag = this%dis%con%ia(n)
478  this%flowja(idiag) = this%flowja(idiag) + rate
479  end do
Here is the call graph for this function:

◆ prt_cr()

subroutine, public prtmodule::prt_cr ( character(len=*), intent(in)  filename,
integer(i4b), intent(in)  id,
character(len=*), intent(in)  modelname 
)

Definition at line 122 of file prt.f90.

123  ! modules
124  use listsmodule, only: basemodellist
127  use compilerversion
132  ! dummy
133  character(len=*), intent(in) :: filename
134  integer(I4B), intent(in) :: id
135  character(len=*), intent(in) :: modelname
136  ! local
137  type(PrtModelType), pointer :: this
138  class(BaseModelType), pointer :: model
139  character(len=LENMEMPATH) :: input_mempath
140  character(len=LINELENGTH) :: lst_fname
141  type(PrtNamParamFoundType) :: found
142 
143  ! Allocate a new PRT Model (this)
144  allocate (this)
145 
146  ! Set this before any allocs in the memory manager can be done
147  this%memoryPath = create_mem_path(modelname)
148 
149  ! Allocate event system and track output manager
150  allocate (this%events)
151  allocate (this%tracks)
152 
153  ! Allocate scalars and add model to basemodellist
154  call this%allocate_scalars(modelname)
155  model => this
156  call addbasemodeltolist(basemodellist, model)
157 
158  ! Assign variables
159  this%filename = filename
160  this%name = modelname
161  this%macronym = 'PRT'
162  this%id = id
163 
164  ! Set input model namfile memory path
165  input_mempath = create_mem_path(modelname, 'NAM', idm_context)
166 
167  ! Copy options from input context
168  call mem_set_value(this%iprpak, 'PRINT_INPUT', input_mempath, &
169  found%print_input)
170  call mem_set_value(this%iprflow, 'PRINT_FLOWS', input_mempath, &
171  found%print_flows)
172  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', input_mempath, &
173  found%save_flows)
174 
175  ! Create the list file
176  call this%create_lstfile(lst_fname, filename, found%list, &
177  'PARTICLE TRACKING MODEL (PRT)')
178 
179  ! Activate save_flows if found
180  if (found%save_flows) then
181  this%ipakcb = -1
182  end if
183 
184  ! Create model packages
185  call this%create_packages()
186 
187  ! Log options
188  if (this%iout > 0) then
189  call this%log_namfile_options(found)
190  end if
191 
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:161
type(listtype), public basemodellist
Definition: mf6lists.f90:16
Here is the call graph for this function:
Here is the caller graph for this function:

◆ prt_da()

subroutine prtmodule::prt_da ( class(prtmodeltype this)

Definition at line 716 of file prt.f90.

717  ! modules
724  ! dummy
725  class(PrtModelType) :: this
726  ! local
727  integer(I4B) :: ip
728  class(BndType), pointer :: packobj
729 
730  ! Deallocate idm memory
731  call memorystore_remove(this%name, 'NAM', idm_context)
732  call memorystore_remove(component=this%name, context=idm_context)
733 
734  ! Internal packages
735  call this%dis%dis_da()
736  call this%fmi%fmi_da()
737  call this%mip%mip_da()
738  call this%budget%budget_da()
739  call this%oc%oc_da()
740  deallocate (this%dis)
741  deallocate (this%fmi)
742  deallocate (this%mip)
743  deallocate (this%budget)
744  deallocate (this%oc)
745 
746  ! Method objects
749  call destroy_method_pool()
750 
751  ! Boundary packages
752  do ip = 1, this%bndlist%Count()
753  packobj => getbndfromlist(this%bndlist, ip)
754  call packobj%bnd_da()
755  deallocate (packobj)
756  end do
757 
758  ! Scalars
759  call mem_deallocate(this%infmi)
760  call mem_deallocate(this%inmip)
761  call mem_deallocate(this%inadv)
762  call mem_deallocate(this%indsp)
763  call mem_deallocate(this%inssm)
764  call mem_deallocate(this%inmst)
765  call mem_deallocate(this%inmvt)
766  call mem_deallocate(this%inoc)
767 
768  ! Arrays
769  call mem_deallocate(this%masssto)
770  call mem_deallocate(this%massstoold)
771  call mem_deallocate(this%ratesto)
772 
773  call this%tracks%destroy()
774  deallocate (this%events)
775  deallocate (this%tracks)
776 
777  call this%NumericalModelType%model_da()
subroutine, public memorystore_remove(component, subcomponent, context)
subroutine, public destroy_method_cell_pool()
Destroy the cell method pool.
subroutine, public destroy_method_pool()
Destroy the method pool.
Definition: MethodPool.f90:24
subroutine, public destroy_method_subcell_pool()
Destroy the subcell method pool.
Here is the call graph for this function:

◆ prt_df()

subroutine prtmodule::prt_df ( class(prtmodeltype this)

(1) call df routines for each package (2) set variables and pointers

Definition at line 199 of file prt.f90.

200  ! modules
201  use prtprpmodule, only: prtprptype
202  ! dummy
203  class(PrtModelType) :: this
204  ! local
205  integer(I4B) :: ip
206  class(BndType), pointer :: packobj
207 
208  ! Define packages and utility objects
209  call this%dis%dis_df()
210  call this%fmi%fmi_df(this%dis, 1)
211  call this%oc%oc_df()
212  call this%budget%budget_df(niunit_prt, 'MASS', 'M')
213 
214  ! Define packages and assign iout for time series managers
215  do ip = 1, this%bndlist%Count()
216  packobj => getbndfromlist(this%bndlist, ip)
217  call packobj%bnd_df(this%dis%nodes, this%dis)
218  packobj%TsManager%iout = this%iout
219  packobj%TasManager%iout = this%iout
220  end do
221 
222  ! Allocate model arrays
223  call this%allocate_arrays()
224 
Here is the call graph for this function:

◆ prt_ot()

subroutine prtmodule::prt_ot ( class(prtmodeltype this)

Definition at line 517 of file prt.f90.

518  use tdismodule, only: tdis_ot, endofperiod
519  ! dummy
520  class(PrtModelType) :: this
521  ! local
522  integer(I4B) :: idvsave
523  integer(I4B) :: idvprint
524  integer(I4B) :: icbcfl
525  integer(I4B) :: icbcun
526  integer(I4B) :: ibudfl
527  integer(I4B) :: ipflag
528 
529  ! Note: particle tracking output is handled elsewhere
530 
531  ! Set write and print flags
532  idvsave = 0
533  idvprint = 0
534  icbcfl = 0
535  ibudfl = 0
536  if (this%oc%oc_save('CONCENTRATION')) idvsave = 1
537  if (this%oc%oc_print('CONCENTRATION')) idvprint = 1
538  if (this%oc%oc_save('BUDGET')) icbcfl = 1
539  if (this%oc%oc_print('BUDGET')) ibudfl = 1
540  icbcun = this%oc%oc_save_unit('BUDGET')
541 
542  ! Override ibudfl and idvprint flags for nonconvergence
543  ! and end of period
544  ibudfl = this%oc%set_print_flag('BUDGET', 1, endofperiod)
545  idvprint = this%oc%set_print_flag('CONCENTRATION', 1, endofperiod)
546 
547  ! Save and print flows
548  call this%prt_ot_flow(icbcfl, ibudfl, icbcun)
549 
550  ! Save and print dependent variables
551  call this%prt_ot_dv(idvsave, idvprint, ipflag)
552 
553  ! Print budget summaries
554  call this%prt_ot_bdsummary(ibudfl, ipflag)
555 
556  ! Timing Output; if any dependent variables or budgets
557  ! are printed, then ipflag is set to 1.
558  if (ipflag == 1) call tdis_ot(this%iout)
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
subroutine, public tdis_ot(iout)
Print simulation time.
Definition: tdis.f90:274
Here is the call graph for this function:

◆ prt_ot_bdsummary()

subroutine prtmodule::prt_ot_bdsummary ( class(prtmodeltype this,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(inout)  ipflag 
)
private

Definition at line 686 of file prt.f90.

687  ! modules
688  use tdismodule, only: kstp, kper, totim, delt
689  ! dummy
690  class(PrtModelType) :: this
691  integer(I4B), intent(in) :: ibudfl
692  integer(I4B), intent(inout) :: ipflag
693  ! local
694  class(BndType), pointer :: packobj
695  integer(I4B) :: ip
696 
697  ! Package budget summary
698  do ip = 1, this%bndlist%Count()
699  packobj => getbndfromlist(this%bndlist, ip)
700  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
701  end do
702 
703  ! model budget summary
704  call this%budget%finalize_step(delt)
705  if (ibudfl /= 0) then
706  ipflag = 1
707  ! model budget summary
708  call this%budget%budget_ot(kstp, kper, this%iout)
709  end if
710 
711  ! Write to budget csv
712  call this%budget%writecsv(totim)
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
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:

◆ prt_ot_dv()

subroutine prtmodule::prt_ot_dv ( class(prtmodeltype this,
integer(i4b), intent(in)  idvsave,
integer(i4b), intent(in)  idvprint,
integer(i4b), intent(inout)  ipflag 
)

Definition at line 665 of file prt.f90.

666  ! dummy
667  class(PrtModelType) :: this
668  integer(I4B), intent(in) :: idvsave
669  integer(I4B), intent(in) :: idvprint
670  integer(I4B), intent(inout) :: ipflag
671  ! local
672  class(BndType), pointer :: packobj
673  integer(I4B) :: ip
674 
675  ! Print advanced package dependent variables
676  do ip = 1, this%bndlist%Count()
677  packobj => getbndfromlist(this%bndlist, ip)
678  call packobj%bnd_ot_dv(idvsave, idvprint)
679  end do
680 
681  ! save head and print head
682  call this%oc%oc_ot(ipflag)
Here is the call graph for this function:

◆ prt_ot_flow()

subroutine prtmodule::prt_ot_flow ( class(prtmodeltype this,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  ibudfl,
integer(i4b), intent(in)  icbcun 
)

Definition at line 562 of file prt.f90.

563  use prtprpmodule, only: prtprptype
564  class(PrtModelType) :: this
565  integer(I4B), intent(in) :: icbcfl
566  integer(I4B), intent(in) :: ibudfl
567  integer(I4B), intent(in) :: icbcun
568  class(BndType), pointer :: packobj
569  integer(I4B) :: ip
570 
571  ! Save PRT flows
572  call this%prt_ot_saveflow(this%dis%nja, this%flowja, icbcfl, icbcun)
573  do ip = 1, this%bndlist%Count()
574  packobj => getbndfromlist(this%bndlist, ip)
575  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
576  end do
577 
578  ! Save advanced package flows
579  do ip = 1, this%bndlist%Count()
580  packobj => getbndfromlist(this%bndlist, ip)
581  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
582  end do
583 
584  ! Print PRT flows
585  call this%prt_ot_printflow(ibudfl, this%flowja)
586  do ip = 1, this%bndlist%Count()
587  packobj => getbndfromlist(this%bndlist, ip)
588  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
589  end do
590 
591  ! Print advanced package flows
592  do ip = 1, this%bndlist%Count()
593  packobj => getbndfromlist(this%bndlist, ip)
594  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
595  end do
Here is the call graph for this function:

◆ prt_ot_printflow()

subroutine prtmodule::prt_ot_printflow ( class(prtmodeltype this,
integer(i4b), intent(in)  ibudfl,
real(dp), dimension(:), intent(inout)  flowja 
)
private

Definition at line 626 of file prt.f90.

627  ! modules
628  use tdismodule, only: kper, kstp
629  use constantsmodule, only: lenbigline
630  ! dummy
631  class(PrtModelType) :: this
632  integer(I4B), intent(in) :: ibudfl
633  real(DP), intent(inout), dimension(:) :: flowja
634  ! local
635  character(len=LENBIGLINE) :: line
636  character(len=30) :: tempstr
637  integer(I4B) :: n, ipos, m
638  real(DP) :: qnm
639  ! formats
640  character(len=*), parameter :: fmtiprflow = &
641  "(/,4x,'CALCULATED INTERCELL FLOW &
642  &FOR PERIOD ', i0, ' STEP ', i0)"
643 
644  ! Write flowja to list file if requested
645  if (ibudfl /= 0 .and. this%iprflow > 0) then
646  write (this%iout, fmtiprflow) kper, kstp
647  do n = 1, this%dis%nodes
648  line = ''
649  call this%dis%noder_to_string(n, tempstr)
650  line = trim(tempstr)//':'
651  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
652  m = this%dis%con%ja(ipos)
653  call this%dis%noder_to_string(m, tempstr)
654  line = trim(line)//' '//trim(tempstr)
655  qnm = flowja(ipos)
656  write (tempstr, '(1pg15.6)') qnm
657  line = trim(line)//' '//trim(adjustl(tempstr))
658  end do
659  write (this%iout, '(a)') trim(line)
660  end do
661  end if
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15

◆ prt_ot_saveflow()

subroutine prtmodule::prt_ot_saveflow ( class(prtmodeltype this,
integer(i4b), intent(in)  nja,
real(dp), dimension(nja), intent(in)  flowja,
integer(i4b), intent(in)  icbcfl,
integer(i4b), intent(in)  icbcun 
)

Definition at line 599 of file prt.f90.

600  ! dummy
601  class(PrtModelType) :: this
602  integer(I4B), intent(in) :: nja
603  real(DP), dimension(nja), intent(in) :: flowja
604  integer(I4B), intent(in) :: icbcfl
605  integer(I4B), intent(in) :: icbcun
606  ! local
607  integer(I4B) :: ibinun
608 
609  ! Set unit number for binary output
610  if (this%ipakcb < 0) then
611  ibinun = icbcun
612  elseif (this%ipakcb == 0) then
613  ibinun = 0
614  else
615  ibinun = this%ipakcb
616  end if
617  if (icbcfl == 0) ibinun = 0
618 
619  ! Write the face flows if requested
620  if (ibinun /= 0) then
621  call this%dis%record_connection_array(flowja, ibinun, this%iout)
622  end if

◆ prt_rp()

subroutine prtmodule::prt_rp ( class(prtmodeltype this)

Definition at line 326 of file prt.f90.

327  use tdismodule, only: readnewdata
328  ! dummy
329  class(PrtModelType) :: this
330  ! local
331  class(BndType), pointer :: packobj
332  integer(I4B) :: ip
333 
334  ! Check with TDIS on whether or not it is time to RP
335  if (.not. readnewdata) return
336 
337  ! Read and prepare
338  if (this%inoc > 0) call this%oc%oc_rp()
339  do ip = 1, this%bndlist%Count()
340  packobj => getbndfromlist(this%bndlist, ip)
341  call packobj%bnd_rp()
342  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:

◆ prt_solve()

subroutine prtmodule::prt_solve ( class(prtmodeltype this)
private

Definition at line 917 of file prt.f90.

919  use prtprpmodule, only: prtprptype
922  ! dummy
923  class(PrtModelType) :: this
924  ! local
925  integer(I4B) :: np, ip
926  class(BndType), pointer :: packobj
927  type(ParticleType), pointer :: particle
928  real(DP) :: tmax
929  integer(I4B) :: iprp
930 
931  ! A single particle is reused in the tracking loops
932  ! to avoid allocating and deallocating it each time.
933  ! get() and put() retrieve and store particle state.
934  call create_particle(particle)
935  ! Loop over PRP packages and particles within them.
936  iprp = 0
937  do ip = 1, this%bndlist%Count()
938  packobj => getbndfromlist(this%bndlist, ip)
939  select type (packobj)
940  type is (prtprptype)
941  iprp = iprp + 1
942  do np = 1, packobj%nparticles
943  ! Get the particle from the store
944  call packobj%particles%get(particle, this%id, iprp, np)
945  ! If particle is permanently unreleased, cycle.
946  ! Raise a termination event if we haven't yet.
947  ! TODO: when we have generic dynamic vectors,
948  ! consider terminating permanently unreleased
949  ! in PRP instead of here. For now, status -8
950  ! indicates the permanently unreleased event
951  ! is not yet recorded, status 8 it has been.
952  if (particle%istatus == (-1 * term_unreleased)) then
953  call this%method%terminate(particle, status=term_unreleased)
954  call packobj%particles%put(particle, np)
955  end if
956  if (particle%istatus > active) cycle ! Skip terminated particles
957  particle%istatus = active ! Set active status in case of release
958  ! If the particle was released this time step, emit a release event
959  if (particle%trelease >= totimc) call this%method%release(particle)
960  ! Maximum time is the end of the time step or the particle
961  ! stop time, whichever comes first, unless it's the final
962  ! time step and the extend option is on, in which case
963  ! it's just the particle stop time.
964  if (endofsimulation .and. particle%iextend > 0) then
965  tmax = particle%tstop
966  else
967  tmax = min(totimc + delt, particle%tstop)
968  end if
969  ! Apply the tracking method until the maximum time.
970  call this%method%apply(particle, tmax)
971  ! Reset cell/zone one-backs, used for cycle detection.
972  ! TODO can remove when we have better cycle detection
973  particle%icp = 0
974  particle%izp = 0
975  ! If the particle timed out, terminate it.
976  ! "Timed out" means it's still active but
977  ! - it reached its stop time, or
978  ! - the simulation is over and not extended.
979  ! We can't detect timeout within the tracking
980  ! method because the method just receives the
981  ! maximum time with no context on what it is.
982  ! TODO maybe think about changing that?
983  if (particle%istatus <= active .and. &
984  (particle%ttrack == particle%tstop .or. &
985  (endofsimulation .and. particle%iextend == 0))) &
986  call this%method%terminate(particle, status=term_timeout)
987  ! Return the particle to the store
988  call packobj%particles%put(particle, np)
989  end do
990  end select
991  end do
992  call particle%destroy()
993  deallocate (particle)
@, public release
particle was released
@, public terminate
particle terminated
@ term_timeout
terminated at stop time or end of simulation
Definition: Particle.f90:37
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:35
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
Here is the call graph for this function:

Variable Documentation

◆ budtxt

character(len=lenbudtxt), dimension(nbditems) prtmodule::budtxt
private

Definition at line 39 of file prt.f90.

39  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt

◆ nbditems

integer(i4b), parameter prtmodule::nbditems = 1
private

Definition at line 38 of file prt.f90.

38  integer(I4B), parameter :: NBDITEMS = 1

◆ niunit_prt

integer(i4b), parameter prtmodule::niunit_prt = PRT_NBASEPKG + PRT_NMULTIPKG
private

Definition at line 117 of file prt.f90.

117  integer(I4B), parameter :: NIUNIT_PRT = prt_nbasepkg + prt_nmultipkg

◆ prt_basepkg

character(len=lenpackagetype), dimension(prt_nbasepkg), public prtmodule::prt_basepkg

Definition at line 98 of file prt.f90.

98  character(len=LENPACKAGETYPE), dimension(PRT_NBASEPKG) :: PRT_BASEPKG

◆ prt_multipkg

character(len=lenpackagetype), dimension(prt_nmultipkg), public prtmodule::prt_multipkg

Definition at line 112 of file prt.f90.

112  character(len=LENPACKAGETYPE), dimension(PRT_NMULTIPKG) :: PRT_MULTIPKG

◆ prt_nbasepkg

integer(i4b), parameter, public prtmodule::prt_nbasepkg = 50

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

Definition at line 97 of file prt.f90.

97  integer(I4B), parameter :: PRT_NBASEPKG = 50

◆ prt_nmultipkg

integer(i4b), parameter, public prtmodule::prt_nmultipkg = 50

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

Definition at line 111 of file prt.f90.

111  integer(I4B), parameter :: PRT_NMULTIPKG = 50