MODFLOW 6  version 6.6.0.dev0
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 non-allocatable members. 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 781 of file prt.f90.

783  class(PrtModelType) :: this
784  integer(I4B) :: n
785 
786  ! -- Allocate arrays in parent type
787  this%nja = this%dis%nja
788  call this%NumericalModelType%allocate_arrays()
789 
790  ! -- Allocate and initialize arrays
791  call mem_allocate(this%masssto, this%dis%nodes, &
792  'MASSSTO', this%memoryPath)
793  call mem_allocate(this%massstoold, this%dis%nodes, &
794  'MASSSTOOLD', this%memoryPath)
795  call mem_allocate(this%ratesto, this%dis%nodes, &
796  'RATESTO', this%memoryPath)
797  ! -- explicit model, so these must be manually allocated
798  call mem_allocate(this%x, this%dis%nodes, 'X', this%memoryPath)
799  call mem_allocate(this%rhs, this%dis%nodes, 'RHS', this%memoryPath)
800  call mem_allocate(this%ibound, this%dis%nodes, 'IBOUND', this%memoryPath)
801  do n = 1, this%dis%nodes
802  this%masssto(n) = dzero
803  this%massstoold(n) = dzero
804  this%ratesto(n) = dzero
805  this%x(n) = dzero
806  this%rhs(n) = dzero
807  this%ibound(n) = 1
808  end do

◆ allocate_scalars()

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

Definition at line 752 of file prt.f90.

753  ! -- dummy
754  class(PrtModelType) :: this
755  character(len=*), intent(in) :: modelname
756 
757  ! -- allocate members from parent class
758  call this%NumericalModelType%allocate_scalars(modelname)
759 
760  ! -- allocate members that are part of model class
761  call mem_allocate(this%infmi, 'INFMI', this%memoryPath)
762  call mem_allocate(this%inmip, 'INMIP', this%memoryPath)
763  call mem_allocate(this%inmvt, 'INMVT', this%memoryPath)
764  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
765  call mem_allocate(this%inadv, 'INADV', this%memoryPath)
766  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
767  call mem_allocate(this%inssm, 'INSSM', this%memoryPath)
768  call mem_allocate(this%inoc, 'INOC ', this%memoryPath)
769 
770  this%infmi = 0
771  this%inmip = 0
772  this%inmvt = 0
773  this%inmst = 0
774  this%inadv = 0
775  this%indsp = 0
776  this%inssm = 0
777  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 970 of file prt.f90.

972  ! -- modules
975  ! -- dummy
976  class(PrtModelType) :: this
977  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
978  type(CharacterStringType), dimension(:), contiguous, &
979  pointer, intent(inout) :: pkgtypes
980  type(CharacterStringType), dimension(:), contiguous, &
981  pointer, intent(inout) :: pkgnames
982  type(CharacterStringType), dimension(:), contiguous, &
983  pointer, intent(inout) :: mempaths
984  integer(I4B), dimension(:), contiguous, &
985  pointer, intent(inout) :: inunits
986  ! -- local
987  integer(I4B) :: ipakid, ipaknum
988  character(len=LENFTYPE) :: pkgtype, bndptype
989  character(len=LENPACKAGENAME) :: pkgname
990  character(len=LENMEMPATH) :: mempath
991  integer(I4B), pointer :: inunit
992  integer(I4B) :: n
993 
994  if (allocated(bndpkgs)) then
995  !
996  ! -- create stress packages
997  ipakid = 1
998  bndptype = ''
999  do n = 1, size(bndpkgs)
1000  !
1001  pkgtype = pkgtypes(bndpkgs(n))
1002  pkgname = pkgnames(bndpkgs(n))
1003  mempath = mempaths(bndpkgs(n))
1004  inunit => inunits(bndpkgs(n))
1005  !
1006  if (bndptype /= pkgtype) then
1007  ipaknum = 1
1008  bndptype = pkgtype
1009  end if
1010  !
1011  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1012  inunit, this%iout)
1013  ipakid = ipakid + 1
1014  ipaknum = ipaknum + 1
1015  end do
1016  !
1017  ! -- cleanup
1018  deallocate (bndpkgs)
1019  end if
1020 
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 1024 of file prt.f90.

1025  ! -- modules
1028  use arrayhandlersmodule, only: expandarray
1029  use memorymanagermodule, only: mem_setptr
1031  use simvariablesmodule, only: idm_context
1032  use budgetmodule, only: budget_cr
1036  use prtmipmodule, only: mip_cr
1037  use prtfmimodule, only: fmi_cr
1038  use prtocmodule, only: oc_cr
1039  ! -- dummy
1040  class(PrtModelType) :: this
1041  ! -- local
1042  type(CharacterStringType), dimension(:), contiguous, &
1043  pointer :: pkgtypes => null()
1044  type(CharacterStringType), dimension(:), contiguous, &
1045  pointer :: pkgnames => null()
1046  type(CharacterStringType), dimension(:), contiguous, &
1047  pointer :: mempaths => null()
1048  integer(I4B), dimension(:), contiguous, &
1049  pointer :: inunits => null()
1050  character(len=LENMEMPATH) :: model_mempath
1051  character(len=LENFTYPE) :: pkgtype
1052  character(len=LENPACKAGENAME) :: pkgname
1053  character(len=LENMEMPATH) :: mempath
1054  integer(I4B), pointer :: inunit
1055  integer(I4B), dimension(:), allocatable :: bndpkgs
1056  integer(I4B) :: n
1057  integer(I4B) :: indis = 0 ! DIS enabled flag
1058  character(len=LENMEMPATH) :: mempathmip = ''
1059 
1060  ! -- set input memory paths, input/model and input/model/namfile
1061  model_mempath = create_mem_path(component=this%name, context=idm_context)
1062 
1063  ! -- set pointers to model path package info
1064  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1065  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1066  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1067  call mem_setptr(inunits, 'INUNITS', model_mempath)
1068 
1069  do n = 1, size(pkgtypes)
1070  ! attributes for this input package
1071  pkgtype = pkgtypes(n)
1072  pkgname = pkgnames(n)
1073  mempath = mempaths(n)
1074  inunit => inunits(n)
1075 
1076  ! -- create dis package first as it is a prerequisite for other packages
1077  select case (pkgtype)
1078  case ('DIS6')
1079  indis = 1
1080  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1081  case ('DISV6')
1082  indis = 1
1083  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1084  case ('DISU6')
1085  indis = 1
1086  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1087  case ('MIP6')
1088  this%inmip = 1
1089  mempathmip = mempath
1090  case ('FMI6')
1091  this%infmi = inunit
1092  case ('OC6')
1093  this%inoc = inunit
1094  case ('PRP6')
1095  call expandarray(bndpkgs)
1096  bndpkgs(size(bndpkgs)) = n
1097  case default
1098  call pstop(1, "Unrecognized package type: "//pkgtype)
1099  end select
1100  end do
1101 
1102  ! -- Create budget manager
1103  call budget_cr(this%budget, this%name)
1104 
1105  ! -- Create tracking method pools
1106  call create_method_pool()
1109 
1110  ! -- Create packages that are tied directly to model
1111  call mip_cr(this%mip, this%name, mempathmip, this%inmip, this%iout, this%dis)
1112  call fmi_cr(this%fmi, this%name, this%infmi, this%iout)
1113  call oc_cr(this%oc, this%name, this%inoc, this%iout)
1114 
1115  ! -- Check to make sure that required ftype's have been specified
1116  call this%ftype_check(indis)
1117 
1118  ! -- Create boundary packages
1119  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, 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:51
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 861 of file prt.f90.

862  ! -- dummy
863  class(PrtModelType) :: this
864  integer(I4B), intent(in) :: indis
865  ! -- local
866  character(len=LINELENGTH) :: errmsg
867 
868  ! -- Check for DIS(u) and MIP. Stop if not present.
869  if (indis == 0) then
870  write (errmsg, '(1x,a)') &
871  'Discretization (DIS6, DISV6, or DISU6) package not specified.'
872  call store_error(errmsg)
873  end if
874  if (this%inmip == 0) then
875  write (errmsg, '(1x,a)') &
876  'Model input (MIP6) package not specified.'
877  call store_error(errmsg)
878  end if
879 
880  if (count_errors() > 0) then
881  write (errmsg, '(1x,a)') 'One or more required package(s) not specified.'
882  call store_error(errmsg)
883  call store_error_filename(this%filename)
884  end if
Here is the call graph for this function:

◆ log_namfile_options()

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

Definition at line 1123 of file prt.f90.

1125  class(PrtModelType) :: this
1126  type(GwfNamParamFoundType), intent(in) :: found
1127 
1128  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1129 
1130  if (found%newton) then
1131  write (this%iout, '(4x,a)') &
1132  'NEWTON-RAPHSON method enabled for the model.'
1133  if (found%under_relaxation) then
1134  write (this%iout, '(4x,a,a)') &
1135  'NEWTON-RAPHSON UNDER-RELAXATION based on the bottom ', &
1136  'elevation of the model will be applied to the model.'
1137  end if
1138  end if
1139 
1140  if (found%print_input) then
1141  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1142  'FOR ALL MODEL STRESS PACKAGES'
1143  end if
1144 
1145  if (found%print_flows) then
1146  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1147  'FOR ALL MODEL PACKAGES'
1148  end if
1149 
1150  if (found%save_flows) then
1151  write (this%iout, '(4x,a)') &
1152  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1153  end if
1154 
1155  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 812 of file prt.f90.

814  ! -- modules
815  use constantsmodule, only: linelength
816  use prtprpmodule, only: prp_create
817  use apimodule, only: api_create
818  ! -- dummy
819  class(PrtModelType) :: this
820  character(len=*), intent(in) :: filtyp
821  character(len=LINELENGTH) :: errmsg
822  integer(I4B), intent(in) :: ipakid
823  integer(I4B), intent(in) :: ipaknum
824  character(len=*), intent(in) :: pakname
825  character(len=*), intent(in) :: mempath
826  integer(I4B), intent(in) :: inunit
827  integer(I4B), intent(in) :: iout
828  ! -- local
829  class(BndType), pointer :: packobj
830  class(BndType), pointer :: packobj2
831  integer(I4B) :: ip
832 
833  ! -- This part creates the package object
834  select case (filtyp)
835  case ('PRP6')
836  call prp_create(packobj, ipakid, ipaknum, inunit, iout, &
837  this%name, pakname, this%fmi)
838  case ('API6')
839  call api_create(packobj, ipakid, ipaknum, inunit, iout, &
840  this%name, pakname)
841  case default
842  write (errmsg, *) 'Invalid package type: ', filtyp
843  call store_error(errmsg, terminate=.true.)
844  end select
845 
846  ! -- Packages is the bndlist that is associated with the parent model
847  ! -- The following statement puts a pointer to this package in the ipakid
848  ! -- position of packages.
849  do ip = 1, this%bndlist%Count()
850  packobj2 => getbndfromlist(this%bndlist, ip)
851  if (packobj2%packName == pakname) then
852  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
853  'already exists: ', trim(pakname)
854  call store_error(errmsg, terminate=.true.)
855  end if
856  end do
857  call addbndtolist(this%bndlist, packobj)
This module contains the API package methods.
Definition: gwf-api.f90:12
subroutine, public api_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
@ brief Create a new package object
Definition: gwf-api.f90:49
subroutine, public prp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, 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 317 of file prt.f90.

318  ! -- modules
320  ! -- dummy
321  class(PrtModelType) :: this
322  class(BndType), pointer :: packobj
323  ! -- local
324  integer(I4B) :: irestore
325  integer(I4B) :: ip, n, i
326 
327  ! -- Reset state variable
328  irestore = 0
329  if (ifailedstepretry > 0) irestore = 1
330 
331  ! -- Copy masssto into massstoold
332  do n = 1, this%dis%nodes
333  this%massstoold(n) = this%masssto(n)
334  end do
335 
336  ! -- Advance fmi
337  call this%fmi%fmi_ad()
338 
339  ! -- Advance
340  do ip = 1, this%bndlist%Count()
341  packobj => getbndfromlist(this%bndlist, ip)
342  call packobj%bnd_ad()
343  if (isimcheck > 0) then
344  call packobj%bnd_ck()
345  end if
346  end do
347  !
348  ! -- Initialize the flowja array. Flowja is calculated each time,
349  ! even if output is suppressed. (Flowja represents flow of particle
350  ! mass and is positive into a cell. Currently, each particle is assigned
351  ! unit mass.) Flowja is updated continually as particles are tracked
352  ! over the time step and at the end of the time step. The diagonal
353  ! position of the flowja array will contain the flow residual.
354  do i = 1, this%dis%nja
355  this%flowja(i) = dzero
356  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 226 of file prt.f90.

227  ! -- modules
228  use constantsmodule, only: dhnoflo
229  use prtprpmodule, only: prtprptype
230  use prtmipmodule, only: prtmiptype
232  ! -- dummy
233  class(PrtModelType) :: this
234  ! -- locals
235  integer(I4B) :: ip
236  class(BndType), pointer :: packobj
237 
238  ! -- Allocate and read modules attached to model
239  call this%fmi%fmi_ar(this%ibound)
240  if (this%inmip > 0) call this%mip%mip_ar()
241 
242  ! -- set up output control
243  call this%oc%oc_ar(this%dis, dhnoflo)
244  call this%budget%set_ibudcsv(this%oc%ibudcsv)
245 
246  ! -- Package input files now open, so allocate and read
247  do ip = 1, this%bndlist%Count()
248  packobj => getbndfromlist(this%bndlist, ip)
249  select type (packobj)
250  type is (prtprptype)
251  call packobj%prp_set_pointers(this%ibound, this%mip%izone, &
252  this%trackfilectl)
253  end select
254  ! -- Read and allocate package
255  call packobj%bnd_ar()
256  end do
257 
258  ! -- Initialize tracking method
259  select type (dis => this%dis)
260  type is (distype)
261  call method_dis%init( &
262  fmi=this%fmi, &
263  trackfilectl=this%trackfilectl, &
264  izone=this%mip%izone, &
265  flowja=this%flowja, &
266  porosity=this%mip%porosity, &
267  retfactor=this%mip%retfactor, &
268  tracktimes=this%oc%tracktimes)
269  this%method => method_dis
270  type is (disvtype)
271  call method_disv%init( &
272  fmi=this%fmi, &
273  trackfilectl=this%trackfilectl, &
274  izone=this%mip%izone, &
275  flowja=this%flowja, &
276  porosity=this%mip%porosity, &
277  retfactor=this%mip%retfactor, &
278  tracktimes=this%oc%tracktimes)
279  this%method => method_disv
280  end select
281 
282  ! -- Initialize track output files and reporting options
283  if (this%oc%itrkout > 0) &
284  call this%trackfilectl%init_track_file(this%oc%itrkout)
285  if (this%oc%itrkcsv > 0) &
286  call this%trackfilectl%init_track_file(this%oc%itrkcsv, csv=.true.)
287  call this%trackfilectl%set_track_events( &
288  this%oc%trackrelease, &
289  this%oc%trackexit, &
290  this%oc%tracktimestep, &
291  this%oc%trackterminate, &
292  this%oc%trackweaksink, &
293  this%oc%trackusertime)
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:41
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 459 of file prt.f90.

460  ! -- modules
461  use tdismodule, only: delt
462  use budgetmodule, only: rate_accumulator
463  ! -- dummy
464  class(PrtModelType) :: this
465  integer(I4B), intent(in) :: icnvg
466  integer(I4B), intent(in) :: isuppress_output
467  ! -- local
468  integer(I4B) :: ip
469  class(BndType), pointer :: packobj
470  real(DP) :: rin
471  real(DP) :: rout
472 
473  ! -- Budget routines (start by resetting). Sole purpose of this section
474  ! is to add in and outs to model budget. All ins and out for a model
475  ! should be added here to this%budget. In a subsequent exchange call,
476  ! exchange flows might also be added.
477  call this%budget%reset()
478  call rate_accumulator(this%ratesto, rin, rout)
479  call this%budget%addentry(rin, rout, delt, budtxt(1), &
480  isuppress_output, ' PRT')
481  do ip = 1, this%bndlist%Count()
482  packobj => getbndfromlist(this%bndlist, ip)
483  call packobj%bnd_bd(this%budget)
484  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 360 of file prt.f90.

361  ! -- modules
362  use sparsemodule, only: csr_diagsum
363  use tdismodule, only: delt
364  use prtprpmodule, only: prtprptype
365  ! -- dummy
366  class(PrtModelType) :: this
367  integer(I4B), intent(in) :: icnvg
368  integer(I4B), intent(in) :: isuppress_output
369  ! -- local
370  integer(I4B) :: i
371  integer(I4B) :: ip
372  class(BndType), pointer :: packobj
373  real(DP) :: tled
374 
375  ! -- Flowja is calculated each time, even if output is suppressed.
376  ! Flowja represents flow of particle mass and is positive into a cell.
377  ! Currently, each particle is assigned unit mass.
378  !
379  ! -- Reciprocal of time step size.
380  tled = done / delt
381  !
382  ! -- Flowja was updated continually as particles were tracked over the
383  ! time step. At this point, flowja contains the net particle mass
384  ! exchanged between cells during the time step. To convert these to
385  ! flow rates (particle mass per time), divide by the time step size.
386  do i = 1, this%dis%nja
387  this%flowja(i) = this%flowja(i) * tled
388  end do
389 
390  ! -- Particle mass storage
391  call this%prt_cq_sto()
392 
393  ! -- Go through packages and call cq routines. Just a formality.
394  do ip = 1, this%bndlist%Count()
395  packobj => getbndfromlist(this%bndlist, ip)
396  call packobj%bnd_cq(this%masssto, this%flowja)
397  end do
398 
399  ! -- Finalize calculation of flowja by adding face flows to the diagonal.
400  ! This results in the flow residual being stored in the diagonal
401  ! position for each cell.
402  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 406 of file prt.f90.

407  ! -- modules
408  use tdismodule, only: delt
409  use prtprpmodule, only: prtprptype
410  ! -- dummy
411  class(PrtModelType) :: this
412  ! -- local
413  integer(I4B) :: ip
414  class(BndType), pointer :: packobj
415  integer(I4B) :: n
416  integer(I4B) :: np
417  integer(I4B) :: idiag
418  integer(I4B) :: istatus
419  real(DP) :: tled
420  real(DP) :: rate
421 
422  ! -- Reciprocal of time step size.
423  tled = done / delt
424 
425  ! -- Particle mass storage rate
426  do n = 1, this%dis%nodes
427  this%masssto(n) = dzero
428  this%ratesto(n) = dzero
429  end do
430  do ip = 1, this%bndlist%Count()
431  packobj => getbndfromlist(this%bndlist, ip)
432  select type (packobj)
433  type is (prtprptype)
434  do np = 1, packobj%nparticles
435  istatus = packobj%particles%istatus(np)
436  ! this may need to change if istatus flags change
437  if ((istatus > 0) .and. (istatus /= 8)) then
438  n = packobj%particles%idomain(np, 2)
439  ! -- Each particle currently assigned unit mass
440  this%masssto(n) = this%masssto(n) + done
441  end if
442  end do
443  end select
444  end do
445  do n = 1, this%dis%nodes
446  rate = -(this%masssto(n) - this%massstoold(n)) * tled
447  this%ratesto(n) = rate
448  idiag = this%dis%con%ia(n)
449  this%flowja(idiag) = this%flowja(idiag) + rate
450  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 118 of file prt.f90.

119  ! -- modules
120  use listsmodule, only: basemodellist
123  use compilerversion
128  ! -- dummy
129  character(len=*), intent(in) :: filename
130  integer(I4B), intent(in) :: id
131  character(len=*), intent(in) :: modelname
132  ! -- local
133  type(PrtModelType), pointer :: this
134  class(BaseModelType), pointer :: model
135  character(len=LENMEMPATH) :: input_mempath
136  character(len=LINELENGTH) :: lst_fname
137  type(GwfNamParamFoundType) :: found
138 
139  ! -- Allocate a new PRT Model (this)
140  allocate (this)
141 
142  ! -- Set this before any allocs in the memory manager can be done
143  this%memoryPath = create_mem_path(modelname)
144 
145  ! -- Allocate track control object
146  allocate (this%trackfilectl)
147 
148  ! -- Allocate scalars and add model to basemodellist
149  call this%allocate_scalars(modelname)
150  model => this
151  call addbasemodeltolist(basemodellist, model)
152 
153  ! -- Assign variables
154  this%filename = filename
155  this%name = modelname
156  this%macronym = 'PRT'
157  this%id = id
158 
159  ! -- Set input model namfile memory path
160  input_mempath = create_mem_path(modelname, 'NAM', idm_context)
161 
162  ! -- Copy options from input context
163  call mem_set_value(this%iprpak, 'PRINT_INPUT', input_mempath, &
164  found%print_input)
165  call mem_set_value(this%iprflow, 'PRINT_FLOWS', input_mempath, &
166  found%print_flows)
167  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', input_mempath, &
168  found%save_flows)
169 
170  ! -- Create the list file
171  call this%create_lstfile(lst_fname, filename, found%list, &
172  'PARTICLE TRACKING MODEL (PRT)')
173 
174  ! -- Activate save_flows if found
175  if (found%save_flows) then
176  this%ipakcb = -1
177  end if
178 
179  ! -- Log options
180  if (this%iout > 0) then
181  call this%log_namfile_options(found)
182  end if
183 
184  ! -- Create model packages
185  call this%create_packages()
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 687 of file prt.f90.

688  ! -- modules
695  ! -- dummy
696  class(PrtModelType) :: this
697  ! -- local
698  integer(I4B) :: ip
699  class(BndType), pointer :: packobj
700 
701  ! -- Deallocate idm memory
702  call memorystore_remove(this%name, 'NAM', idm_context)
703  call memorystore_remove(component=this%name, context=idm_context)
704 
705  ! -- Internal packages
706  call this%dis%dis_da()
707  call this%fmi%fmi_da()
708  call this%mip%mip_da()
709  call this%budget%budget_da()
710  call this%oc%oc_da()
711  deallocate (this%dis)
712  deallocate (this%fmi)
713  deallocate (this%mip)
714  deallocate (this%budget)
715  deallocate (this%oc)
716 
717  ! -- Method objects
720  call destroy_method_pool()
721 
722  ! -- Boundary packages
723  do ip = 1, this%bndlist%Count()
724  packobj => getbndfromlist(this%bndlist, ip)
725  call packobj%bnd_da()
726  deallocate (packobj)
727  end do
728 
729  ! -- Scalars
730  call mem_deallocate(this%infmi)
731  call mem_deallocate(this%inmip)
732  call mem_deallocate(this%inadv)
733  call mem_deallocate(this%indsp)
734  call mem_deallocate(this%inssm)
735  call mem_deallocate(this%inmst)
736  call mem_deallocate(this%inmvt)
737  call mem_deallocate(this%inoc)
738 
739  ! -- Arrays
740  call mem_deallocate(this%masssto)
741  call mem_deallocate(this%massstoold)
742  call mem_deallocate(this%ratesto)
743 
744  ! -- Track file control
745  deallocate (this%trackfilectl)
746 
747  ! -- Parent type
748  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 193 of file prt.f90.

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

◆ prt_ot()

subroutine prtmodule::prt_ot ( class(prtmodeltype this)

Definition at line 488 of file prt.f90.

489  use tdismodule, only: tdis_ot, endofperiod
490  ! -- dummy
491  class(PrtModelType) :: this
492  ! -- local
493  integer(I4B) :: idvsave
494  integer(I4B) :: idvprint
495  integer(I4B) :: icbcfl
496  integer(I4B) :: icbcun
497  integer(I4B) :: ibudfl
498  integer(I4B) :: ipflag
499 
500  ! -- Note: particle tracking output is handled elsewhere
501 
502  ! -- Set write and print flags
503  idvsave = 0
504  idvprint = 0
505  icbcfl = 0
506  ibudfl = 0
507  if (this%oc%oc_save('CONCENTRATION')) idvsave = 1
508  if (this%oc%oc_print('CONCENTRATION')) idvprint = 1
509  if (this%oc%oc_save('BUDGET')) icbcfl = 1
510  if (this%oc%oc_print('BUDGET')) ibudfl = 1
511  icbcun = this%oc%oc_save_unit('BUDGET')
512 
513  ! -- Override ibudfl and idvprint flags for nonconvergence
514  ! and end of period
515  ibudfl = this%oc%set_print_flag('BUDGET', 1, endofperiod)
516  idvprint = this%oc%set_print_flag('CONCENTRATION', 1, endofperiod)
517 
518  ! -- Save and print flows
519  call this%prt_ot_flow(icbcfl, ibudfl, icbcun)
520 
521  ! -- Save and print dependent variables
522  call this%prt_ot_dv(idvsave, idvprint, ipflag)
523 
524  ! -- Print budget summaries
525  call this%prt_ot_bdsummary(ibudfl, ipflag)
526 
527  ! -- Timing Output; if any dependent variables or budgets
528  ! are printed, then ipflag is set to 1.
529  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 657 of file prt.f90.

658  ! -- modules
659  use tdismodule, only: kstp, kper, totim, delt
660  ! -- dummy
661  class(PrtModelType) :: this
662  integer(I4B), intent(in) :: ibudfl
663  integer(I4B), intent(inout) :: ipflag
664  ! -- local
665  class(BndType), pointer :: packobj
666  integer(I4B) :: ip
667 
668  ! -- Package budget summary
669  do ip = 1, this%bndlist%Count()
670  packobj => getbndfromlist(this%bndlist, ip)
671  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
672  end do
673 
674  ! -- model budget summary
675  call this%budget%finalize_step(delt)
676  if (ibudfl /= 0) then
677  ipflag = 1
678  ! -- model budget summary
679  call this%budget%budget_ot(kstp, kper, this%iout)
680  end if
681 
682  ! -- Write to budget csv
683  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 636 of file prt.f90.

637  ! -- dummy
638  class(PrtModelType) :: this
639  integer(I4B), intent(in) :: idvsave
640  integer(I4B), intent(in) :: idvprint
641  integer(I4B), intent(inout) :: ipflag
642  ! -- local
643  class(BndType), pointer :: packobj
644  integer(I4B) :: ip
645 
646  ! -- Print advanced package dependent variables
647  do ip = 1, this%bndlist%Count()
648  packobj => getbndfromlist(this%bndlist, ip)
649  call packobj%bnd_ot_dv(idvsave, idvprint)
650  end do
651 
652  ! -- save head and print head
653  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 533 of file prt.f90.

534  use prtprpmodule, only: prtprptype
535  class(PrtModelType) :: this
536  integer(I4B), intent(in) :: icbcfl
537  integer(I4B), intent(in) :: ibudfl
538  integer(I4B), intent(in) :: icbcun
539  class(BndType), pointer :: packobj
540  integer(I4B) :: ip
541 
542  ! -- Save PRT flows
543  call this%prt_ot_saveflow(this%dis%nja, this%flowja, icbcfl, icbcun)
544  do ip = 1, this%bndlist%Count()
545  packobj => getbndfromlist(this%bndlist, ip)
546  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
547  end do
548 
549  ! -- Save advanced package flows
550  do ip = 1, this%bndlist%Count()
551  packobj => getbndfromlist(this%bndlist, ip)
552  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
553  end do
554 
555  ! -- Print PRT flows
556  call this%prt_ot_printflow(ibudfl, this%flowja)
557  do ip = 1, this%bndlist%Count()
558  packobj => getbndfromlist(this%bndlist, ip)
559  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
560  end do
561 
562  ! -- Print advanced package flows
563  do ip = 1, this%bndlist%Count()
564  packobj => getbndfromlist(this%bndlist, ip)
565  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
566  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 597 of file prt.f90.

598  ! -- modules
599  use tdismodule, only: kper, kstp
600  use constantsmodule, only: lenbigline
601  ! -- dummy
602  class(PrtModelType) :: this
603  integer(I4B), intent(in) :: ibudfl
604  real(DP), intent(inout), dimension(:) :: flowja
605  ! -- local
606  character(len=LENBIGLINE) :: line
607  character(len=30) :: tempstr
608  integer(I4B) :: n, ipos, m
609  real(DP) :: qnm
610  ! -- formats
611  character(len=*), parameter :: fmtiprflow = &
612  "(/,4x,'CALCULATED INTERCELL FLOW &
613  &FOR PERIOD ', i0, ' STEP ', i0)"
614 
615  ! -- Write flowja to list file if requested
616  if (ibudfl /= 0 .and. this%iprflow > 0) then
617  write (this%iout, fmtiprflow) kper, kstp
618  do n = 1, this%dis%nodes
619  line = ''
620  call this%dis%noder_to_string(n, tempstr)
621  line = trim(tempstr)//':'
622  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
623  m = this%dis%con%ja(ipos)
624  call this%dis%noder_to_string(m, tempstr)
625  line = trim(line)//' '//trim(tempstr)
626  qnm = flowja(ipos)
627  write (tempstr, '(1pg15.6)') qnm
628  line = trim(line)//' '//trim(adjustl(tempstr))
629  end do
630  write (this%iout, '(a)') trim(line)
631  end do
632  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 570 of file prt.f90.

571  ! -- dummy
572  class(PrtModelType) :: this
573  integer(I4B), intent(in) :: nja
574  real(DP), dimension(nja), intent(in) :: flowja
575  integer(I4B), intent(in) :: icbcfl
576  integer(I4B), intent(in) :: icbcun
577  ! -- local
578  integer(I4B) :: ibinun
579 
580  ! -- Set unit number for binary output
581  if (this%ipakcb < 0) then
582  ibinun = icbcun
583  elseif (this%ipakcb == 0) then
584  ibinun = 0
585  else
586  ibinun = this%ipakcb
587  end if
588  if (icbcfl == 0) ibinun = 0
589 
590  ! -- Write the face flows if requested
591  if (ibinun /= 0) then
592  call this%dis%record_connection_array(flowja, ibinun, this%iout)
593  end if

◆ prt_rp()

subroutine prtmodule::prt_rp ( class(prtmodeltype this)

Definition at line 297 of file prt.f90.

298  use tdismodule, only: readnewdata
299  ! -- dummy
300  class(PrtModelType) :: this
301  ! -- local
302  class(BndType), pointer :: packobj
303  integer(I4B) :: ip
304 
305  ! -- Check with TDIS on whether or not it is time to RP
306  if (.not. readnewdata) return
307 
308  ! -- Read and prepare
309  if (this%inoc > 0) call this%oc%oc_rp()
310  do ip = 1, this%bndlist%Count()
311  packobj => getbndfromlist(this%bndlist, ip)
312  call packobj%bnd_rp()
313  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 888 of file prt.f90.

889  ! -- modules
890  use tdismodule, only: kper, kstp, totimc, nper, nstp, delt
891  use prtprpmodule, only: prtprptype
892  ! -- dummy variables
893  class(PrtModelType) :: this
894  ! -- local variables
895  integer(I4B) :: np, ip
896  class(BndType), pointer :: packobj
897  type(ParticleType), pointer :: particle
898  real(DP) :: tmax
899  integer(I4B) :: iprp
900 
901  ! -- Initialize particle
902  call create_particle(particle)
903 
904  ! -- Loop over PRP packages
905  iprp = 0
906  do ip = 1, this%bndlist%Count()
907  packobj => getbndfromlist(this%bndlist, ip)
908  select type (packobj)
909  type is (prtprptype)
910  ! -- Update PRP index
911  iprp = iprp + 1
912 
913  ! -- Initialize PRP-specific track files, if enabled
914  if (packobj%itrkout > 0) then
915  call this%trackfilectl%init_track_file( &
916  packobj%itrkout, &
917  iprp=iprp)
918  end if
919  if (packobj%itrkcsv > 0) then
920  call this%trackfilectl%init_track_file( &
921  packobj%itrkcsv, &
922  csv=.true., &
923  iprp=iprp)
924  end if
925 
926  ! -- Loop over particles in package
927  do np = 1, packobj%nparticles
928  ! -- Load particle from storage
929  call particle%load_particle(packobj%particles, &
930  this%id, iprp, np)
931 
932  ! -- If particle is permanently unreleased, record its initial/terminal state
933  if (particle%istatus == 8) &
934  call this%method%save(particle, reason=3) ! reason=3: termination
935 
936  ! If particle is inactive or not yet to be released, cycle
937  if (particle%istatus > 1) cycle
938 
939  ! If particle released this time step, record its initial state
940  particle%istatus = 1
941  if (particle%trelease >= totimc) &
942  call this%method%save(particle, reason=0) ! reason=0: release
943 
944  ! Maximum time is the end of the time step or the particle
945  ! stop time, whichever comes first, unless it's the final
946  ! time step and the extend option is on, in which case
947  ! it's just the particle stop time.
948  if (nper == kper .and. &
949  nstp(kper) == kstp .and. &
950  particle%iextend > 0) then
951  tmax = particle%tstop
952  else
953  tmax = min(totimc + delt, particle%tstop)
954  end if
955 
956  ! Get and apply the tracking method
957  call this%method%apply(particle, tmax)
958 
959  ! Update particle storage
960  call packobj%particles%save_particle(particle, np)
961  end do
962  end select
963  end do
964 
965  ! -- Deallocate particle
966  deallocate (particle)
integer(i4b), dimension(:), pointer, public, contiguous nstp
number of time steps in each stress period
Definition: tdis.f90:39
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21
Here is the call graph for this function:

Variable Documentation

◆ budtxt

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

Definition at line 36 of file prt.f90.

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

◆ nbditems

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

Definition at line 35 of file prt.f90.

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

◆ niunit_prt

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

Definition at line 113 of file prt.f90.

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

◆ prt_basepkg

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

Definition at line 94 of file prt.f90.

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

◆ prt_multipkg

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

Definition at line 108 of file prt.f90.

108  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 93 of file prt.f90.

93  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 107 of file prt.f90.

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