MODFLOW 6  version 6.8.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_budterms (this)
 Calculate particle mass budget terms. 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 = 2
 
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 916 of file prt.f90.

918  class(PrtModelType) :: this
919  integer(I4B) :: n
920 
921  ! Allocate arrays in parent type (ibound, flowja, nja)
922  call this%ExplicitModelType%allocate_arrays()
923 
924  ! Allocate and initialize PRT-specific arrays
925  call mem_allocate(this%masssto, this%dis%nodes, &
926  'MASSSTO', this%memoryPath)
927  call mem_allocate(this%massstoold, this%dis%nodes, &
928  'MASSSTOOLD', this%memoryPath)
929  call mem_allocate(this%ratesto, this%dis%nodes, &
930  'RATESTO', this%memoryPath)
931  call mem_allocate(this%masstrm, this%dis%nodes, &
932  'MASSTRM', this%memoryPath)
933  call mem_allocate(this%ratetrm, this%dis%nodes, &
934  'RATETRM', this%memoryPath)
935  do n = 1, this%dis%nodes
936  this%masssto(n) = dzero
937  this%massstoold(n) = dzero
938  this%ratesto(n) = dzero
939  this%masstrm(n) = dzero
940  this%ratetrm(n) = dzero
941  end do

◆ allocate_scalars()

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

Definition at line 887 of file prt.f90.

888  ! dummy
889  class(PrtModelType) :: this
890  character(len=*), intent(in) :: modelname
891 
892  ! allocate members from parent class
893  call this%ExplicitModelType%allocate_scalars(modelname)
894 
895  ! allocate members that are part of model class
896  call mem_allocate(this%infmi, 'INFMI', this%memoryPath)
897  call mem_allocate(this%inmip, 'INMIP', this%memoryPath)
898  call mem_allocate(this%inmvt, 'INMVT', this%memoryPath)
899  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
900  call mem_allocate(this%inadv, 'INADV', this%memoryPath)
901  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
902  call mem_allocate(this%inssm, 'INSSM', this%memoryPath)
903  call mem_allocate(this%inoc, 'INOC ', this%memoryPath)
904 
905  this%infmi = 0
906  this%inmip = 0
907  this%inmvt = 0
908  this%inmst = 0
909  this%inadv = 0
910  this%indsp = 0
911  this%inssm = 0
912  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 1096 of file prt.f90.

1098  ! modules
1101  ! dummy
1102  class(PrtModelType) :: this
1103  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1104  type(CharacterStringType), dimension(:), contiguous, &
1105  pointer, intent(inout) :: pkgtypes
1106  type(CharacterStringType), dimension(:), contiguous, &
1107  pointer, intent(inout) :: pkgnames
1108  type(CharacterStringType), dimension(:), contiguous, &
1109  pointer, intent(inout) :: mempaths
1110  integer(I4B), dimension(:), contiguous, &
1111  pointer, intent(inout) :: inunits
1112  ! local
1113  integer(I4B) :: ipakid, ipaknum
1114  character(len=LENFTYPE) :: pkgtype, bndptype
1115  character(len=LENPACKAGENAME) :: pkgname
1116  character(len=LENMEMPATH) :: mempath
1117  integer(I4B), pointer :: inunit
1118  integer(I4B) :: n
1119 
1120  if (allocated(bndpkgs)) then
1121  ! create stress packages
1122  ipakid = 1
1123  bndptype = ''
1124  do n = 1, size(bndpkgs)
1125  pkgtype = pkgtypes(bndpkgs(n))
1126  pkgname = pkgnames(bndpkgs(n))
1127  mempath = mempaths(bndpkgs(n))
1128  inunit => inunits(bndpkgs(n))
1129 
1130  if (bndptype /= pkgtype) then
1131  ipaknum = 1
1132  bndptype = pkgtype
1133  end if
1134 
1135  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1136  inunit, this%iout)
1137  ipakid = ipakid + 1
1138  ipaknum = ipaknum + 1
1139  end do
1140 
1141  ! cleanup
1142  deallocate (bndpkgs)
1143  end if
1144 
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 1148 of file prt.f90.

1149  ! modules
1152  use arrayhandlersmodule, only: expandarray
1153  use memorymanagermodule, only: mem_setptr
1155  use simvariablesmodule, only: idm_context
1156  use budgetmodule, only: budget_cr
1160  use prtmipmodule, only: mip_cr
1161  use prtfmimodule, only: fmi_cr
1162  use prtocmodule, only: oc_cr
1163  ! dummy
1164  class(PrtModelType) :: this
1165  ! local
1166  type(CharacterStringType), dimension(:), contiguous, &
1167  pointer :: pkgtypes => null()
1168  type(CharacterStringType), dimension(:), contiguous, &
1169  pointer :: pkgnames => null()
1170  type(CharacterStringType), dimension(:), contiguous, &
1171  pointer :: mempaths => null()
1172  integer(I4B), dimension(:), contiguous, &
1173  pointer :: inunits => null()
1174  character(len=LENMEMPATH) :: model_mempath
1175  character(len=LENFTYPE) :: pkgtype
1176  character(len=LENPACKAGENAME) :: pkgname
1177  character(len=LENMEMPATH) :: mempath
1178  integer(I4B), pointer :: inunit
1179  integer(I4B), dimension(:), allocatable :: bndpkgs
1180  integer(I4B) :: n
1181  integer(I4B) :: indis = 0 ! DIS enabled flag
1182  character(len=LENMEMPATH) :: mempathmip = ''
1183  character(len=LENMEMPATH) :: mempathfmi = ''
1184  character(len=LENMEMPATH) :: mempathoc = ''
1185 
1186  ! set input memory paths, input/model and input/model/namfile
1187  model_mempath = create_mem_path(component=this%name, context=idm_context)
1188 
1189  ! set pointers to model path package info
1190  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1191  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1192  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1193  call mem_setptr(inunits, 'INUNITS', model_mempath)
1194 
1195  do n = 1, size(pkgtypes)
1196  ! attributes for this input package
1197  pkgtype = pkgtypes(n)
1198  pkgname = pkgnames(n)
1199  mempath = mempaths(n)
1200  inunit => inunits(n)
1201 
1202  ! create dis package first as it is a prerequisite for other packages
1203  select case (pkgtype)
1204  case ('DIS6')
1205  indis = 1
1206  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1207  case ('DISV6')
1208  indis = 1
1209  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1210  case ('DISU6')
1211  indis = 1
1212  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1213  case ('MIP6')
1214  this%inmip = 1
1215  mempathmip = mempath
1216  case ('FMI6')
1217  this%infmi = 1
1218  mempathfmi = mempath
1219  case ('OC6')
1220  this%inoc = 1
1221  mempathoc = mempath
1222  case ('PRP6')
1223  call expandarray(bndpkgs)
1224  bndpkgs(size(bndpkgs)) = n
1225  case default
1226  call pstop(1, "Unrecognized package type: "//pkgtype)
1227  end select
1228  end do
1229 
1230  ! Create budget manager
1231  call budget_cr(this%budget, this%name)
1232 
1233  ! Create tracking method pools
1234  call create_method_pool()
1237 
1238  ! Create packages that are tied directly to model
1239  call mip_cr(this%mip, this%name, mempathmip, this%inmip, this%iout, this%dis)
1240  call fmi_cr(this%fmi, this%name, mempathfmi, this%infmi, this%iout)
1241  call oc_cr(this%oc, this%name, mempathoc, this%inoc, this%iout)
1242 
1243  ! Check to make sure that required ftype's have been specified
1244  call this%ftype_check(indis)
1245 
1246  ! Create boundary packages
1247  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:51
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, input_mempath, 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 994 of file prt.f90.

995  ! dummy
996  class(PrtModelType) :: this
997  integer(I4B), intent(in) :: indis
998  ! local
999  character(len=LINELENGTH) :: errmsg
1000 
1001  ! Check for DIS(u) and MIP. Stop if not present.
1002  if (indis == 0) then
1003  write (errmsg, '(1x,a)') &
1004  'Discretization (DIS6, DISV6, or DISU6) package not specified.'
1005  call store_error(errmsg)
1006  end if
1007  if (this%inmip == 0) then
1008  write (errmsg, '(1x,a)') &
1009  'Model input (MIP6) package not specified.'
1010  call store_error(errmsg)
1011  end if
1012 
1013  if (count_errors() > 0) then
1014  write (errmsg, '(1x,a)') 'One or more required package(s) not specified.'
1015  call store_error(errmsg)
1016  call store_error_filename(this%filename)
1017  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 1251 of file prt.f90.

1253  class(PrtModelType) :: this
1254  type(PrtNamParamFoundType), intent(in) :: found
1255 
1256  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1257 
1258  if (found%print_input) then
1259  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1260  'FOR ALL MODEL STRESS PACKAGES'
1261  end if
1262 
1263  if (found%print_flows) then
1264  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1265  'FOR ALL MODEL PACKAGES'
1266  end if
1267 
1268  if (found%save_flows) then
1269  write (this%iout, '(4x,a)') &
1270  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1271  end if
1272 
1273  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 945 of file prt.f90.

947  ! modules
948  use constantsmodule, only: linelength
949  use prtprpmodule, only: prp_create
950  use apimodule, only: api_create
951  ! dummy
952  class(PrtModelType) :: this
953  character(len=*), intent(in) :: filtyp
954  character(len=LINELENGTH) :: errmsg
955  integer(I4B), intent(in) :: ipakid
956  integer(I4B), intent(in) :: ipaknum
957  character(len=*), intent(in) :: pakname
958  character(len=*), intent(in) :: mempath
959  integer(I4B), intent(in) :: inunit
960  integer(I4B), intent(in) :: iout
961  ! local
962  class(BndType), pointer :: packobj
963  class(BndType), pointer :: packobj2
964  integer(I4B) :: ip
965 
966  ! This part creates the package object
967  select case (filtyp)
968  case ('PRP6')
969  call prp_create(packobj, ipakid, ipaknum, inunit, iout, &
970  this%name, pakname, mempath, this%fmi)
971  case ('API6')
972  call api_create(packobj, ipakid, ipaknum, inunit, iout, &
973  this%name, pakname, mempath)
974  case default
975  write (errmsg, *) 'Invalid package type: ', filtyp
976  call store_error(errmsg, terminate=.true.)
977  end select
978 
979  ! Packages is the bndlist that is associated with the parent model
980  ! The following statement puts a pointer to this package in the ipakid
981  ! position of packages.
982  do ip = 1, this%bndlist%Count()
983  packobj2 => getbndfromlist(this%bndlist, ip)
984  if (packobj2%packName == pakname) then
985  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
986  'already exists: ', trim(pakname)
987  call store_error(errmsg, terminate=.true.)
988  end if
989  end do
990  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:106
Here is the call graph for this function:

◆ prt_ad()

subroutine prtmodule::prt_ad ( class(prtmodeltype this)

Definition at line 355 of file prt.f90.

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

241  ! modules
242  use constantsmodule, only: dhnoflo
243  use prtprpmodule, only: prtprptype
244  use prtmipmodule, only: prtmiptype
246  ! dummy
247  class(PrtModelType) :: this
248  ! locals
249  integer(I4B) :: ip, nprp
250  class(BndType), pointer :: packobj
251 
252  ! Set up basic packages
253  call this%fmi%fmi_ar(this%ibound)
254  if (this%inmip > 0) call this%mip%mip_ar()
255 
256  ! Set up output control and budget
257  call this%oc%oc_ar(this%dis, dhnoflo)
258  call this%budget%set_ibudcsv(this%oc%ibudcsv)
259 
260  ! Select tracking events
261  call this%tracks%select_events( &
262  this%oc%trackrelease, &
263  this%oc%trackfeatexit, &
264  this%oc%tracktimestep, &
265  this%oc%trackterminate, &
266  this%oc%trackweaksink, &
267  this%oc%trackusertime, &
268  this%oc%tracksubfexit, &
269  this%oc%trackdropped)
270 
271  ! Set up boundary pkgs and pkg-scoped track files
272  nprp = 0
273  do ip = 1, this%bndlist%Count()
274  packobj => getbndfromlist(this%bndlist, ip)
275  select type (packobj)
276  type is (prtprptype)
277  nprp = nprp + 1
278  call packobj%prp_set_pointers(this%ibound, this%mip%izone)
279  call packobj%bnd_ar()
280  call packobj%bnd_ar()
281  if (packobj%itrkout > 0) then
282  call this%tracks%init_file( &
283  packobj%itrkout, &
284  iprp=nprp)
285  end if
286  if (packobj%itrkcsv > 0) then
287  call this%tracks%init_file( &
288  packobj%itrkcsv, &
289  csv=.true., &
290  iprp=nprp)
291  end if
292  class default
293  call packobj%bnd_ar()
294  end select
295  end do
296 
297  ! Set up model-scoped track files
298  if (this%oc%itrkout > 0) &
299  call this%tracks%init_file(this%oc%itrkout)
300  if (this%oc%itrkcsv > 0) &
301  call this%tracks%init_file(this%oc%itrkcsv, csv=.true.)
302 
303  ! Set up the tracking method
304  select type (dis => this%dis)
305  type is (distype)
306  call method_dis%init( &
307  fmi=this%fmi, &
308  events=this%events, &
309  izone=this%mip%izone, &
310  flowja=this%flowja, &
311  porosity=this%mip%porosity, &
312  retfactor=this%mip%retfactor, &
313  tracktimes=this%oc%tracktimes)
314  this%method => method_dis
315  type is (disvtype)
316  call method_disv%init( &
317  fmi=this%fmi, &
318  events=this%events, &
319  izone=this%mip%izone, &
320  flowja=this%flowja, &
321  porosity=this%mip%porosity, &
322  retfactor=this%mip%retfactor, &
323  tracktimes=this%oc%tracktimes)
324  this%method => method_disv
325  end select
326 
327  ! Subscribe track output manager to events
328  call this%events%subscribe(this%tracks)
329 
330  ! Set verbose tracing if requested
331  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:39
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 523 of file prt.f90.

524  ! modules
525  use tdismodule, only: delt
526  use budgetmodule, only: rate_accumulator
527  ! dummy
528  class(PrtModelType) :: this
529  integer(I4B), intent(in) :: icnvg
530  integer(I4B), intent(in) :: isuppress_output
531  ! local
532  integer(I4B) :: ip
533  class(BndType), pointer :: packobj
534  real(DP) :: rin
535  real(DP) :: rout
536 
537  ! Budget routines (start by resetting). Sole purpose of this section
538  ! is to add in and outs to model budget. All ins and out for a model
539  ! should be added here to this%budget. In a subsequent exchange call,
540  ! exchange flows might also be added.
541  call this%budget%reset()
542  ! storage term
543  call rate_accumulator(this%ratesto, rin, rout)
544  call this%budget%addentry(rin, rout, delt, budtxt(1), &
545  isuppress_output, ' PRT')
546  ! termination term
547  call rate_accumulator(this%ratetrm, rin, rout)
548  call this%budget%addentry(rin, rout, delt, budtxt(2), &
549  isuppress_output, ' PRT')
550  ! boundary packages
551  do ip = 1, this%bndlist%Count()
552  packobj => getbndfromlist(this%bndlist, ip)
553  call packobj%bnd_bd(this%budget)
554  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 398 of file prt.f90.

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

subroutine prtmodule::prt_cq_budterms ( class(prtmodeltype this)

Definition at line 444 of file prt.f90.

445  ! modules
446  use tdismodule, only: delt
447  use prtprpmodule, only: prtprptype
448  ! dummy
449  class(PrtModelType) :: this
450  ! local
451  integer(I4B) :: ip
452  class(BndType), pointer :: packobj
453  integer(I4B) :: n
454  integer(I4B) :: np
455  integer(I4B) :: idiag
456  integer(I4B) :: iprp
457  integer(I4B) :: istatus
458  real(DP) :: tled
459  real(DP) :: ratesto, ratetrm
460  character(len=:), allocatable :: particle_id
461  type(ParticleType), pointer :: particle
462 
463  call create_particle(particle)
464 
465  ! Reciprocal of time step size.
466  tled = done / delt
467 
468  ! Reset mass and rate arrays
469  do n = 1, this%dis%nodes
470  this%masssto(n) = dzero
471  this%masstrm(n) = dzero
472  this%ratesto(n) = dzero
473  this%ratetrm(n) = dzero
474  end do
475 
476  ! Loop over PRP packages and assign particle mass to the
477  ! appropriate budget term based on the particle status.
478  iprp = 0
479  do ip = 1, this%bndlist%Count()
480  packobj => getbndfromlist(this%bndlist, ip)
481  select type (packobj)
482  type is (prtprptype)
483  do np = 1, packobj%nparticles
484  call packobj%particles%get(particle, this%id, iprp, np)
485  istatus = packobj%particles%istatus(np)
486  particle_id = particle%get_id()
487  if (istatus == active) then
488  ! calculate storage mass
489  n = packobj%particles%itrdomain(np, level_feature)
490  this%masssto(n) = this%masssto(n) + done ! unit mass
491  else if (istatus > active) then
492  if (this%trm_ids%get(particle_id) /= 0) cycle
493  ! calculate terminating mass
494  n = packobj%particles%itrdomain(np, level_feature)
495  this%masstrm(n) = this%masstrm(n) + done ! unit mass
496  call this%trm_ids%add(particle_id, 1) ! mark id terminated
497  end if
498  end do
499  end select
500  end do
501 
502  ! Calculate rates and update flowja
503  do n = 1, this%dis%nodes
504  ratesto = -(this%masssto(n) - this%massstoold(n)) * tled
505  ratetrm = -this%masstrm(n) * tled
506  this%ratesto(n) = ratesto
507  this%ratetrm(n) = ratetrm
508  idiag = this%dis%con%ia(n)
509  this%flowja(idiag) = this%flowja(idiag) + ratesto
510  end do
511 
512  call particle%destroy()
513  deallocate (particle)
514 
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 127 of file prt.f90.

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

821  ! modules
828  ! dummy
829  class(PrtModelType) :: this
830  ! local
831  integer(I4B) :: ip
832  class(BndType), pointer :: packobj
833 
834  ! Deallocate idm memory
835  call memorystore_remove(this%name, 'NAM', idm_context)
836  call memorystore_remove(component=this%name, context=idm_context)
837 
838  ! Internal packages
839  call this%dis%dis_da()
840  call this%fmi%fmi_da()
841  call this%mip%mip_da()
842  call this%budget%budget_da()
843  call this%oc%oc_da()
844  deallocate (this%dis)
845  deallocate (this%fmi)
846  deallocate (this%mip)
847  deallocate (this%budget)
848  deallocate (this%oc)
849 
850  ! Method objects
853  call destroy_method_pool()
854 
855  ! Boundary packages
856  do ip = 1, this%bndlist%Count()
857  packobj => getbndfromlist(this%bndlist, ip)
858  call packobj%bnd_da()
859  deallocate (packobj)
860  end do
861 
862  ! Scalars
863  call mem_deallocate(this%infmi)
864  call mem_deallocate(this%inmip)
865  call mem_deallocate(this%inadv)
866  call mem_deallocate(this%indsp)
867  call mem_deallocate(this%inssm)
868  call mem_deallocate(this%inmst)
869  call mem_deallocate(this%inmvt)
870  call mem_deallocate(this%inoc)
871 
872  ! Arrays
873  call mem_deallocate(this%masssto)
874  call mem_deallocate(this%massstoold)
875  call mem_deallocate(this%ratesto)
876  call mem_deallocate(this%masstrm)
877  call mem_deallocate(this%ratetrm)
878 
879  call this%tracks%destroy()
880  deallocate (this%events)
881  deallocate (this%tracks)
882 
883  call this%ExplicitModelType%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 207 of file prt.f90.

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

◆ prt_ot()

subroutine prtmodule::prt_ot ( class(prtmodeltype this)

Definition at line 558 of file prt.f90.

559  use tdismodule, only: tdis_ot, endofperiod
560  ! dummy
561  class(PrtModelType) :: this
562  ! local
563  integer(I4B) :: idvsave
564  integer(I4B) :: idvprint
565  integer(I4B) :: icbcfl
566  integer(I4B) :: icbcun
567  integer(I4B) :: ibudfl
568  integer(I4B) :: ipflag
569 
570  ! Note: particle tracking output is handled elsewhere
571 
572  ! Set write and print flags
573  idvsave = 0
574  idvprint = 0
575  icbcfl = 0
576  ibudfl = 0
577  if (this%oc%oc_save('CONCENTRATION')) idvsave = 1
578  if (this%oc%oc_print('CONCENTRATION')) idvprint = 1
579  if (this%oc%oc_save('BUDGET')) icbcfl = 1
580  if (this%oc%oc_print('BUDGET')) ibudfl = 1
581  icbcun = this%oc%oc_save_unit('BUDGET')
582 
583  ! Override ibudfl and idvprint flags for nonconvergence
584  ! and end of period
585  ibudfl = this%oc%set_print_flag('BUDGET', 1, endofperiod)
586  idvprint = this%oc%set_print_flag('CONCENTRATION', 1, endofperiod)
587 
588  ! Save and print flows
589  call this%prt_ot_flow(icbcfl, ibudfl, icbcun)
590 
591  ! Save and print dependent variables
592  call this%prt_ot_dv(idvsave, idvprint, ipflag)
593 
594  ! Print budget summaries
595  call this%prt_ot_bdsummary(ibudfl, ipflag)
596 
597  ! Timing Output; if any dependent variables or budgets
598  ! are printed, then ipflag is set to 1.
599  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 790 of file prt.f90.

791  ! modules
792  use tdismodule, only: kstp, kper, totim, delt
793  ! dummy
794  class(PrtModelType) :: this
795  integer(I4B), intent(in) :: ibudfl
796  integer(I4B), intent(inout) :: ipflag
797  ! local
798  class(BndType), pointer :: packobj
799  integer(I4B) :: ip
800 
801  ! Package budget summary
802  do ip = 1, this%bndlist%Count()
803  packobj => getbndfromlist(this%bndlist, ip)
804  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
805  end do
806 
807  ! model budget summary
808  call this%budget%finalize_step(delt)
809  if (ibudfl /= 0) then
810  ipflag = 1
811  ! model budget summary
812  call this%budget%budget_ot(kstp, kper, this%iout)
813  end if
814 
815  ! Write to budget csv
816  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 769 of file prt.f90.

770  ! dummy
771  class(PrtModelType) :: this
772  integer(I4B), intent(in) :: idvsave
773  integer(I4B), intent(in) :: idvprint
774  integer(I4B), intent(inout) :: ipflag
775  ! local
776  class(BndType), pointer :: packobj
777  integer(I4B) :: ip
778 
779  ! Print advanced package dependent variables
780  do ip = 1, this%bndlist%Count()
781  packobj => getbndfromlist(this%bndlist, ip)
782  call packobj%bnd_ot_dv(idvsave, idvprint)
783  end do
784 
785  ! save head and print head
786  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 603 of file prt.f90.

604  use prtprpmodule, only: prtprptype
605  class(PrtModelType) :: this
606  integer(I4B), intent(in) :: icbcfl
607  integer(I4B), intent(in) :: ibudfl
608  integer(I4B), intent(in) :: icbcun
609  class(BndType), pointer :: packobj
610  integer(I4B) :: ip
611 
612  ! Save PRT flows
613  call this%prt_ot_saveflow(this%dis%nja, this%flowja, icbcfl, icbcun)
614  do ip = 1, this%bndlist%Count()
615  packobj => getbndfromlist(this%bndlist, ip)
616  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
617  end do
618 
619  ! Save advanced package flows
620  do ip = 1, this%bndlist%Count()
621  packobj => getbndfromlist(this%bndlist, ip)
622  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
623  end do
624 
625  ! Print PRT flows
626  call this%prt_ot_printflow(ibudfl, this%flowja)
627  do ip = 1, this%bndlist%Count()
628  packobj => getbndfromlist(this%bndlist, ip)
629  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
630  end do
631 
632  ! Print advanced package flows
633  do ip = 1, this%bndlist%Count()
634  packobj => getbndfromlist(this%bndlist, ip)
635  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
636  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 730 of file prt.f90.

731  ! modules
732  use tdismodule, only: kper, kstp
733  use constantsmodule, only: lenbigline
734  ! dummy
735  class(PrtModelType) :: this
736  integer(I4B), intent(in) :: ibudfl
737  real(DP), intent(inout), dimension(:) :: flowja
738  ! local
739  character(len=LENBIGLINE) :: line
740  character(len=30) :: tempstr
741  integer(I4B) :: n, ipos, m
742  real(DP) :: qnm
743  ! formats
744  character(len=*), parameter :: fmtiprflow = &
745  "(/,4x,'CALCULATED INTERCELL FLOW &
746  &FOR PERIOD ', i0, ' STEP ', i0)"
747 
748  ! Write flowja to list file if requested
749  if (ibudfl /= 0 .and. this%iprflow > 0) then
750  write (this%iout, fmtiprflow) kper, kstp
751  do n = 1, this%dis%nodes
752  line = ''
753  call this%dis%noder_to_string(n, tempstr)
754  line = trim(tempstr)//':'
755  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
756  m = this%dis%con%ja(ipos)
757  call this%dis%noder_to_string(m, tempstr)
758  line = trim(line)//' '//trim(tempstr)
759  qnm = flowja(ipos)
760  write (tempstr, '(1pg15.6)') qnm
761  line = trim(line)//' '//trim(adjustl(tempstr))
762  end do
763  write (this%iout, '(a)') trim(line)
764  end do
765  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 640 of file prt.f90.

641  ! dummy
642  class(PrtModelType) :: this
643  integer(I4B), intent(in) :: nja
644  real(DP), dimension(nja), intent(in) :: flowja
645  integer(I4B), intent(in) :: icbcfl
646  integer(I4B), intent(in) :: icbcun
647  ! local
648  integer(I4B) :: ibinun
649  integer(I4B) :: naux
650  real(DP), dimension(0) :: auxrow
651  character(len=LENAUXNAME), dimension(0) :: auxname
652  logical(LGP) :: header_written
653  integer(I4B) :: i, nn
654  real(DP) :: m
655  integer(I4B) :: nsto, ntrm
656  logical(LGP), allocatable :: msto_mask(:), mtrm_mask(:)
657  integer(I4B), allocatable :: msto_nns(:), mtrm_nns(:)
658  real(DP), allocatable :: msto_vals(:), mtrm_vals(:)
659 
660  ! Set unit number for binary output
661  if (this%ipakcb < 0) then
662  ibinun = icbcun
663  elseif (this%ipakcb == 0) then
664  ibinun = 0
665  else
666  ibinun = this%ipakcb
667  end if
668  if (icbcfl == 0) ibinun = 0
669 
670  ! Return if nothing to do
671  if (ibinun == 0) return
672 
673  ! Write mass face flows
674  call this%dis%record_connection_array(flowja, ibinun, this%iout)
675 
676  ! Write mass storage term
677  naux = 0
678  header_written = .false.
679  msto_mask = this%masssto > dzero
680  msto_vals = pack(this%masssto, msto_mask)
681  msto_nns = [(i, i=1, size(this%masssto))]
682  msto_nns = pack(msto_nns, msto_mask)
683  nsto = size(msto_nns)
684  do i = 1, nsto
685  nn = msto_nns(i)
686  m = msto_vals(i)
687  if (.not. header_written) then
688  call this%dis%record_srcdst_list_header(budtxt(1), &
689  'PRT ', &
690  'PRT ', &
691  'PRT ', &
692  'STORAGE ', &
693  naux, auxname, ibinun, &
694  nsto, this%iout)
695  header_written = .true.
696  end if
697  call this%dis%record_mf6_list_entry(ibinun, nn, nn, m, &
698  0, auxrow, &
699  olconv2=.false.)
700  end do
701 
702  ! Write mass termination term
703  header_written = .false.
704  mtrm_mask = this%masstrm > dzero
705  mtrm_vals = pack(this%masstrm, mtrm_mask)
706  mtrm_nns = [(i, i=1, size(this%masstrm))]
707  mtrm_nns = pack(mtrm_nns, mtrm_mask)
708  ntrm = size(mtrm_nns)
709  do i = 1, ntrm
710  nn = mtrm_nns(i)
711  m = mtrm_vals(i)
712  if (.not. header_written) then
713  call this%dis%record_srcdst_list_header(budtxt(2), &
714  'PRT ', &
715  'PRT ', &
716  'PRT ', &
717  'TERMINATION ', &
718  naux, auxname, ibinun, &
719  ntrm, this%iout)
720  header_written = .true.
721  end if
722  call this%dis%record_mf6_list_entry(ibinun, nn, nn, m, &
723  0, auxrow, &
724  olconv2=.false.)
725  end do
726 

◆ prt_rp()

subroutine prtmodule::prt_rp ( class(prtmodeltype this)

Definition at line 335 of file prt.f90.

336  use tdismodule, only: readnewdata
337  ! dummy
338  class(PrtModelType) :: this
339  ! local
340  class(BndType), pointer :: packobj
341  integer(I4B) :: ip
342 
343  ! Check with TDIS on whether or not it is time to RP
344  if (.not. readnewdata) return
345 
346  ! Read and prepare
347  if (this%inoc > 0) call this%oc%oc_rp()
348  do ip = 1, this%bndlist%Count()
349  packobj => getbndfromlist(this%bndlist, ip)
350  call packobj%bnd_rp()
351  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 1021 of file prt.f90.

1022  use tdismodule, only: totimc, delt, endofsimulation
1023  use prtprpmodule, only: prtprptype
1026  ! dummy
1027  class(PrtModelType) :: this
1028  ! local
1029  integer(I4B) :: np, ip
1030  class(BndType), pointer :: packobj
1031  type(ParticleType), pointer :: particle
1032  real(DP) :: tmax
1033  integer(I4B) :: iprp
1034 
1035  ! A single particle is reused in the tracking loops
1036  ! to avoid allocating and deallocating it each time.
1037  ! get() and put() retrieve and store particle state.
1038  call create_particle(particle)
1039  ! Loop over PRP packages and particles within them.
1040  iprp = 0
1041  do ip = 1, this%bndlist%Count()
1042  packobj => getbndfromlist(this%bndlist, ip)
1043  select type (packobj)
1044  type is (prtprptype)
1045  iprp = iprp + 1
1046  do np = 1, packobj%nparticles
1047  ! Get the particle from the store
1048  call packobj%particles%get(particle, this%id, iprp, np)
1049  ! If particle is permanently unreleased, cycle.
1050  ! Raise a termination event if we haven't yet.
1051  ! TODO: when we have generic dynamic vectors,
1052  ! consider terminating permanently unreleased
1053  ! in PRP instead of here. For now, status -8
1054  ! indicates the permanently unreleased event
1055  ! is not yet recorded, status 8 it has been.
1056  if (particle%istatus == (-1 * term_unreleased)) then
1057  call this%method%terminate(particle, status=term_unreleased)
1058  call packobj%particles%put(particle, np)
1059  end if
1060  if (particle%istatus > active) cycle ! Skip terminated particles
1061  particle%istatus = active ! Set active status in case of release
1062  ! If the particle was released this time step, emit a release event
1063  if (particle%trelease >= totimc) call this%method%release(particle)
1064  ! Maximum time is the end of the time step or the particle
1065  ! stop time, whichever comes first, unless it's the final
1066  ! time step and the extend option is on, in which case
1067  ! it's just the particle stop time.
1068  if (endofsimulation .and. particle%extend) then
1069  tmax = particle%tstop
1070  else
1071  tmax = min(totimc + delt, particle%tstop)
1072  end if
1073  ! Apply the tracking method until the maximum time.
1074  call this%method%apply(particle, tmax)
1075  ! If the particle timed out, terminate it.
1076  ! "Timed out" means it's still active but
1077  ! - it reached its stop time, or
1078  ! - the simulation is over.
1079  ! We can't detect timeout within the tracking
1080  ! method because the method just receives the
1081  ! maximum time with no context on what it is.
1082  ! TODO maybe think about changing that?
1083  if (particle%istatus <= active .and. &
1084  (particle%ttrack == particle%tstop .or. endofsimulation)) &
1085  call this%method%terminate(particle, status=term_timeout)
1086  ! Return the particle to the store
1087  call packobj%particles%put(particle, np)
1088  end do
1089  end select
1090  end do
1091  call particle%destroy()
1092  deallocate (particle)
@, public release
particle was released
@, public terminate
particle terminated
@ term_timeout
terminated at stop time or end of simulation
Definition: Particle.f90:38
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:36
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 41 of file prt.f90.

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

◆ nbditems

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

Definition at line 40 of file prt.f90.

40  integer(I4B), parameter :: NBDITEMS = 2

◆ niunit_prt

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

Definition at line 122 of file prt.f90.

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

◆ prt_basepkg

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

Definition at line 103 of file prt.f90.

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

◆ prt_multipkg

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

Definition at line 117 of file prt.f90.

117  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 102 of file prt.f90.

102  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 116 of file prt.f90.

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