MODFLOW 6  version 6.7.0.dev1
USGS Modular Hydrologic Model
prt.f90
Go to the documentation of this file.
1 module prtmodule
2  use kindmodule, only: dp, i4b, lgp
3  use errorutilmodule, only: pstop
12  use dismodule, only: distype, dis_cr
13  use disvmodule, only: disvtype, disv_cr
14  use disumodule, only: disutype, disu_cr
15  use prtprpmodule, only: prtprptype
16  use prtfmimodule, only: prtfmitype
17  use prtmipmodule, only: prtmiptype
18  use prtocmodule, only: prtoctype
19  use budgetmodule, only: budgettype
20  use listmodule, only: listtype
26  use methodmodule, only: methodtype
27 
28  implicit none
29 
30  private
31  public :: prt_cr
32  public :: prtmodeltype
33  public :: prt_nbasepkg, prt_nmultipkg
34  public :: prt_basepkg, prt_multipkg
35 
36  integer(I4B), parameter :: nbditems = 1
37  character(len=LENBUDTXT), dimension(NBDITEMS) :: budtxt
38  data budtxt/' STORAGE'/
39 
40  !> @brief Particle tracking (PRT) model
41  type, extends(numericalmodeltype) :: prtmodeltype
42  type(prtfmitype), pointer :: fmi => null() ! flow model interface
43  type(prtmiptype), pointer :: mip => null() ! model input package
44  type(prtoctype), pointer :: oc => null() ! output control package
45  type(budgettype), pointer :: budget => null() ! budget object
46  class(methodtype), pointer :: method => null() ! tracking method
47  type(trackcontroltype), pointer :: trackctl ! track control
48  integer(I4B), pointer :: infmi => null() ! unit number FMI
49  integer(I4B), pointer :: inmip => null() ! unit number MIP
50  integer(I4B), pointer :: inmvt => null() ! unit number MVT
51  integer(I4B), pointer :: inmst => null() ! unit number MST
52  integer(I4B), pointer :: inadv => null() ! unit number ADV
53  integer(I4B), pointer :: indsp => null() ! unit number DSP
54  integer(I4B), pointer :: inssm => null() ! unit number SSM
55  integer(I4B), pointer :: inoc => null() ! unit number OC
56  integer(I4B), pointer :: nprp => null() ! number of PRP packages in the model
57  real(dp), dimension(:), pointer, contiguous :: masssto => null() !< particle mass storage in cells, new value
58  real(dp), dimension(:), pointer, contiguous :: massstoold => null() !< particle mass storage in cells, old value
59  real(dp), dimension(:), pointer, contiguous :: ratesto => null() !< particle mass storage rate in cells
60  contains
61  ! Override BaseModelType procs
62  procedure :: model_df => prt_df
63  procedure :: model_ar => prt_ar
64  procedure :: model_rp => prt_rp
65  procedure :: model_ad => prt_ad
66  procedure :: model_cq => prt_cq
67  procedure :: model_bd => prt_bd
68  procedure :: model_ot => prt_ot
69  procedure :: model_da => prt_da
70  procedure :: model_solve => prt_solve
71 
72  ! Private utilities
73  procedure :: allocate_scalars
74  procedure :: allocate_arrays
75  procedure, private :: package_create
76  procedure, private :: ftype_check
77  procedure, private :: prt_ot_flow
78  procedure, private :: prt_ot_saveflow
79  procedure, private :: prt_ot_printflow
80  procedure, private :: prt_ot_dv
81  procedure, private :: prt_ot_bdsummary
82  procedure, private :: prt_cq_sto
83  procedure, private :: create_packages
84  procedure, private :: create_bndpkgs
85  procedure, private :: log_namfile_options
86 
87  end type prtmodeltype
88 
89  !> @brief PRT base package array descriptors
90  !!
91  !! PRT6 model base package types. Only listed packages are candidates
92  !! for input and these will be loaded in the order specified.
93  !<
94  integer(I4B), parameter :: prt_nbasepkg = 50
95  character(len=LENPACKAGETYPE), dimension(PRT_NBASEPKG) :: prt_basepkg
96  data prt_basepkg/'DIS6 ', 'DISV6', 'DISU6', 'IC6 ', 'MST6 ', & ! 5
97  &'ADV6 ', 'DSP6 ', 'SSM6 ', 'MIP6 ', 'CNC6 ', & ! 10
98  &'OC6 ', ' ', 'FMI6 ', ' ', 'IST6 ', & ! 15
99  &'LKT6 ', 'SFT6 ', 'MWT6 ', 'UZT6 ', 'MVT6 ', & ! 20
100  &'API6 ', ' ', ' ', ' ', ' ', & ! 25
101  25*' '/ ! 50
102 
103  !> @brief PRT multi package array descriptors
104  !!
105  !! PRT6 model multi-instance package types. Only listed packages are
106  !! candidates for input and these will be loaded in the order specified.
107  !<
108  integer(I4B), parameter :: prt_nmultipkg = 50
109  character(len=LENPACKAGETYPE), dimension(PRT_NMULTIPKG) :: prt_multipkg
110  data prt_multipkg/'PRP6 ', ' ', ' ', ' ', ' ', & ! 5
111  &45*' '/ ! 50
112 
113  ! size of supported model package arrays
114  integer(I4B), parameter :: niunit_prt = prt_nbasepkg + prt_nmultipkg
115 
116 contains
117 
118  !> @brief Create a new particle tracking model object
119  subroutine prt_cr(filename, id, modelname)
120  ! modules
121  use listsmodule, only: basemodellist
124  use compilerversion
129  ! dummy
130  character(len=*), intent(in) :: filename
131  integer(I4B), intent(in) :: id
132  character(len=*), intent(in) :: modelname
133  ! local
134  type(prtmodeltype), pointer :: this
135  class(basemodeltype), pointer :: model
136  character(len=LENMEMPATH) :: input_mempath
137  character(len=LINELENGTH) :: lst_fname
138  type(gwfnamparamfoundtype) :: found
139 
140  ! Allocate a new PRT Model (this)
141  allocate (this)
142 
143  ! Set this before any allocs in the memory manager can be done
144  this%memoryPath = create_mem_path(modelname)
145 
146  ! Allocate track control object
147  allocate (this%trackctl)
148 
149  ! Allocate scalars and add model to basemodellist
150  call this%allocate_scalars(modelname)
151  model => this
152  call addbasemodeltolist(basemodellist, model)
153 
154  ! Assign variables
155  this%filename = filename
156  this%name = modelname
157  this%macronym = 'PRT'
158  this%id = id
159 
160  ! Set input model namfile memory path
161  input_mempath = create_mem_path(modelname, 'NAM', idm_context)
162 
163  ! Copy options from input context
164  call mem_set_value(this%iprpak, 'PRINT_INPUT', input_mempath, &
165  found%print_input)
166  call mem_set_value(this%iprflow, 'PRINT_FLOWS', input_mempath, &
167  found%print_flows)
168  call mem_set_value(this%ipakcb, 'SAVE_FLOWS', input_mempath, &
169  found%save_flows)
170 
171  ! Create the list file
172  call this%create_lstfile(lst_fname, filename, found%list, &
173  'PARTICLE TRACKING MODEL (PRT)')
174 
175  ! Activate save_flows if found
176  if (found%save_flows) then
177  this%ipakcb = -1
178  end if
179 
180  ! Log options
181  if (this%iout > 0) then
182  call this%log_namfile_options(found)
183  end if
184 
185  ! Create model packages
186  call this%create_packages()
187  end subroutine prt_cr
188 
189  !> @brief Define packages
190  !!
191  !! (1) call df routines for each package
192  !! (2) set variables and pointers
193  !<
194  subroutine prt_df(this)
195  ! modules
196  use prtprpmodule, only: prtprptype
197  ! dummy
198  class(prtmodeltype) :: this
199  ! local
200  integer(I4B) :: ip
201  class(bndtype), pointer :: packobj
202 
203  ! Define packages and utility objects
204  call this%dis%dis_df()
205  call this%fmi%fmi_df(this%dis, 1)
206  call this%oc%oc_df()
207  call this%budget%budget_df(niunit_prt, 'MASS', 'M')
208 
209  ! Define packages and assign iout for time series managers
210  do ip = 1, this%bndlist%Count()
211  packobj => getbndfromlist(this%bndlist, ip)
212  call packobj%bnd_df(this%dis%nodes, this%dis)
213  packobj%TsManager%iout = this%iout
214  packobj%TasManager%iout = this%iout
215  end do
216 
217  ! Allocate model arrays
218  call this%allocate_arrays()
219 
220  end subroutine prt_df
221 
222  !> @brief Allocate and read
223  !!
224  !! (1) allocates and reads packages part of this model,
225  !! (2) allocates memory for arrays part of this model object
226  !<
227  subroutine prt_ar(this)
228  ! modules
229  use constantsmodule, only: dhnoflo
230  use prtprpmodule, only: prtprptype
231  use prtmipmodule, only: prtmiptype
233  ! dummy
234  class(prtmodeltype) :: this
235  ! locals
236  integer(I4B) :: ip
237  class(bndtype), pointer :: packobj
238 
239  ! Allocate and read modules attached to model
240  call this%fmi%fmi_ar(this%ibound)
241  if (this%inmip > 0) call this%mip%mip_ar()
242 
243  ! set up output control
244  call this%oc%oc_ar(this%dis, dhnoflo)
245  call this%budget%set_ibudcsv(this%oc%ibudcsv)
246 
247  ! Package input files now open, so allocate and read
248  do ip = 1, this%bndlist%Count()
249  packobj => getbndfromlist(this%bndlist, ip)
250  select type (packobj)
251  type is (prtprptype)
252  call packobj%prp_set_pointers(this%ibound, this%mip%izone, &
253  this%trackctl)
254  end select
255  ! Read and allocate package
256  call packobj%bnd_ar()
257  end do
258 
259  ! Initialize tracking method
260  select type (dis => this%dis)
261  type is (distype)
262  call method_dis%init( &
263  fmi=this%fmi, &
264  trackctl=this%trackctl, &
265  izone=this%mip%izone, &
266  flowja=this%flowja, &
267  porosity=this%mip%porosity, &
268  retfactor=this%mip%retfactor, &
269  tracktimes=this%oc%tracktimes)
270  this%method => method_dis
271  type is (disvtype)
272  call method_disv%init( &
273  fmi=this%fmi, &
274  trackctl=this%trackctl, &
275  izone=this%mip%izone, &
276  flowja=this%flowja, &
277  porosity=this%mip%porosity, &
278  retfactor=this%mip%retfactor, &
279  tracktimes=this%oc%tracktimes)
280  this%method => method_disv
281  end select
282 
283  ! Initialize track output files and reporting options
284  if (this%oc%itrkout > 0) &
285  call this%trackctl%init_track_file(this%oc%itrkout)
286  if (this%oc%itrkcsv > 0) &
287  call this%trackctl%init_track_file(this%oc%itrkcsv, csv=.true.)
288  call this%trackctl%set_track_events( &
289  this%oc%trackrelease, &
290  this%oc%trackexit, &
291  this%oc%tracktimestep, &
292  this%oc%trackterminate, &
293  this%oc%trackweaksink, &
294  this%oc%trackusertime)
295  end subroutine prt_ar
296 
297  !> @brief Read and prepare (calls package read and prepare routines)
298  subroutine prt_rp(this)
299  use tdismodule, only: readnewdata
300  ! dummy
301  class(prtmodeltype) :: this
302  ! local
303  class(bndtype), pointer :: packobj
304  integer(I4B) :: ip
305 
306  ! Check with TDIS on whether or not it is time to RP
307  if (.not. readnewdata) return
308 
309  ! Read and prepare
310  if (this%inoc > 0) call this%oc%oc_rp()
311  do ip = 1, this%bndlist%Count()
312  packobj => getbndfromlist(this%bndlist, ip)
313  call packobj%bnd_rp()
314  end do
315  end subroutine prt_rp
316 
317  !> @brief Time step advance (calls package advance subroutines)
318  subroutine prt_ad(this)
319  ! modules
321  ! dummy
322  class(prtmodeltype) :: this
323  class(bndtype), pointer :: packobj
324  ! local
325  integer(I4B) :: irestore
326  integer(I4B) :: ip, n, i
327 
328  ! Reset state variable
329  irestore = 0
330  if (ifailedstepretry > 0) irestore = 1
331 
332  ! Copy masssto into massstoold
333  do n = 1, this%dis%nodes
334  this%massstoold(n) = this%masssto(n)
335  end do
336 
337  ! Advance fmi
338  call this%fmi%fmi_ad()
339 
340  ! Advance
341  do ip = 1, this%bndlist%Count()
342  packobj => getbndfromlist(this%bndlist, ip)
343  call packobj%bnd_ad()
344  if (isimcheck > 0) then
345  call packobj%bnd_ck()
346  end if
347  end do
348  !
349  ! Initialize the flowja array. Flowja is calculated each time,
350  ! even if output is suppressed. (Flowja represents flow of particle
351  ! mass and is positive into a cell. Currently, each particle is assigned
352  ! unit mass.) Flowja is updated continually as particles are tracked
353  ! over the time step and at the end of the time step. The diagonal
354  ! position of the flowja array will contain the flow residual.
355  do i = 1, this%dis%nja
356  this%flowja(i) = dzero
357  end do
358  end subroutine prt_ad
359 
360  !> @brief Calculate intercell flow (flowja)
361  subroutine prt_cq(this, icnvg, isuppress_output)
362  ! modules
363  use sparsemodule, only: csr_diagsum
364  use tdismodule, only: delt
365  use prtprpmodule, only: prtprptype
366  ! dummy
367  class(prtmodeltype) :: this
368  integer(I4B), intent(in) :: icnvg
369  integer(I4B), intent(in) :: isuppress_output
370  ! local
371  integer(I4B) :: i
372  integer(I4B) :: ip
373  class(bndtype), pointer :: packobj
374  real(DP) :: tled
375 
376  ! Flowja is calculated each time, even if output is suppressed.
377  ! Flowja represents flow of particle mass and is positive into a cell.
378  ! Currently, each particle is assigned unit mass.
379  !
380  ! Reciprocal of time step size.
381  tled = done / delt
382  !
383  ! Flowja was updated continually as particles were tracked over the
384  ! time step. At this point, flowja contains the net particle mass
385  ! exchanged between cells during the time step. To convert these to
386  ! flow rates (particle mass per time), divide by the time step size.
387  do i = 1, this%dis%nja
388  this%flowja(i) = this%flowja(i) * tled
389  end do
390 
391  ! Particle mass storage
392  call this%prt_cq_sto()
393 
394  ! Go through packages and call cq routines. Just a formality.
395  do ip = 1, this%bndlist%Count()
396  packobj => getbndfromlist(this%bndlist, ip)
397  call packobj%bnd_cq(this%masssto, this%flowja)
398  end do
399 
400  ! Finalize calculation of flowja by adding face flows to the diagonal.
401  ! This results in the flow residual being stored in the diagonal
402  ! position for each cell.
403  call csr_diagsum(this%dis%con%ia, this%flowja)
404  end subroutine prt_cq
405 
406  !> @brief Calculate particle mass storage
407  subroutine prt_cq_sto(this)
408  ! modules
409  use tdismodule, only: delt
410  use prtprpmodule, only: prtprptype
411  ! dummy
412  class(prtmodeltype) :: this
413  ! local
414  integer(I4B) :: ip
415  class(bndtype), pointer :: packobj
416  integer(I4B) :: n
417  integer(I4B) :: np
418  integer(I4B) :: idiag
419  integer(I4B) :: istatus
420  real(DP) :: tled
421  real(DP) :: rate
422 
423  ! Reciprocal of time step size.
424  tled = done / delt
425 
426  ! Particle mass storage rate
427  do n = 1, this%dis%nodes
428  this%masssto(n) = dzero
429  this%ratesto(n) = dzero
430  end do
431  do ip = 1, this%bndlist%Count()
432  packobj => getbndfromlist(this%bndlist, ip)
433  select type (packobj)
434  type is (prtprptype)
435  do np = 1, packobj%nparticles
436  istatus = packobj%particles%istatus(np)
437  ! this may need to change if istatus flags change
438  if ((istatus > 0) .and. (istatus /= 8)) then
439  n = packobj%particles%idomain(np, 2)
440  ! Each particle currently assigned unit mass
441  this%masssto(n) = this%masssto(n) + done
442  end if
443  end do
444  end select
445  end do
446  do n = 1, this%dis%nodes
447  rate = -(this%masssto(n) - this%massstoold(n)) * tled
448  this%ratesto(n) = rate
449  idiag = this%dis%con%ia(n)
450  this%flowja(idiag) = this%flowja(idiag) + rate
451  end do
452  end subroutine prt_cq_sto
453 
454  !> @brief Calculate flows and budget
455  !!
456  !! (1) Calculate intercell flows (flowja)
457  !! (2) Calculate package contributions to model budget
458  !!
459  !<
460  subroutine prt_bd(this, icnvg, isuppress_output)
461  ! modules
462  use tdismodule, only: delt
463  use budgetmodule, only: rate_accumulator
464  ! dummy
465  class(prtmodeltype) :: this
466  integer(I4B), intent(in) :: icnvg
467  integer(I4B), intent(in) :: isuppress_output
468  ! local
469  integer(I4B) :: ip
470  class(bndtype), pointer :: packobj
471  real(DP) :: rin
472  real(DP) :: rout
473 
474  ! Budget routines (start by resetting). Sole purpose of this section
475  ! is to add in and outs to model budget. All ins and out for a model
476  ! should be added here to this%budget. In a subsequent exchange call,
477  ! exchange flows might also be added.
478  call this%budget%reset()
479  call rate_accumulator(this%ratesto, rin, rout)
480  call this%budget%addentry(rin, rout, delt, budtxt(1), &
481  isuppress_output, ' PRT')
482  do ip = 1, this%bndlist%Count()
483  packobj => getbndfromlist(this%bndlist, ip)
484  call packobj%bnd_bd(this%budget)
485  end do
486  end subroutine prt_bd
487 
488  !> @brief Print and/or save model output
489  subroutine prt_ot(this)
490  use tdismodule, only: tdis_ot, endofperiod
491  ! dummy
492  class(prtmodeltype) :: this
493  ! local
494  integer(I4B) :: idvsave
495  integer(I4B) :: idvprint
496  integer(I4B) :: icbcfl
497  integer(I4B) :: icbcun
498  integer(I4B) :: ibudfl
499  integer(I4B) :: ipflag
500 
501  ! Note: particle tracking output is handled elsewhere
502 
503  ! Set write and print flags
504  idvsave = 0
505  idvprint = 0
506  icbcfl = 0
507  ibudfl = 0
508  if (this%oc%oc_save('CONCENTRATION')) idvsave = 1
509  if (this%oc%oc_print('CONCENTRATION')) idvprint = 1
510  if (this%oc%oc_save('BUDGET')) icbcfl = 1
511  if (this%oc%oc_print('BUDGET')) ibudfl = 1
512  icbcun = this%oc%oc_save_unit('BUDGET')
513 
514  ! Override ibudfl and idvprint flags for nonconvergence
515  ! and end of period
516  ibudfl = this%oc%set_print_flag('BUDGET', 1, endofperiod)
517  idvprint = this%oc%set_print_flag('CONCENTRATION', 1, endofperiod)
518 
519  ! Save and print flows
520  call this%prt_ot_flow(icbcfl, ibudfl, icbcun)
521 
522  ! Save and print dependent variables
523  call this%prt_ot_dv(idvsave, idvprint, ipflag)
524 
525  ! Print budget summaries
526  call this%prt_ot_bdsummary(ibudfl, ipflag)
527 
528  ! Timing Output; if any dependent variables or budgets
529  ! are printed, then ipflag is set to 1.
530  if (ipflag == 1) call tdis_ot(this%iout)
531  end subroutine prt_ot
532 
533  !> @brief Save flows
534  subroutine prt_ot_flow(this, icbcfl, ibudfl, icbcun)
535  use prtprpmodule, only: prtprptype
536  class(prtmodeltype) :: this
537  integer(I4B), intent(in) :: icbcfl
538  integer(I4B), intent(in) :: ibudfl
539  integer(I4B), intent(in) :: icbcun
540  class(bndtype), pointer :: packobj
541  integer(I4B) :: ip
542 
543  ! Save PRT flows
544  call this%prt_ot_saveflow(this%dis%nja, this%flowja, icbcfl, icbcun)
545  do ip = 1, this%bndlist%Count()
546  packobj => getbndfromlist(this%bndlist, ip)
547  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=0, icbcun=icbcun)
548  end do
549 
550  ! Save advanced package flows
551  do ip = 1, this%bndlist%Count()
552  packobj => getbndfromlist(this%bndlist, ip)
553  call packobj%bnd_ot_package_flows(icbcfl=icbcfl, ibudfl=0)
554  end do
555 
556  ! Print PRT flows
557  call this%prt_ot_printflow(ibudfl, this%flowja)
558  do ip = 1, this%bndlist%Count()
559  packobj => getbndfromlist(this%bndlist, ip)
560  call packobj%bnd_ot_model_flows(icbcfl=icbcfl, ibudfl=ibudfl, icbcun=0)
561  end do
562 
563  ! Print advanced package flows
564  do ip = 1, this%bndlist%Count()
565  packobj => getbndfromlist(this%bndlist, ip)
566  call packobj%bnd_ot_package_flows(icbcfl=0, ibudfl=ibudfl)
567  end do
568  end subroutine prt_ot_flow
569 
570  !> @brief Save intercell flows
571  subroutine prt_ot_saveflow(this, nja, flowja, icbcfl, icbcun)
572  ! dummy
573  class(prtmodeltype) :: this
574  integer(I4B), intent(in) :: nja
575  real(DP), dimension(nja), intent(in) :: flowja
576  integer(I4B), intent(in) :: icbcfl
577  integer(I4B), intent(in) :: icbcun
578  ! local
579  integer(I4B) :: ibinun
580 
581  ! Set unit number for binary output
582  if (this%ipakcb < 0) then
583  ibinun = icbcun
584  elseif (this%ipakcb == 0) then
585  ibinun = 0
586  else
587  ibinun = this%ipakcb
588  end if
589  if (icbcfl == 0) ibinun = 0
590 
591  ! Write the face flows if requested
592  if (ibinun /= 0) then
593  call this%dis%record_connection_array(flowja, ibinun, this%iout)
594  end if
595  end subroutine prt_ot_saveflow
596 
597  !> @brief Print intercell flows
598  subroutine prt_ot_printflow(this, ibudfl, flowja)
599  ! modules
600  use tdismodule, only: kper, kstp
601  use constantsmodule, only: lenbigline
602  ! dummy
603  class(prtmodeltype) :: this
604  integer(I4B), intent(in) :: ibudfl
605  real(DP), intent(inout), dimension(:) :: flowja
606  ! local
607  character(len=LENBIGLINE) :: line
608  character(len=30) :: tempstr
609  integer(I4B) :: n, ipos, m
610  real(DP) :: qnm
611  ! formats
612  character(len=*), parameter :: fmtiprflow = &
613  "(/,4x,'CALCULATED INTERCELL FLOW &
614  &FOR PERIOD ', i0, ' STEP ', i0)"
615 
616  ! Write flowja to list file if requested
617  if (ibudfl /= 0 .and. this%iprflow > 0) then
618  write (this%iout, fmtiprflow) kper, kstp
619  do n = 1, this%dis%nodes
620  line = ''
621  call this%dis%noder_to_string(n, tempstr)
622  line = trim(tempstr)//':'
623  do ipos = this%dis%con%ia(n) + 1, this%dis%con%ia(n + 1) - 1
624  m = this%dis%con%ja(ipos)
625  call this%dis%noder_to_string(m, tempstr)
626  line = trim(line)//' '//trim(tempstr)
627  qnm = flowja(ipos)
628  write (tempstr, '(1pg15.6)') qnm
629  line = trim(line)//' '//trim(adjustl(tempstr))
630  end do
631  write (this%iout, '(a)') trim(line)
632  end do
633  end if
634  end subroutine prt_ot_printflow
635 
636  !> @brief Print dependent variables
637  subroutine prt_ot_dv(this, idvsave, idvprint, ipflag)
638  ! dummy
639  class(prtmodeltype) :: this
640  integer(I4B), intent(in) :: idvsave
641  integer(I4B), intent(in) :: idvprint
642  integer(I4B), intent(inout) :: ipflag
643  ! local
644  class(bndtype), pointer :: packobj
645  integer(I4B) :: ip
646 
647  ! Print advanced package dependent variables
648  do ip = 1, this%bndlist%Count()
649  packobj => getbndfromlist(this%bndlist, ip)
650  call packobj%bnd_ot_dv(idvsave, idvprint)
651  end do
652 
653  ! save head and print head
654  call this%oc%oc_ot(ipflag)
655  end subroutine prt_ot_dv
656 
657  !> @brief Print budget summary
658  subroutine prt_ot_bdsummary(this, ibudfl, ipflag)
659  ! modules
660  use tdismodule, only: kstp, kper, totim, delt
661  ! dummy
662  class(prtmodeltype) :: this
663  integer(I4B), intent(in) :: ibudfl
664  integer(I4B), intent(inout) :: ipflag
665  ! local
666  class(bndtype), pointer :: packobj
667  integer(I4B) :: ip
668 
669  ! Package budget summary
670  do ip = 1, this%bndlist%Count()
671  packobj => getbndfromlist(this%bndlist, ip)
672  call packobj%bnd_ot_bdsummary(kstp, kper, this%iout, ibudfl)
673  end do
674 
675  ! model budget summary
676  call this%budget%finalize_step(delt)
677  if (ibudfl /= 0) then
678  ipflag = 1
679  ! model budget summary
680  call this%budget%budget_ot(kstp, kper, this%iout)
681  end if
682 
683  ! Write to budget csv
684  call this%budget%writecsv(totim)
685  end subroutine prt_ot_bdsummary
686 
687  !> @brief Deallocate
688  subroutine prt_da(this)
689  ! modules
696  ! dummy
697  class(prtmodeltype) :: this
698  ! local
699  integer(I4B) :: ip
700  class(bndtype), pointer :: packobj
701 
702  ! Deallocate idm memory
703  call memorystore_remove(this%name, 'NAM', idm_context)
704  call memorystore_remove(component=this%name, context=idm_context)
705 
706  ! Internal packages
707  call this%dis%dis_da()
708  call this%fmi%fmi_da()
709  call this%mip%mip_da()
710  call this%budget%budget_da()
711  call this%oc%oc_da()
712  deallocate (this%dis)
713  deallocate (this%fmi)
714  deallocate (this%mip)
715  deallocate (this%budget)
716  deallocate (this%oc)
717 
718  ! Method objects
721  call destroy_method_pool()
722 
723  ! Boundary packages
724  do ip = 1, this%bndlist%Count()
725  packobj => getbndfromlist(this%bndlist, ip)
726  call packobj%bnd_da()
727  deallocate (packobj)
728  end do
729 
730  ! Scalars
731  call mem_deallocate(this%infmi)
732  call mem_deallocate(this%inmip)
733  call mem_deallocate(this%inadv)
734  call mem_deallocate(this%indsp)
735  call mem_deallocate(this%inssm)
736  call mem_deallocate(this%inmst)
737  call mem_deallocate(this%inmvt)
738  call mem_deallocate(this%inoc)
739 
740  ! Arrays
741  call mem_deallocate(this%masssto)
742  call mem_deallocate(this%massstoold)
743  call mem_deallocate(this%ratesto)
744 
745  deallocate (this%trackctl)
746 
747  call this%NumericalModelType%model_da()
748  end subroutine prt_da
749 
750  !> @brief Allocate memory for scalars
751  subroutine allocate_scalars(this, modelname)
752  ! dummy
753  class(prtmodeltype) :: this
754  character(len=*), intent(in) :: modelname
755 
756  ! allocate members from parent class
757  call this%NumericalModelType%allocate_scalars(modelname)
758 
759  ! allocate members that are part of model class
760  call mem_allocate(this%infmi, 'INFMI', this%memoryPath)
761  call mem_allocate(this%inmip, 'INMIP', this%memoryPath)
762  call mem_allocate(this%inmvt, 'INMVT', this%memoryPath)
763  call mem_allocate(this%inmst, 'INMST', this%memoryPath)
764  call mem_allocate(this%inadv, 'INADV', this%memoryPath)
765  call mem_allocate(this%indsp, 'INDSP', this%memoryPath)
766  call mem_allocate(this%inssm, 'INSSM', this%memoryPath)
767  call mem_allocate(this%inoc, 'INOC ', this%memoryPath)
768 
769  this%infmi = 0
770  this%inmip = 0
771  this%inmvt = 0
772  this%inmst = 0
773  this%inadv = 0
774  this%indsp = 0
775  this%inssm = 0
776  this%inoc = 0
777  end subroutine allocate_scalars
778 
779  !> @brief Allocate arrays
780  subroutine allocate_arrays(this)
782  class(prtmodeltype) :: this
783  integer(I4B) :: n
784 
785  ! Allocate arrays in parent type
786  this%nja = this%dis%nja
787  call this%NumericalModelType%allocate_arrays()
788 
789  ! Allocate and initialize arrays
790  call mem_allocate(this%masssto, this%dis%nodes, &
791  'MASSSTO', this%memoryPath)
792  call mem_allocate(this%massstoold, this%dis%nodes, &
793  'MASSSTOOLD', this%memoryPath)
794  call mem_allocate(this%ratesto, this%dis%nodes, &
795  'RATESTO', this%memoryPath)
796  ! explicit model, so these must be manually allocated
797  call mem_allocate(this%x, this%dis%nodes, 'X', this%memoryPath)
798  call mem_allocate(this%rhs, this%dis%nodes, 'RHS', this%memoryPath)
799  call mem_allocate(this%ibound, this%dis%nodes, 'IBOUND', this%memoryPath)
800  do n = 1, this%dis%nodes
801  this%masssto(n) = dzero
802  this%massstoold(n) = dzero
803  this%ratesto(n) = dzero
804  this%x(n) = dzero
805  this%rhs(n) = dzero
806  this%ibound(n) = 1
807  end do
808  end subroutine allocate_arrays
809 
810  !> @brief Create boundary condition packages for this model
811  subroutine package_create(this, filtyp, ipakid, ipaknum, pakname, mempath, &
812  inunit, iout)
813  ! modules
814  use constantsmodule, only: linelength
815  use prtprpmodule, only: prp_create
816  use apimodule, only: api_create
817  ! dummy
818  class(prtmodeltype) :: this
819  character(len=*), intent(in) :: filtyp
820  character(len=LINELENGTH) :: errmsg
821  integer(I4B), intent(in) :: ipakid
822  integer(I4B), intent(in) :: ipaknum
823  character(len=*), intent(in) :: pakname
824  character(len=*), intent(in) :: mempath
825  integer(I4B), intent(in) :: inunit
826  integer(I4B), intent(in) :: iout
827  ! local
828  class(bndtype), pointer :: packobj
829  class(bndtype), pointer :: packobj2
830  integer(I4B) :: ip
831 
832  ! This part creates the package object
833  select case (filtyp)
834  case ('PRP6')
835  call prp_create(packobj, ipakid, ipaknum, inunit, iout, &
836  this%name, pakname, this%fmi)
837  case ('API6')
838  call api_create(packobj, ipakid, ipaknum, inunit, iout, &
839  this%name, pakname)
840  case default
841  write (errmsg, *) 'Invalid package type: ', filtyp
842  call store_error(errmsg, terminate=.true.)
843  end select
844 
845  ! Packages is the bndlist that is associated with the parent model
846  ! The following statement puts a pointer to this package in the ipakid
847  ! position of packages.
848  do ip = 1, this%bndlist%Count()
849  packobj2 => getbndfromlist(this%bndlist, ip)
850  if (packobj2%packName == pakname) then
851  write (errmsg, '(a,a)') 'Cannot create package. Package name '// &
852  'already exists: ', trim(pakname)
853  call store_error(errmsg, terminate=.true.)
854  end if
855  end do
856  call addbndtolist(this%bndlist, packobj)
857  end subroutine package_create
858 
859  !> @brief Check to make sure required input files have been specified
860  subroutine ftype_check(this, indis)
861  ! dummy
862  class(prtmodeltype) :: this
863  integer(I4B), intent(in) :: indis
864  ! local
865  character(len=LINELENGTH) :: errmsg
866 
867  ! Check for DIS(u) and MIP. Stop if not present.
868  if (indis == 0) then
869  write (errmsg, '(1x,a)') &
870  'Discretization (DIS6, DISV6, or DISU6) package not specified.'
871  call store_error(errmsg)
872  end if
873  if (this%inmip == 0) then
874  write (errmsg, '(1x,a)') &
875  'Model input (MIP6) package not specified.'
876  call store_error(errmsg)
877  end if
878 
879  if (count_errors() > 0) then
880  write (errmsg, '(1x,a)') 'One or more required package(s) not specified.'
881  call store_error(errmsg)
882  call store_error_filename(this%filename)
883  end if
884  end subroutine ftype_check
885 
886  !> @brief Solve the model
887  subroutine prt_solve(this)
888  ! modules
890  use prtprpmodule, only: prtprptype
892  ! dummy variables
893  class(prtmodeltype) :: this
894  ! local variables
895  integer(I4B) :: np, ip
896  class(bndtype), pointer :: packobj
897  type(particletype), pointer :: particle
898  real(DP) :: tmax
899  integer(I4B) :: iprp
900 
901  call create_particle(particle)
902 
903  ! Apply tracking solution to PRP packages
904  iprp = 0
905  do ip = 1, this%bndlist%Count()
906  packobj => getbndfromlist(this%bndlist, ip)
907  select type (packobj)
908  type is (prtprptype)
909  ! Update PRP index
910  iprp = iprp + 1
911 
912  ! Initialize PRP track files
913  if (kper == 1 .and. kstp == 1) then
914  if (packobj%itrkout > 0) then
915  call this%trackctl%init_track_file( &
916  packobj%itrkout, &
917  iprp=iprp)
918  end if
919  if (packobj%itrkcsv > 0) then
920  call this%trackctl%init_track_file( &
921  packobj%itrkcsv, &
922  csv=.true., &
923  iprp=iprp)
924  end if
925  end if
926 
927  ! Track particles
928  do np = 1, packobj%nparticles
929  ! Load particle from storage
930  call packobj%particles%get(particle, this%id, iprp, np)
931 
932  ! If particle is permanently unreleased, cycle.
933  ! Save a termination record if we haven't yet.
934  ! TODO: when we have generic dynamic vectors,
935  ! consider terminating permanently unreleased
936  ! in PRP instead of here. For now, status -8
937  ! indicates the permanently unreleased event
938  ! is not yet recorded, status 8 it has been.
939  if (particle%istatus == (-1 * term_unreleased)) then
940  particle%istatus = term_unreleased
941  call this%method%save(particle, reason=3)
942  call packobj%particles%put(particle, np)
943  end if
944 
945  ! Skip terminated particles
946  if (particle%istatus > active) cycle
947 
948  ! If particle was released this time step, record release
949  particle%istatus = active
950  if (particle%trelease >= totimc) &
951  call this%method%save(particle, reason=0)
952 
953  ! Maximum time is the end of the time step or the particle
954  ! stop time, whichever comes first, unless it's the final
955  ! time step and the extend option is on, in which case
956  ! it's just the particle stop time.
957  if (endofsimulation .and. particle%iextend > 0) then
958  tmax = particle%tstop
959  else
960  tmax = min(totimc + delt, particle%tstop)
961  end if
962 
963  ! Get and apply the tracking method
964  call this%method%apply(particle, tmax)
965 
966  ! Reset previous cell and zone numbers
967  particle%icp = 0
968  particle%izp = 0
969 
970  ! If the particle timed out, terminate it.
971  ! "Timeout" means it remains active, but
972  ! - it reached its stop time, or
973  ! - the simulation is over and no extending.
974  ! We can't detect timeout within the tracking
975  ! method because the method just receives the
976  ! maximum time with no context on what it is.
977  ! TODO maybe think about changing that?
978  if (particle%istatus <= active .and. &
979  (particle%ttrack == particle%tstop .or. &
980  (endofsimulation .and. particle%iextend == 0))) then
981  particle%istatus = term_timeout
982  call this%method%save(particle, reason=3)
983  end if
984 
985  ! Update particle storage
986  call packobj%particles%put(particle, np)
987  end do
988  end select
989  end do
990 
991  deallocate (particle)
992  end subroutine prt_solve
993 
994  !> @brief Source package info and begin to process
995  subroutine create_bndpkgs(this, bndpkgs, pkgtypes, pkgnames, &
996  mempaths, inunits)
997  ! modules
1000  ! dummy
1001  class(prtmodeltype) :: this
1002  integer(I4B), dimension(:), allocatable, intent(inout) :: bndpkgs
1003  type(characterstringtype), dimension(:), contiguous, &
1004  pointer, intent(inout) :: pkgtypes
1005  type(characterstringtype), dimension(:), contiguous, &
1006  pointer, intent(inout) :: pkgnames
1007  type(characterstringtype), dimension(:), contiguous, &
1008  pointer, intent(inout) :: mempaths
1009  integer(I4B), dimension(:), contiguous, &
1010  pointer, intent(inout) :: inunits
1011  ! local
1012  integer(I4B) :: ipakid, ipaknum
1013  character(len=LENFTYPE) :: pkgtype, bndptype
1014  character(len=LENPACKAGENAME) :: pkgname
1015  character(len=LENMEMPATH) :: mempath
1016  integer(I4B), pointer :: inunit
1017  integer(I4B) :: n
1018 
1019  if (allocated(bndpkgs)) then
1020  !
1021  ! create stress packages
1022  ipakid = 1
1023  bndptype = ''
1024  do n = 1, size(bndpkgs)
1025  !
1026  pkgtype = pkgtypes(bndpkgs(n))
1027  pkgname = pkgnames(bndpkgs(n))
1028  mempath = mempaths(bndpkgs(n))
1029  inunit => inunits(bndpkgs(n))
1030  !
1031  if (bndptype /= pkgtype) then
1032  ipaknum = 1
1033  bndptype = pkgtype
1034  end if
1035  !
1036  call this%package_create(pkgtype, ipakid, ipaknum, pkgname, mempath, &
1037  inunit, this%iout)
1038  ipakid = ipakid + 1
1039  ipaknum = ipaknum + 1
1040  end do
1041  !
1042  ! cleanup
1043  deallocate (bndpkgs)
1044  end if
1045 
1046  end subroutine create_bndpkgs
1047 
1048  !> @brief Source package info and begin to process
1049  subroutine create_packages(this)
1050  ! modules
1053  use arrayhandlersmodule, only: expandarray
1054  use memorymanagermodule, only: mem_setptr
1056  use simvariablesmodule, only: idm_context
1057  use budgetmodule, only: budget_cr
1061  use prtmipmodule, only: mip_cr
1062  use prtfmimodule, only: fmi_cr
1063  use prtocmodule, only: oc_cr
1064  ! dummy
1065  class(prtmodeltype) :: this
1066  ! local
1067  type(characterstringtype), dimension(:), contiguous, &
1068  pointer :: pkgtypes => null()
1069  type(characterstringtype), dimension(:), contiguous, &
1070  pointer :: pkgnames => null()
1071  type(characterstringtype), dimension(:), contiguous, &
1072  pointer :: mempaths => null()
1073  integer(I4B), dimension(:), contiguous, &
1074  pointer :: inunits => null()
1075  character(len=LENMEMPATH) :: model_mempath
1076  character(len=LENFTYPE) :: pkgtype
1077  character(len=LENPACKAGENAME) :: pkgname
1078  character(len=LENMEMPATH) :: mempath
1079  integer(I4B), pointer :: inunit
1080  integer(I4B), dimension(:), allocatable :: bndpkgs
1081  integer(I4B) :: n
1082  integer(I4B) :: indis = 0 ! DIS enabled flag
1083  character(len=LENMEMPATH) :: mempathmip = ''
1084 
1085  ! set input memory paths, input/model and input/model/namfile
1086  model_mempath = create_mem_path(component=this%name, context=idm_context)
1087 
1088  ! set pointers to model path package info
1089  call mem_setptr(pkgtypes, 'PKGTYPES', model_mempath)
1090  call mem_setptr(pkgnames, 'PKGNAMES', model_mempath)
1091  call mem_setptr(mempaths, 'MEMPATHS', model_mempath)
1092  call mem_setptr(inunits, 'INUNITS', model_mempath)
1093 
1094  do n = 1, size(pkgtypes)
1095  ! attributes for this input package
1096  pkgtype = pkgtypes(n)
1097  pkgname = pkgnames(n)
1098  mempath = mempaths(n)
1099  inunit => inunits(n)
1100 
1101  ! create dis package first as it is a prerequisite for other packages
1102  select case (pkgtype)
1103  case ('DIS6')
1104  indis = 1
1105  call dis_cr(this%dis, this%name, mempath, indis, this%iout)
1106  case ('DISV6')
1107  indis = 1
1108  call disv_cr(this%dis, this%name, mempath, indis, this%iout)
1109  case ('DISU6')
1110  indis = 1
1111  call disu_cr(this%dis, this%name, mempath, indis, this%iout)
1112  case ('MIP6')
1113  this%inmip = 1
1114  mempathmip = mempath
1115  case ('FMI6')
1116  this%infmi = inunit
1117  case ('OC6')
1118  this%inoc = inunit
1119  case ('PRP6')
1120  call expandarray(bndpkgs)
1121  bndpkgs(size(bndpkgs)) = n
1122  case default
1123  call pstop(1, "Unrecognized package type: "//pkgtype)
1124  end select
1125  end do
1126 
1127  ! Create budget manager
1128  call budget_cr(this%budget, this%name)
1129 
1130  ! Create tracking method pools
1131  call create_method_pool()
1134 
1135  ! Create packages that are tied directly to model
1136  call mip_cr(this%mip, this%name, mempathmip, this%inmip, this%iout, this%dis)
1137  call fmi_cr(this%fmi, this%name, this%infmi, this%iout)
1138  call oc_cr(this%oc, this%name, this%inoc, this%iout)
1139 
1140  ! Check to make sure that required ftype's have been specified
1141  call this%ftype_check(indis)
1142 
1143  ! Create boundary packages
1144  call this%create_bndpkgs(bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
1145  end subroutine create_packages
1146 
1147  !> @brief Write model namfile options to list file
1148  subroutine log_namfile_options(this, found)
1150  class(prtmodeltype) :: this
1151  type(gwfnamparamfoundtype), intent(in) :: found
1152 
1153  write (this%iout, '(1x,a)') 'NAMEFILE OPTIONS:'
1154 
1155  if (found%newton) then
1156  write (this%iout, '(4x,a)') &
1157  'NEWTON-RAPHSON method enabled for the model.'
1158  if (found%under_relaxation) then
1159  write (this%iout, '(4x,a,a)') &
1160  'NEWTON-RAPHSON UNDER-RELAXATION based on the bottom ', &
1161  'elevation of the model will be applied to the model.'
1162  end if
1163  end if
1164 
1165  if (found%print_input) then
1166  write (this%iout, '(4x,a)') 'STRESS PACKAGE INPUT WILL BE PRINTED '// &
1167  'FOR ALL MODEL STRESS PACKAGES'
1168  end if
1169 
1170  if (found%print_flows) then
1171  write (this%iout, '(4x,a)') 'PACKAGE FLOWS WILL BE PRINTED '// &
1172  'FOR ALL MODEL PACKAGES'
1173  end if
1174 
1175  if (found%save_flows) then
1176  write (this%iout, '(4x,a)') &
1177  'FLOWS WILL BE SAVED TO BUDGET FILE SPECIFIED IN OUTPUT CONTROL'
1178  end if
1179 
1180  write (this%iout, '(1x,a)') 'END NAMEFILE OPTIONS:'
1181  end subroutine log_namfile_options
1182 
1183 end module prtmodule
This module contains the API package methods.
Definition: gwf-api.f90:12
subroutine, public api_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname)
@ brief Create a new package object
Definition: gwf-api.f90:49
subroutine, public addbasemodeltolist(list, model)
Definition: BaseModel.f90:161
This module contains the base boundary package.
subroutine, public addbndtolist(list, bnd)
Add boundary to package list.
class(bndtype) function, pointer, public getbndfromlist(list, idx)
Get boundary from package list.
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
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:632
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
@ mnormal
normal output mode
Definition: Constants.f90:206
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter lenpackagetype
maximum length of a package type (DIS6, SFR6, CSUB6, etc.)
Definition: Constants.f90:38
integer(i4b), parameter lenbigline
maximum length of a big line
Definition: Constants.f90:15
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
integer(i4b), parameter lenpakloc
maximum length of a package location
Definition: Constants.f90:50
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:37
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Definition: Dis.f90:1
subroutine, public dis_cr(dis, name_model, input_mempath, inunit, iout)
Create a new structured discretization object.
Definition: Dis.f90:97
subroutine, public disu_cr(dis, name_model, input_mempath, inunit, iout)
Create a new unstructured discretization object.
Definition: Disu.f90:127
subroutine, public disv_cr(dis, name_model, input_mempath, inunit, iout)
Create a new discretization by vertices object.
Definition: Disv.f90:109
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
subroutine, public lowcase(word)
Convert to lower case.
subroutine, public parseline(line, nwords, words, inunit, filename)
Parse a line into words.
subroutine, public upcase(word)
Convert to upper case.
This module defines variable data types.
Definition: kind.f90:8
type(listtype), public basemodellist
Definition: mf6lists.f90:16
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public memorystore_remove(component, subcomponent, context)
Cell-level tracking methods.
subroutine, public create_method_cell_pool()
Create the cell method pool.
subroutine, public destroy_method_cell_pool()
Destroy the cell method pool.
Particle tracking strategies.
Definition: Method.f90:2
Model-level tracking methods.
Definition: MethodPool.f90:2
type(methoddisvtype), pointer, public method_disv
Definition: MethodPool.f90:12
type(methoddistype), pointer, public method_dis
Definition: MethodPool.f90:11
subroutine, public create_method_pool()
Create the method pool.
Definition: MethodPool.f90:18
subroutine, public destroy_method_pool()
Destroy the method pool.
Definition: MethodPool.f90:24
Subcell-level tracking methods.
subroutine, public create_method_subcell_pool()
Create the subcell method pool.
subroutine, public destroy_method_subcell_pool()
Destroy the subcell method pool.
@ term_timeout
terminated at stop time or end of simulation
Definition: Particle.f90:40
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:38
subroutine create_particle(particle)
Create a new particle.
Definition: Particle.f90:164
subroutine, public fmi_cr(fmiobj, name_model, inunit, iout)
Create a new PrtFmi object.
Definition: prt-fmi.f90:38
subroutine, public mip_cr(mip, name_model, input_mempath, inunit, iout, dis)
Create a model input object.
Definition: prt-mip.f90:34
Definition: prt.f90:1
integer(i4b), parameter niunit_prt
Definition: prt.f90:114
subroutine prt_ot(this)
Print and/or save model output.
Definition: prt.f90:490
subroutine prt_rp(this)
Read and prepare (calls package read and prepare routines)
Definition: prt.f90:299
subroutine create_bndpkgs(this, bndpkgs, pkgtypes, pkgnames, mempaths, inunits)
Source package info and begin to process.
Definition: prt.f90:997
subroutine prt_ar(this)
Allocate and read.
Definition: prt.f90:228
subroutine ftype_check(this, indis)
Check to make sure required input files have been specified.
Definition: prt.f90:861
subroutine prt_ot_saveflow(this, nja, flowja, icbcfl, icbcun)
Save intercell flows.
Definition: prt.f90:572
subroutine prt_ad(this)
Time step advance (calls package advance subroutines)
Definition: prt.f90:319
subroutine prt_cq(this, icnvg, isuppress_output)
Calculate intercell flow (flowja)
Definition: prt.f90:362
subroutine package_create(this, filtyp, ipakid, ipaknum, pakname, mempath, inunit, iout)
Create boundary condition packages for this model.
Definition: prt.f90:813
subroutine prt_ot_flow(this, icbcfl, ibudfl, icbcun)
Save flows.
Definition: prt.f90:535
subroutine allocate_scalars(this, modelname)
Allocate memory for scalars.
Definition: prt.f90:752
subroutine prt_ot_bdsummary(this, ibudfl, ipflag)
Print budget summary.
Definition: prt.f90:659
character(len=lenpackagetype), dimension(prt_nmultipkg), public prt_multipkg
Definition: prt.f90:109
subroutine create_packages(this)
Source package info and begin to process.
Definition: prt.f90:1050
character(len=lenpackagetype), dimension(prt_nbasepkg), public prt_basepkg
Definition: prt.f90:95
integer(i4b), parameter, public prt_nmultipkg
PRT multi package array descriptors.
Definition: prt.f90:108
character(len=lenbudtxt), dimension(nbditems) budtxt
Definition: prt.f90:37
subroutine prt_da(this)
Deallocate.
Definition: prt.f90:689
subroutine prt_cq_sto(this)
Calculate particle mass storage.
Definition: prt.f90:408
subroutine, public prt_cr(filename, id, modelname)
Create a new particle tracking model object.
Definition: prt.f90:120
subroutine prt_ot_printflow(this, ibudfl, flowja)
Print intercell flows.
Definition: prt.f90:599
subroutine prt_bd(this, icnvg, isuppress_output)
Calculate flows and budget.
Definition: prt.f90:461
subroutine prt_df(this)
Define packages.
Definition: prt.f90:195
integer(i4b), parameter, public prt_nbasepkg
PRT base package array descriptors.
Definition: prt.f90:94
integer(i4b), parameter nbditems
Definition: prt.f90:36
subroutine allocate_arrays(this)
Allocate arrays.
Definition: prt.f90:781
subroutine log_namfile_options(this, found)
Write model namfile options to list file.
Definition: prt.f90:1149
subroutine prt_ot_dv(this, idvsave, idvprint, ipflag)
Print dependent variables.
Definition: prt.f90:638
subroutine prt_solve(this)
Solve the model.
Definition: prt.f90:888
subroutine, public oc_cr(ocobj, name_model, inunit, iout)
@ brief Create an output control object
Definition: prt-oc.f90:51
subroutine, public prp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi)
Create a new particle release point package.
Definition: prt-prp.f90:102
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
integer(i4b) isimcheck
simulation input check flag (1) to check input, (0) to ignore checks
integer(i4b) ifailedstepretry
current retry for this time step
subroutine csr_diagsum(ia, flowja)
Definition: Sparse.f90:263
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
subroutine, public tdis_ot(iout)
Print simulation time.
Definition: tdis.f90:274
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:26
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
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
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
This module contains version information.
Definition: version.f90:7
subroutine write_listfile_header(iout, cmodel_type, write_sys_command, write_kind_info)
@ brief Write program header
Definition: version.f90:98
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
@ brief BndType
Derived type for the Budget object.
Definition: Budget.f90:39
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Structured grid discretization.
Definition: Dis.f90:23
Unstructured grid discretization.
Definition: Disu.f90:28
Vertex grid discretization.
Definition: Disv.f90:24
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
Base type for particle tracking methods.
Definition: Method.f90:31
Particle tracked by the PRT model.
Definition: Particle.f90:78
Particle tracking (PRT) model.
Definition: prt.f90:41
@ brief Output control for particle tracking models
Definition: prt-oc.f90:22
Particle release point (PRP) package.
Definition: prt-prp.f90:41
Manages particle track (i.e. pathline) files.
Output file containing all or some particle pathlines.
Definition: TrackFile.f90:28