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, isuppress_output)
 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 create_exg_prp (this)
 Create an exchange PRP package for particles entering this model from other model. 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 940 of file prt.f90.

942  class(PrtModelType) :: this
943  integer(I4B) :: n
944 
945  ! Allocate arrays in parent type (ibound, flowja, nja)
946  call this%ExplicitModelType%allocate_arrays()
947 
948  ! Allocate and initialize PRT-specific arrays
949  call mem_allocate(this%masssto, this%dis%nodes, &
950  'MASSSTO', this%memoryPath)
951  call mem_allocate(this%massstoold, this%dis%nodes, &
952  'MASSSTOOLD', this%memoryPath)
953  call mem_allocate(this%ratesto, this%dis%nodes, &
954  'RATESTO', this%memoryPath)
955  call mem_allocate(this%masstrm, this%dis%nodes, &
956  'MASSTRM', this%memoryPath)
957  call mem_allocate(this%ratetrm, this%dis%nodes, &
958  'RATETRM', this%memoryPath)
959  do n = 1, this%dis%nodes
960  this%masssto(n) = dzero
961  this%massstoold(n) = dzero
962  this%ratesto(n) = dzero
963  this%masstrm(n) = dzero
964  this%ratetrm(n) = dzero
965  end do

◆ allocate_scalars()

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

Definition at line 911 of file prt.f90.

912  ! dummy
913  class(PrtModelType) :: this
914  character(len=*), intent(in) :: modelname
915 
916  ! allocate members from parent class
917  call this%ExplicitModelType%allocate_scalars(modelname)
918 
919  ! allocate members that are part of model class
920  call mem_allocate(this%infmi, 'INFMI', this%memoryPath)
921  call mem_allocate(this%inmip, 'INMIP', this%memoryPath)
922  call mem_allocate(this%inmvt, 'INMVT', this%memoryPath)
923  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
924  call mem_allocate(this%inadv, 'INADV', this%memoryPath)
925  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
926  call mem_allocate(this%inssm, 'INSSM', this%memoryPath)
927  call mem_allocate(this%inoc, 'INOC ', this%memoryPath)
928 
929  this%infmi = 0
930  this%inmip = 0
931  this%inmvt = 0
932  this%inmst = 0
933  this%inadv = 0
934  this%indsp = 0
935  this%inssm = 0
936  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 1125 of file prt.f90.

1127  ! modules
1130  ! dummy
1131  class(PrtModelType) :: this
1132  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1133  type(CharacterStringType), dimension(:), contiguous, &
1134  pointer, intent(inout) :: pkgtypes
1135  type(CharacterStringType), dimension(:), contiguous, &
1136  pointer, intent(inout) :: pkgnames
1137  type(CharacterStringType), dimension(:), contiguous, &
1138  pointer, intent(inout) :: mempaths
1139  integer(I4B), dimension(:), contiguous, &
1140  pointer, intent(inout) :: inunits
1141  ! local
1142  integer(I4B) :: ipakid, ipaknum
1143  character(len=LENFTYPE) :: pkgtype, bndptype
1144  character(len=LENPACKAGENAME) :: pkgname
1145  character(len=LENMEMPATH) :: mempath
1146  integer(I4B), pointer :: inunit
1147  integer(I4B) :: n
1148 
1149  if (allocated(bndpkgs)) then
1150  ! create stress packages
1151  ipakid = 1
1152  bndptype = ''
1153  do n = 1, size(bndpkgs)
1154  pkgtype = pkgtypes(bndpkgs(n))
1155  pkgname = pkgnames(bndpkgs(n))
1156  mempath = mempaths(bndpkgs(n))
1157  inunit => inunits(bndpkgs(n))
1158 
1159  if (bndptype /= pkgtype) then
1160  ipaknum = 1
1161  bndptype = pkgtype
1162  end if
1163 
1164  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1165  inunit, this%iout)
1166  ipakid = ipakid + 1
1167  ipaknum = ipaknum + 1
1168  end do
1169 
1170  ! cleanup
1171  deallocate (bndpkgs)
1172  end if
1173 
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_exg_prp()

subroutine prtmodule::create_exg_prp ( class(prtmodeltype this)

Definition at line 1278 of file prt.f90.

1279  class(PrtModelType) :: this
1280  ! local
1281  class(BndType), pointer :: packobj
1282  character(len=LENPACKAGENAME) :: exgprp_name
1283 
1284  exgprp_name = 'EXGPRP'
1285 
1286  call prp_create(packobj, &
1287  id=0, &
1288  ibcnum=0, &
1289  inunit=-1, &
1290  iout=this%iout, &
1291  namemodel=this%name, &
1292  pakname=exgprp_name, &
1293  fmi=this%fmi)
1294  call addbndtolist(this%bndlist, packobj)
Here is the call graph for this function:

◆ create_packages()

subroutine prtmodule::create_packages ( class(prtmodeltype this)

Definition at line 1177 of file prt.f90.

1178  ! modules
1181  use arrayhandlersmodule, only: expandarray
1182  use memorymanagermodule, only: mem_setptr
1184  use simvariablesmodule, only: idm_context
1185  use budgetmodule, only: budget_cr
1186  use prtmipmodule, only: mip_cr
1187  use prtfmimodule, only: fmi_cr
1188  use prtocmodule, only: oc_cr
1189  ! dummy
1190  class(PrtModelType) :: this
1191  ! local
1192  type(CharacterStringType), dimension(:), contiguous, &
1193  pointer :: pkgtypes => null()
1194  type(CharacterStringType), dimension(:), contiguous, &
1195  pointer :: pkgnames => null()
1196  type(CharacterStringType), dimension(:), contiguous, &
1197  pointer :: mempaths => null()
1198  integer(I4B), dimension(:), contiguous, &
1199  pointer :: inunits => null()
1200  character(len=LENMEMPATH) :: model_mempath
1201  character(len=LENFTYPE) :: pkgtype
1202  character(len=LENPACKAGENAME) :: pkgname
1203  character(len=LENMEMPATH) :: mempath
1204  integer(I4B), pointer :: inunit
1205  integer(I4B), dimension(:), allocatable :: bndpkgs
1206  integer(I4B) :: n
1207  integer(I4B) :: indis = 0 ! DIS enabled flag
1208  character(len=LENMEMPATH) :: mempathmip = ''
1209  character(len=LENMEMPATH) :: mempathfmi = ''
1210  character(len=LENMEMPATH) :: mempathoc = ''
1211 
1212  ! set input memory paths, input/model and input/model/namfile
1213  model_mempath = create_mem_path(component=this%name, context=idm_context)
1214 
1215  ! set pointers to model path package info
1216  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1217  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1218  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1219  call mem_setptr(inunits, 'INUNITS', model_mempath)
1220 
1221  ! determine which packages we have. create
1222  ! dis up front as the others depend on it.
1223  do n = 1, size(pkgtypes)
1224  pkgtype = pkgtypes(n)
1225  pkgname = pkgnames(n)
1226  mempath = mempaths(n)
1227  inunit => inunits(n)
1228 
1229  select case (pkgtype)
1230  case ('DIS6')
1231  indis = 1
1232  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1233  case ('DISV6')
1234  indis = 1
1235  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1236  case ('DISU6')
1237  indis = 1
1238  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1239  case ('MIP6')
1240  this%inmip = 1
1241  mempathmip = mempath
1242  case ('FMI6')
1243  this%infmi = 1
1244  mempathfmi = mempath
1245  case ('OC6')
1246  this%inoc = 1
1247  mempathoc = mempath
1248  case ('PRP6')
1249  call expandarray(bndpkgs)
1250  bndpkgs(size(bndpkgs)) = n
1251  case default
1252  call pstop(1, "Unrecognized package type: "//pkgtype)
1253  end select
1254  end do
1255 
1256  ! Create budget manager
1257  call budget_cr(this%budget, this%name)
1258 
1259  ! Create tracking methods
1260  call create_method_dis(this%method_dis)
1261  call create_method_disv(this%method_disv)
1262 
1263  ! Create non-boundary packages
1264  call mip_cr(this%mip, this%name, mempathmip, this%inmip, this%iout, this%dis)
1265  call fmi_cr(this%fmi, this%name, mempathfmi, this%infmi, this%iout)
1266  call oc_cr(this%oc, this%name, mempathoc, this%inoc, this%iout)
1267 
1268  ! Check required input files
1269  call this%ftype_check(indis)
1270 
1271  ! Create boundary packages
1272  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
1273  call this%create_exg_prp()
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
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:53
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 1017 of file prt.f90.

1018  ! dummy
1019  class(PrtModelType) :: this
1020  integer(I4B), intent(in) :: indis
1021  ! local
1022  character(len=LINELENGTH) :: errmsg
1023 
1024  ! Check for DIS(u) and MIP. Stop if not present.
1025  if (indis == 0) then
1026  write (errmsg, '(1x,a)') &
1027  'Discretization (DIS6, DISV6, or DISU6) package not specified.'
1028  call store_error(errmsg)
1029  end if
1030  if (this%inmip == 0) then
1031  write (errmsg, '(1x,a)') &
1032  'Model input (MIP6) package not specified.'
1033  call store_error(errmsg)
1034  end if
1035 
1036  if (count_errors() > 0) then
1037  write (errmsg, '(1x,a)') 'One or more required package(s) not specified.'
1038  call store_error(errmsg)
1039  call store_error_filename(this%filename)
1040  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 
)
private

Definition at line 1298 of file prt.f90.

1300  class(PrtModelType) :: this
1301  type(PrtNamParamFoundType), intent(in) :: found
1302 
1303  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1304 
1305  if (found%print_input) then
1306  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1307  'FOR ALL MODEL STRESS PACKAGES'
1308  end if
1309 
1310  if (found%print_flows) then
1311  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1312  'FOR ALL MODEL PACKAGES'
1313  end if
1314 
1315  if (found%save_flows) then
1316  write (this%iout, '(4x,a)') &
1317  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1318  end if
1319 
1320  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 969 of file prt.f90.

971  ! modules
972  use constantsmodule, only: linelength
973  use apimodule, only: api_create
974  ! dummy
975  class(PrtModelType) :: this
976  character(len=*), intent(in) :: filtyp
977  character(len=LINELENGTH) :: errmsg
978  integer(I4B), intent(in) :: ipakid
979  integer(I4B), intent(in) :: ipaknum
980  character(len=*), intent(in) :: pakname
981  character(len=*), intent(in) :: mempath
982  integer(I4B), intent(in) :: inunit
983  integer(I4B), intent(in) :: iout
984  ! local
985  class(BndType), pointer :: packobj
986  class(BndType), pointer :: packobj2
987  integer(I4B) :: ip
988 
989  ! This part creates the package object
990  select case (filtyp)
991  case ('PRP6')
992  call prp_create(packobj, ipakid, ipaknum, inunit, iout, &
993  this%name, pakname, this%fmi, mempath)
994  case ('API6')
995  call api_create(packobj, ipakid, ipaknum, inunit, iout, &
996  this%name, pakname, mempath)
997  case default
998  write (errmsg, *) 'Invalid package type: ', filtyp
999  call store_error(errmsg, terminate=.true.)
1000  end select
1001 
1002  ! Packages is the bndlist that is associated with the parent model
1003  ! The following statement puts a pointer to this package in the ipakid
1004  ! position of packages.
1005  do ip = 1, this%bndlist%Count()
1006  packobj2 => getbndfromlist(this%bndlist, ip)
1007  if (packobj2%packName == pakname) then
1008  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
1009  'already exists: ', trim(pakname)
1010  call store_error(errmsg, terminate=.true.)
1011  end if
1012  end do
1013  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
Here is the call graph for this function:

◆ prt_ad()

subroutine prtmodule::prt_ad ( class(prtmodeltype this)

Definition at line 364 of file prt.f90.

365  ! modules
366  use simvariablesmodule, only: isimcheck
367  ! dummy
368  class(PrtModelType) :: this
369  class(BndType), pointer :: packobj
370  ! local
371  integer(I4B) :: ip, n, i
372 
373  ! Discard buffered events from previous time step solve attempts.
374  ! prt_advance() is called on every sln_ca(): once per ATS retry,
375  ! once per Picard iteration, and once for the post-Picard output
376  ! rerun if mxiter > 1. Tracking is skipped during Picard iterations
377  ! (isuppress_output=1) in prt_solve, so the buffer is empty here on
378  ! Picard calls and the discard is a no-op. On the output rerun
379  ! (isuppress_output=0), tracking runs once and the buffer is cleared
380  ! before it is filled, so only events from that run reach disk.
381  call this%tracks%discard_buffer()
382 
383  ! Update look-behind mass
384  do n = 1, this%dis%nodes
385  this%massstoold(n) = this%masssto(n)
386  end do
387 
388  ! Advance fmi
389  call this%fmi%fmi_ad()
390 
391  ! Advance release packages
392  do ip = 1, this%bndlist%Count()
393  packobj => getbndfromlist(this%bndlist, ip)
394  call packobj%bnd_ad()
395  if (isimcheck > 0) &
396  call packobj%bnd_ck()
397  end do
398 
399  ! Initialize the flowja array. Flowja is calculated each time,
400  ! even if output is suppressed. (Flowja represents flow of particle
401  ! mass and is positive into a cell. Currently, each particle is assigned
402  ! unit mass.) Flowja is updated continually as particles are tracked
403  ! over the time step and at the end of the time step. The diagonal
404  ! position of the flowja array will contain the flow residual.
405  do i = 1, this%dis%nja
406  this%flowja(i) = dzero
407  end do
integer(i4b) isimcheck
simulation input check flag (1) to check input, (0) to ignore checks
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 245 of file prt.f90.

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

537  ! modules
538  use tdismodule, only: delt
539  use budgetmodule, only: rate_accumulator
540  ! dummy
541  class(PrtModelType) :: this
542  integer(I4B), intent(in) :: icnvg
543  integer(I4B), intent(in) :: isuppress_output
544  ! local
545  integer(I4B) :: ip
546  class(BndType), pointer :: packobj
547  real(DP) :: rin
548  real(DP) :: rout
549 
550  ! Budget routines (start by resetting). Sole purpose of this section
551  ! is to add in and outs to model budget. All ins and out for a model
552  ! should be added here to this%budget. In a subsequent exchange call,
553  ! exchange flows might also be added.
554  call this%budget%reset()
555  ! storage term
556  call rate_accumulator(this%ratesto, rin, rout)
557  call this%budget%addentry(rin, rout, delt, budtxt(1), &
558  isuppress_output, ' PRT')
559  ! termination term
560  call rate_accumulator(this%ratetrm, rin, rout)
561  call this%budget%addentry(rin, rout, delt, budtxt(2), &
562  isuppress_output, ' PRT')
563  ! boundary packages
564  do ip = 1, this%bndlist%Count()
565  packobj => getbndfromlist(this%bndlist, ip)
566  call packobj%bnd_bd(this%budget)
567  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:32
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 411 of file prt.f90.

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

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

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

847  ! modules
851  ! dummy
852  class(PrtModelType) :: this
853  ! local
854  integer(I4B) :: ip
855  class(BndType), pointer :: packobj
856 
857  ! Deallocate idm memory
858  call memorystore_remove(this%name, 'NAM', idm_context)
859  call memorystore_remove(component=this%name, context=idm_context)
860 
861  ! Internal packages
862  call this%dis%dis_da()
863  call this%fmi%fmi_da()
864  call this%mip%mip_da()
865  call this%budget%budget_da()
866  call this%oc%oc_da()
867  deallocate (this%dis)
868  deallocate (this%fmi)
869  deallocate (this%mip)
870  deallocate (this%budget)
871  deallocate (this%oc)
872 
873  ! Method objects
874  call this%method_dis%deallocate()
875  deallocate (this%method_dis)
876  call this%method_disv%deallocate()
877  deallocate (this%method_disv)
878 
879  ! Boundary packages
880  do ip = 1, this%bndlist%Count()
881  packobj => getbndfromlist(this%bndlist, ip)
882  call packobj%bnd_da()
883  deallocate (packobj)
884  end do
885 
886  ! Scalars
887  call mem_deallocate(this%infmi)
888  call mem_deallocate(this%inmip)
889  call mem_deallocate(this%inadv)
890  call mem_deallocate(this%indsp)
891  call mem_deallocate(this%inssm)
892  call mem_deallocate(this%inmst)
893  call mem_deallocate(this%inmvt)
894  call mem_deallocate(this%inoc)
895 
896  ! Arrays
897  call mem_deallocate(this%masssto)
898  call mem_deallocate(this%massstoold)
899  call mem_deallocate(this%ratesto)
900  call mem_deallocate(this%masstrm)
901  call mem_deallocate(this%ratetrm)
902 
903  call this%tracks%destroy()
904  deallocate (this%events)
905  deallocate (this%tracks)
906 
907  call this%ExplicitModelType%model_da()
subroutine, public memorystore_remove(component, subcomponent, context)
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 212 of file prt.f90.

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

◆ prt_ot()

subroutine prtmodule::prt_ot ( class(prtmodeltype this)

Definition at line 571 of file prt.f90.

572  use tdismodule, only: tdis_ot, endofperiod
573  use prtprpmodule, only: prtprptype
574  ! dummy
575  class(PrtModelType) :: this
576  ! local
577  integer(I4B) :: idvsave
578  integer(I4B) :: idvprint
579  integer(I4B) :: icbcfl
580  integer(I4B) :: icbcun
581  integer(I4B) :: ibudfl
582  integer(I4B) :: ipflag
583  integer(I4B) :: ip
584  class(BndType), pointer :: packobj
585 
586  ! Flush buffered events to disk
587  call this%tracks%flush_buffer()
588 
589  ! Commit each PRP's staged state
590  do ip = 1, this%bndlist%Count()
591  packobj => getbndfromlist(this%bndlist, ip)
592  select type (packobj)
593  type is (prtprptype)
594  call packobj%prp_commit()
595  end select
596  end do
597 
598  ! Set write and print flags
599  idvsave = 0
600  idvprint = 0
601  icbcfl = 0
602  ibudfl = 0
603  if (this%oc%oc_save('CONCENTRATION')) idvsave = 1
604  if (this%oc%oc_print('CONCENTRATION')) idvprint = 1
605  if (this%oc%oc_save('BUDGET')) icbcfl = 1
606  if (this%oc%oc_print('BUDGET')) ibudfl = 1
607  icbcun = this%oc%oc_save_unit('BUDGET')
608 
609  ! Override ibudfl and idvprint flags for nonconvergence
610  ! and end of period
611  ibudfl = this%oc%set_print_flag('BUDGET', 1, endofperiod)
612  idvprint = this%oc%set_print_flag('CONCENTRATION', 1, endofperiod)
613 
614  ! Save and print flows
615  call this%prt_ot_flow(icbcfl, ibudfl, icbcun)
616 
617  ! Save and print dependent variables
618  call this%prt_ot_dv(idvsave, idvprint, ipflag)
619 
620  ! Print budget summaries
621  call this%prt_ot_bdsummary(ibudfl, ipflag)
622 
623  ! Timing Output; if any dependent variables or budgets
624  ! are printed, then ipflag is set to 1.
625  if (ipflag == 1) call tdis_ot(this%iout)
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:30
subroutine, public tdis_ot(iout)
Print simulation time.
Definition: tdis.f90:271
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 816 of file prt.f90.

817  ! modules
818  use tdismodule, only: kstp, kper, totim, delt
819  ! dummy
820  class(PrtModelType) :: this
821  integer(I4B), intent(in) :: ibudfl
822  integer(I4B), intent(inout) :: ipflag
823  ! local
824  class(BndType), pointer :: packobj
825  integer(I4B) :: ip
826 
827  ! Package budget summary
828  do ip = 1, this%bndlist%Count()
829  packobj => getbndfromlist(this%bndlist, ip)
830  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
831  end do
832 
833  ! model budget summary
834  call this%budget%finalize_step(delt)
835  if (ibudfl /= 0) then
836  ipflag = 1
837  ! model budget summary
838  call this%budget%budget_ot(kstp, kper, this%iout)
839  end if
840 
841  ! Write to budget csv
842  call this%budget%writecsv(totim)
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:35
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:27
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:26
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 795 of file prt.f90.

796  ! dummy
797  class(PrtModelType) :: this
798  integer(I4B), intent(in) :: idvsave
799  integer(I4B), intent(in) :: idvprint
800  integer(I4B), intent(inout) :: ipflag
801  ! local
802  class(BndType), pointer :: packobj
803  integer(I4B) :: ip
804 
805  ! Print advanced package dependent variables
806  do ip = 1, this%bndlist%Count()
807  packobj => getbndfromlist(this%bndlist, ip)
808  call packobj%bnd_ot_dv(idvsave, idvprint)
809  end do
810 
811  ! save head and print head
812  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 629 of file prt.f90.

630  use prtprpmodule, only: prtprptype
631  class(PrtModelType) :: this
632  integer(I4B), intent(in) :: icbcfl
633  integer(I4B), intent(in) :: ibudfl
634  integer(I4B), intent(in) :: icbcun
635  class(BndType), pointer :: packobj
636  integer(I4B) :: ip
637 
638  ! Save PRT flows
639  call this%prt_ot_saveflow(this%dis%nja, this%flowja, icbcfl, icbcun)
640  do ip = 1, this%bndlist%Count()
641  packobj => getbndfromlist(this%bndlist, ip)
642  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
643  end do
644 
645  ! Save advanced package flows
646  do ip = 1, this%bndlist%Count()
647  packobj => getbndfromlist(this%bndlist, ip)
648  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
649  end do
650 
651  ! Print PRT flows
652  call this%prt_ot_printflow(ibudfl, this%flowja)
653  do ip = 1, this%bndlist%Count()
654  packobj => getbndfromlist(this%bndlist, ip)
655  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
656  end do
657 
658  ! Print advanced package flows
659  do ip = 1, this%bndlist%Count()
660  packobj => getbndfromlist(this%bndlist, ip)
661  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
662  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 756 of file prt.f90.

757  ! modules
758  use tdismodule, only: kper, kstp
759  use constantsmodule, only: lenbigline
760  ! dummy
761  class(PrtModelType) :: this
762  integer(I4B), intent(in) :: ibudfl
763  real(DP), intent(inout), dimension(:) :: flowja
764  ! local
765  character(len=LENBIGLINE) :: line
766  character(len=30) :: tempstr
767  integer(I4B) :: n, ipos, m
768  real(DP) :: qnm
769  ! formats
770  character(len=*), parameter :: fmtiprflow = &
771  "(/,4x,'CALCULATED INTERCELL FLOW &
772  &FOR PERIOD ', i0, ' STEP ', i0)"
773 
774  ! Write flowja to list file if requested
775  if (ibudfl /= 0 .and. this%iprflow > 0) then
776  write (this%iout, fmtiprflow) kper, kstp
777  do n = 1, this%dis%nodes
778  line = ''
779  call this%dis%noder_to_string(n, tempstr)
780  line = trim(tempstr)//':'
781  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
782  m = this%dis%con%ja(ipos)
783  call this%dis%noder_to_string(m, tempstr)
784  line = trim(line)//' '//trim(tempstr)
785  qnm = flowja(ipos)
786  write (tempstr, '(1pg15.6)') qnm
787  line = trim(line)//' '//trim(adjustl(tempstr))
788  end do
789  write (this%iout, '(a)') trim(line)
790  end do
791  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 666 of file prt.f90.

667  ! dummy
668  class(PrtModelType) :: this
669  integer(I4B), intent(in) :: nja
670  real(DP), dimension(nja), intent(in) :: flowja
671  integer(I4B), intent(in) :: icbcfl
672  integer(I4B), intent(in) :: icbcun
673  ! local
674  integer(I4B) :: ibinun
675  integer(I4B) :: naux
676  real(DP), dimension(0) :: auxrow
677  character(len=LENAUXNAME), dimension(0) :: auxname
678  logical(LGP) :: header_written
679  integer(I4B) :: i, nn
680  real(DP) :: m
681  integer(I4B) :: nsto, ntrm
682  logical(LGP), allocatable :: msto_mask(:), mtrm_mask(:)
683  integer(I4B), allocatable :: msto_nns(:), mtrm_nns(:)
684  real(DP), allocatable :: msto_vals(:), mtrm_vals(:)
685 
686  ! Set unit number for binary output
687  if (this%ipakcb < 0) then
688  ibinun = icbcun
689  elseif (this%ipakcb == 0) then
690  ibinun = 0
691  else
692  ibinun = this%ipakcb
693  end if
694  if (icbcfl == 0) ibinun = 0
695 
696  ! Return if nothing to do
697  if (ibinun == 0) return
698 
699  ! Write mass face flows
700  call this%dis%record_connection_array(flowja, ibinun, this%iout)
701 
702  ! Write mass storage term
703  naux = 0
704  header_written = .false.
705  msto_mask = this%masssto > dzero
706  msto_vals = pack(this%masssto, msto_mask)
707  msto_nns = [(i, i=1, size(this%masssto))]
708  msto_nns = pack(msto_nns, msto_mask)
709  nsto = size(msto_nns)
710  do i = 1, nsto
711  nn = msto_nns(i)
712  m = msto_vals(i)
713  if (.not. header_written) then
714  call this%dis%record_srcdst_list_header(budtxt(1), &
715  'PRT ', &
716  'PRT ', &
717  'PRT ', &
718  'STORAGE ', &
719  naux, auxname, ibinun, &
720  nsto, this%iout)
721  header_written = .true.
722  end if
723  call this%dis%record_mf6_list_entry(ibinun, nn, nn, m, &
724  0, auxrow, &
725  olconv2=.false.)
726  end do
727 
728  ! Write mass termination term
729  header_written = .false.
730  mtrm_mask = this%masstrm > dzero
731  mtrm_vals = pack(this%masstrm, mtrm_mask)
732  mtrm_nns = [(i, i=1, size(this%masstrm))]
733  mtrm_nns = pack(mtrm_nns, mtrm_mask)
734  ntrm = size(mtrm_nns)
735  do i = 1, ntrm
736  nn = mtrm_nns(i)
737  m = mtrm_vals(i)
738  if (.not. header_written) then
739  call this%dis%record_srcdst_list_header(budtxt(2), &
740  'PRT ', &
741  'PRT ', &
742  'PRT ', &
743  'TERMINATION ', &
744  naux, auxname, ibinun, &
745  ntrm, this%iout)
746  header_written = .true.
747  end if
748  call this%dis%record_mf6_list_entry(ibinun, nn, nn, m, &
749  0, auxrow, &
750  olconv2=.false.)
751  end do
752 

◆ prt_rp()

subroutine prtmodule::prt_rp ( class(prtmodeltype this)

Definition at line 344 of file prt.f90.

345  use tdismodule, only: readnewdata
346  ! dummy
347  class(PrtModelType) :: this
348  ! local
349  class(BndType), pointer :: packobj
350  integer(I4B) :: ip
351 
352  ! Check with TDIS on whether or not it is time to RP
353  if (.not. readnewdata) return
354 
355  ! Read and prepare
356  if (this%inoc > 0) call this%oc%oc_rp()
357  do ip = 1, this%bndlist%Count()
358  packobj => getbndfromlist(this%bndlist, ip)
359  call packobj%bnd_rp()
360  end do
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:29
Here is the call graph for this function:

◆ prt_solve()

subroutine prtmodule::prt_solve ( class(prtmodeltype this,
integer(i4b), intent(in)  isuppress_output 
)
private

Definition at line 1044 of file prt.f90.

1045  use tdismodule, only: totimc, delt, endofsimulation
1046  use prtprpmodule, only: prtprptype
1049  ! dummy
1050  class(PrtModelType) :: this
1051  integer(I4B), intent(in) :: isuppress_output
1052  ! local
1053  integer(I4B) :: np, ip
1054  class(BndType), pointer :: packobj
1055  type(ParticleType), pointer :: particle
1056  real(DP) :: tmax
1057  integer(I4B) :: iprp
1058 
1059  ! Skip tracking during Picard iterations: PRT doesn't affect the flow
1060  ! solution, so tracking is only needed once on the final output rerun
1061  ! (or on the single call when mxiter=1, which also has isuppress_output=0).
1062  if (isuppress_output /= 0) return
1063 
1064  ! A single particle is reused in the tracking loops
1065  ! to avoid allocating and deallocating it each time.
1066  ! get() and put() retrieve and store particle state.
1067  call create_particle(particle)
1068  ! Loop over PRP packages and particles within them.
1069  iprp = 0
1070  do ip = 1, this%bndlist%Count()
1071  packobj => getbndfromlist(this%bndlist, ip)
1072  select type (packobj)
1073  type is (prtprptype)
1074  iprp = iprp + 1
1075  do np = 1, packobj%nparticles
1076  ! Get the particle from the staging store
1077  call packobj%particles_staging%get(particle, this%id, iprp, np)
1078  ! If particle is permanently unreleased, cycle.
1079  ! Raise a termination event if we haven't yet.
1080  ! TODO: when we have generic dynamic vectors,
1081  ! consider terminating permanently unreleased
1082  ! in PRP instead of here. For now, status -8
1083  ! indicates the permanently unreleased event
1084  ! is not yet recorded, status 8 it has been.
1085  if (particle%istatus == (-1 * term_unreleased)) then
1086  call this%method%terminate(particle, status=term_unreleased)
1087  call packobj%particles_staging%put(particle, np)
1088  end if
1089  if (particle%istatus > active) cycle ! Skip terminated particles
1090  particle%istatus = active ! Set active status in case of release
1091  ! If the particle was released this time step, emit a release event
1092  if (particle%trelease >= totimc) call this%method%release(particle)
1093  ! Maximum time is the end of the time step or the particle
1094  ! stop time, whichever comes first, unless it's the final
1095  ! time step and the extend option is on, in which case
1096  ! it's just the particle stop time.
1097  if (endofsimulation .and. particle%extend) then
1098  tmax = particle%tstop
1099  else
1100  tmax = min(totimc + delt, particle%tstop)
1101  end if
1102  ! Apply the tracking method until the maximum time.
1103  call this%method%apply(particle, tmax)
1104  ! If the particle timed out, terminate it.
1105  ! "Timed out" means it's still active but
1106  ! - it reached its stop time, or
1107  ! - the simulation is over.
1108  ! We can't detect timeout within the tracking
1109  ! method because the method just receives the
1110  ! maximum time with no context on what it is.
1111  ! TODO maybe think about changing that?
1112  if (particle%istatus <= active .and. &
1113  (particle%ttrack == particle%tstop .or. endofsimulation)) &
1114  call this%method%terminate(particle, status=term_timeout)
1115  ! Return the particle to the staging store
1116  call packobj%particles_staging%put(particle, np)
1117  end do
1118  end select
1119  end do
1120  call particle%destroy()
1121  deallocate (particle)
@, public release
particle was released
@, public terminate
particle terminated
@ term_timeout
terminated at stop time or end of simulation
Definition: Particle.f90:40
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:38
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:31
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:36
Here is the call graph for this function:

Variable Documentation

◆ budtxt

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

Definition at line 43 of file prt.f90.

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

◆ nbditems

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

Definition at line 42 of file prt.f90.

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

◆ niunit_prt

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

Definition at line 127 of file prt.f90.

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

◆ prt_basepkg

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

Definition at line 108 of file prt.f90.

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

◆ prt_multipkg

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

Definition at line 122 of file prt.f90.

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

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

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