MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
prt-prp.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
6  use bndmodule, only: bndtype
7  use bndextmodule, only: bndexttype
9  use prtfmimodule, only: prtfmitype
22  use dismodule, only: distype
23  use disvmodule, only: disvtype
24  use errorutilmodule, only: pstop
25  use mathutilmodule, only: arange, is_close
27 
28  implicit none
29 
30  private
31  public :: prtprptype
32  public :: prp_create
33 
34  character(len=LENFTYPE) :: ftype = 'PRP'
35  character(len=16) :: text = ' PRP'
36  real(dp), parameter :: default_exit_solve_tolerance = dem5
37 
38  !> @brief Particle release point (PRP) package
39  type, extends(bndexttype) :: prtprptype
40  ! options
41  logical(LGP), pointer :: extend => null() !< extend tracking beyond simulation's end
42  logical(LGP), pointer :: frctrn => null() !< force ternary solution for quad grids
43  logical(LGP), pointer :: drape => null() !< whether to drape particle to topmost active cell
44  logical(LGP), pointer :: localz => null() !< compute z coordinates local to the release cell
45  integer(I4B), pointer :: istopweaksink => null() !< weak sink option: 0 = no stop, 1 = stop
46  integer(I4B), pointer :: istopzone => null() !< optional stop zone number: 0 = no stop zone
47  integer(I4B), pointer :: idrymeth => null() !< dry tracking method: 0 = drop, 1 = stop, 2 = stay
48  integer(I4B), pointer :: itrkout => null() !< binary track file
49  integer(I4B), pointer :: itrkhdr => null() !< track header file
50  integer(I4B), pointer :: itrkcsv => null() !< CSV track file
51  integer(I4B), pointer :: irlstls => null() !< release time file
52  integer(I4B), pointer :: iexmeth => null() !< method for iterative solution of particle exit location and time in generalized Pollock's method
53  integer(I4B), pointer :: ichkmeth => null() !< method for checking particle release coordinates are in the specified cells, 0 = none, 1 = eager
54  integer(I4B), pointer :: icycwin => null() !< cycle detection window size
55  real(dp), pointer :: extol => null() !< tolerance for iterative solution of particle exit location and time in generalized Pollock's method
56  real(dp), pointer :: rttol => null() !< tolerance for coincident particle release times
57  real(dp), pointer :: rtfreq => null() !< frequency for regularly spaced release times
58  real(dp), pointer :: offset => null() !< release time offset
59  real(dp), pointer :: stoptime => null() !< stop time for all release points
60  real(dp), pointer :: stoptraveltime => null() !< stop travel time for all points
61  !
62  type(prtfmitype), pointer :: fmi => null() !< flow model interface
63  type(particlestoretype), pointer :: particles => null() !< particle store
64  type(particlereleasescheduletype), pointer :: schedule => null() !< particle release schedule
65  integer(I4B), pointer :: nreleasepoints => null() !< number of release points
66  integer(I4B), pointer :: nreleasetimes => null() !< number of user-specified particle release times
67  integer(I4B), pointer :: nparticles => null() !< number of particles released
68  integer(I4B), pointer, contiguous :: rptnode(:) => null() !< release point reduced nns
69  integer(I4B), pointer, contiguous :: rptzone(:) => null() !< release point zone numbers
70  real(dp), pointer, contiguous :: rptx(:) => null() !< release point x coordinates
71  real(dp), pointer, contiguous :: rpty(:) => null() !< release point y coordinates
72  real(dp), pointer, contiguous :: rptz(:) => null() !< release point z coordinates
73  real(dp), pointer, contiguous :: rptm(:) => null() !< total mass released from point
74  character(len=LENBOUNDNAME), pointer, contiguous :: rptname(:) => null() !< release point names
75  contains
76  procedure :: prp_allocate_arrays
77  procedure :: prp_allocate_scalars
78  procedure :: bnd_ar => prp_ar
79  procedure :: bnd_ad => prp_ad
80  procedure :: bnd_rp => prp_rp
81  procedure :: bnd_cq_simrate => prp_cq_simrate
82  procedure :: bnd_da => prp_da
83  procedure :: define_listlabel
84  procedure :: prp_set_pointers
85  procedure :: source_options => prp_options
86  procedure :: source_dimensions => prp_dimensions
87  procedure :: prp_log_options
88  procedure :: prp_packagedata
89  procedure :: prp_releasetimes
91  procedure :: release
92  procedure :: log_release
94  procedure :: initialize_particle
95  procedure, public :: bnd_obs_supported => prp_obs_supported
96  procedure, public :: bnd_df_obs => prp_df_obs
97  end type prtprptype
98 
99 contains
100 
101  !> @brief Create a new particle release point package
102  subroutine prp_create(packobj, id, ibcnum, inunit, iout, namemodel, &
103  pakname, input_mempath, fmi)
104  ! dummy
105  class(bndtype), pointer :: packobj
106  integer(I4B), intent(in) :: id
107  integer(I4B), intent(in) :: ibcnum
108  integer(I4B), intent(in) :: inunit
109  integer(I4B), intent(in) :: iout
110  character(len=*), intent(in) :: namemodel
111  character(len=*), intent(in) :: pakname
112  character(len=*), intent(in) :: input_mempath
113  type(prtfmitype), pointer :: fmi
114  ! local
115  type(prtprptype), pointer :: prpobj
116  ! formats
117  character(len=*), parameter :: fmtheader = &
118  "(1x, /1x, 'PRP PARTICLE RELEASE POINT PACKAGE', &
119  &' INPUT READ FROM MEMPATH: ', a, /)"
120 
121  ! allocate the object and assign values to object variables
122  allocate (prpobj)
123  packobj => prpobj
124 
125  ! create name and memory path
126  call packobj%set_names(ibcnum, namemodel, pakname, ftype, input_mempath)
127  prpobj%text = text
128 
129  ! allocate scalars
130  call prpobj%prp_allocate_scalars()
131 
132  ! initialize package
133  call packobj%pack_initialize()
134 
135  packobj%inunit = inunit
136  packobj%iout = iout
137  packobj%id = id
138  packobj%ibcnum = ibcnum
139  packobj%ncolbnd = 4
140  packobj%iscloc = 1
141 
142  ! store pointer to flow model interface
143  prpobj%fmi => fmi
144 
145  ! if prp is enabled, print a message identifying it
146  if (inunit > 0) write (iout, fmtheader) input_mempath
147  end subroutine prp_create
148 
149  !> @brief Deallocate memory
150  subroutine prp_da(this)
151  class(prtprptype) :: this
152 
153  ! Deallocate parent
154  call this%BndExtType%bnd_da()
155 
156  ! Deallocate scalars
157  call mem_deallocate(this%localz)
158  call mem_deallocate(this%extend)
159  call mem_deallocate(this%offset)
160  call mem_deallocate(this%stoptime)
161  call mem_deallocate(this%stoptraveltime)
162  call mem_deallocate(this%istopweaksink)
163  call mem_deallocate(this%istopzone)
164  call mem_deallocate(this%drape)
165  call mem_deallocate(this%idrymeth)
166  call mem_deallocate(this%nreleasepoints)
167  call mem_deallocate(this%nreleasetimes)
168  call mem_deallocate(this%nparticles)
169  call mem_deallocate(this%itrkout)
170  call mem_deallocate(this%itrkhdr)
171  call mem_deallocate(this%itrkcsv)
172  call mem_deallocate(this%irlstls)
173  call mem_deallocate(this%frctrn)
174  call mem_deallocate(this%iexmeth)
175  call mem_deallocate(this%ichkmeth)
176  call mem_deallocate(this%icycwin)
177  call mem_deallocate(this%extol)
178  call mem_deallocate(this%rttol)
179  call mem_deallocate(this%rtfreq)
180 
181  ! Deallocate arrays
182  call mem_deallocate(this%rptx)
183  call mem_deallocate(this%rpty)
184  call mem_deallocate(this%rptz)
185  call mem_deallocate(this%rptnode)
186  call mem_deallocate(this%rptm)
187  call mem_deallocate(this%rptname, 'RPTNAME', this%memoryPath)
188 
189  ! Deallocate objects
190  call this%particles%destroy(this%memoryPath)
191  call this%schedule%destroy()
192  deallocate (this%particles)
193  deallocate (this%schedule)
194  end subroutine prp_da
195 
196  !> @ brief Set pointers to model variables
197  subroutine prp_set_pointers(this, ibound, izone)
198  class(prtprptype) :: this
199  integer(I4B), dimension(:), pointer, contiguous :: ibound
200  integer(I4B), dimension(:), pointer, contiguous :: izone
201 
202  this%ibound => ibound
203  this%rptzone => izone
204  end subroutine prp_set_pointers
205 
206  !> @brief Allocate arrays
207  subroutine prp_allocate_arrays(this, nodelist, auxvar)
208  ! dummy
209  class(prtprptype) :: this
210  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
211  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
212  ! local
213  integer(I4B) :: nps
214 
215  call this%BndExtType%allocate_arrays()
216 
217  ! Allocate particle store, starting with the number
218  ! of release points (arrays resized if/when needed)
219  call create_particle_store( &
220  this%particles, &
221  this%nreleasepoints, &
222  this%memoryPath)
223 
224  ! Allocate arrays
225  call mem_allocate(this%rptx, this%nreleasepoints, 'RPTX', this%memoryPath)
226  call mem_allocate(this%rpty, this%nreleasepoints, 'RPTY', this%memoryPath)
227  call mem_allocate(this%rptz, this%nreleasepoints, 'RPTZ', this%memoryPath)
228  call mem_allocate(this%rptm, this%nreleasepoints, 'RPTMASS', &
229  this%memoryPath)
230  call mem_allocate(this%rptnode, this%nreleasepoints, 'RPTNODER', &
231  this%memoryPath)
232  call mem_allocate(this%rptname, lenboundname, this%nreleasepoints, &
233  'RPTNAME', this%memoryPath)
234 
235  ! Initialize arrays
236  do nps = 1, this%nreleasepoints
237  this%rptm(nps) = dzero
238  end do
239  end subroutine prp_allocate_arrays
240 
241  !> @brief Allocate scalars
242  subroutine prp_allocate_scalars(this)
243  class(prtprptype) :: this
244 
245  ! Allocate parent's scalars
246  call this%BndExtType%allocate_scalars()
247 
248  ! Allocate scalars for this type
249  call mem_allocate(this%localz, 'LOCALZ', this%memoryPath)
250  call mem_allocate(this%extend, 'EXTEND', this%memoryPath)
251  call mem_allocate(this%offset, 'OFFSET', this%memoryPath)
252  call mem_allocate(this%stoptime, 'STOPTIME', this%memoryPath)
253  call mem_allocate(this%stoptraveltime, 'STOPTRAVELTIME', this%memoryPath)
254  call mem_allocate(this%istopweaksink, 'ISTOPWEAKSINK', this%memoryPath)
255  call mem_allocate(this%istopzone, 'ISTOPZONE', this%memoryPath)
256  call mem_allocate(this%drape, 'DRAPE', this%memoryPath)
257  call mem_allocate(this%idrymeth, 'IDRYMETH', this%memoryPath)
258  call mem_allocate(this%nreleasepoints, 'NRELEASEPOINTS', this%memoryPath)
259  call mem_allocate(this%nreleasetimes, 'NRELEASETIMES', this%memoryPath)
260  call mem_allocate(this%nparticles, 'NPARTICLES', this%memoryPath)
261  call mem_allocate(this%itrkout, 'ITRKOUT', this%memoryPath)
262  call mem_allocate(this%itrkhdr, 'ITRKHDR', this%memoryPath)
263  call mem_allocate(this%itrkcsv, 'ITRKCSV', this%memoryPath)
264  call mem_allocate(this%irlstls, 'IRLSTLS', this%memoryPath)
265  call mem_allocate(this%frctrn, 'FRCTRN', this%memoryPath)
266  call mem_allocate(this%iexmeth, 'IEXMETH', this%memoryPath)
267  call mem_allocate(this%ichkmeth, 'ICHKMETH', this%memoryPath)
268  call mem_allocate(this%icycwin, 'ICYCWIN', this%memoryPath)
269  call mem_allocate(this%extol, 'EXTOL', this%memoryPath)
270  call mem_allocate(this%rttol, 'RTTOL', this%memoryPath)
271  call mem_allocate(this%rtfreq, 'RTFREQ', this%memoryPath)
272 
273  ! Set values
274  this%localz = .false.
275  this%extend = .false.
276  this%offset = dzero
277  this%stoptime = huge(1d0)
278  this%stoptraveltime = huge(1d0)
279  this%istopweaksink = 0
280  this%istopzone = 0
281  this%drape = .false.
282  this%idrymeth = 0
283  this%nreleasepoints = 0
284  this%nreleasetimes = 0
285  this%nparticles = 0
286  this%itrkout = 0
287  this%itrkhdr = 0
288  this%itrkcsv = 0
289  this%irlstls = 0
290  this%frctrn = .false.
291  this%iexmeth = 0
292  this%ichkmeth = 1
293  this%icycwin = 0
294  this%extol = default_exit_solve_tolerance
295  this%rttol = dsame * dep9
296  this%rtfreq = dzero
297 
298  end subroutine prp_allocate_scalars
299 
300  !> @ brief Allocate and read period data
301  subroutine prp_ar(this)
302  ! dummy variables
303  class(prtprptype), intent(inout) :: this
304  ! local variables
305  integer(I4B) :: n
306 
307  call this%obs%obs_ar()
308 
309  if (this%inamedbound /= 0) then
310  do n = 1, this%nreleasepoints
311  this%boundname(n) = this%rptname(n)
312  end do
313  end if
314  do n = 1, this%nreleasepoints
315  this%nodelist(n) = this%rptnode(n)
316  end do
317  end subroutine prp_ar
318 
319  !> @brief Advance a time step and release particles if scheduled.
320  subroutine prp_ad(this)
321  use tdismodule, only: totalsimtime
322  class(prtprptype) :: this
323  integer(I4B) :: ip, it
324  real(DP) :: t
325 
326  ! Notes
327  ! -----
328  ! Each release point can be thought of as
329  ! a gumball machine with infinite supply:
330  ! a point can release an arbitrary number
331  ! of particles, but only one at any time.
332  ! Coincident release times are merged to
333  ! a single time by the release scheduler.
334 
335  ! Reset mass accumulators for this time step.
336  do ip = 1, this%nreleasepoints
337  this%rptm(ip) = dzero
338  end do
339 
340  ! Advance the release schedule and check if
341  ! any releases will be made this time step.
342  call this%schedule%advance()
343  if (.not. this%schedule%any()) return
344 
345  ! Log the schedule to the list file.
346  call this%log_release()
347 
348  ! Expand the particle store. We know from the
349  ! schedule how many particles will be released.
350  call this%particles%resize( &
351  this%particles%num_stored() + &
352  (this%nreleasepoints * this%schedule%count()), &
353  this%memoryPath)
354 
355  ! Release a particle from each point for
356  ! each release time in the current step.
357  do ip = 1, this%nreleasepoints
358  do it = 1, this%schedule%count()
359  t = this%schedule%times(it)
360  ! Skip the release time if it's before the simulation
361  ! starts, or if no `extend_tracking`, after it ends.
362  if (t < dzero) then
363  write (warnmsg, '(a,g0,a)') &
364  'Skipping negative release time (t=', t, ').'
365  call store_warning(warnmsg)
366  cycle
367  else if (t > totalsimtime .and. .not. this%extend) then
368  write (warnmsg, '(a,g0,a)') &
369  'Skipping release time falling after the end of the &
370  &simulation (t=', t, '). Enable EXTEND_TRACKING to &
371  &release particles after the simulation end time.'
372  call store_warning(warnmsg)
373  cycle
374  end if
375  call this%release(ip, t)
376  end do
377  end do
378  end subroutine prp_ad
379 
380  !> @brief Log the release scheduled for this time step.
381  subroutine log_release(this)
382  class(prtprptype), intent(inout) :: this !< prp
383  if (this%iprpak > 0) then
384  write (this%iout, "(1x,/1x,a,1x,i0)") &
385  'PARTICLE RELEASE FOR PRP', this%ibcnum
386  call this%schedule%log(this%iout)
387  end if
388  end subroutine log_release
389 
390  !> @brief Verify that the release point is in the cell.
391  !!
392  !! Terminate with an error if the release point lies outside the
393  !! given cell, or if the point is above or below the grid top or
394  !! bottom, respectively.
395  !<
396  subroutine validate_release_point(this, ic, x, y, z)
397  class(prtprptype), intent(inout) :: this !< this instance
398  integer(I4B), intent(in) :: ic !< cell index
399  real(DP), intent(in) :: x, y, z !< release point
400  ! local
401  real(DP), allocatable :: polyverts(:, :)
402 
403  call this%fmi%dis%get_polyverts(ic, polyverts)
404  if (.not. point_in_polygon(x, y, polyverts)) then
405  write (errmsg, '(a,g0,a,g0,a,i0)') &
406  'Error: release point (x=', x, ', y=', y, ') is not in cell ', &
407  this%dis%get_nodeuser(ic)
408  call store_error(errmsg, terminate=.false.)
409  call store_error_filename(this%input_fname)
410  end if
411  if (z > maxval(this%dis%top)) then
412  write (errmsg, '(a,g0,a,g0,a,i0)') &
413  'Error: release point (z=', z, ') is above grid top ', &
414  maxval(this%dis%top)
415  call store_error(errmsg, terminate=.false.)
416  call store_error_filename(this%input_fname)
417  else if (z < minval(this%dis%bot)) then
418  write (errmsg, '(a,g0,a,g0,a,i0)') &
419  'Error: release point (z=', z, ') is below grid bottom ', &
420  minval(this%dis%bot)
421  call store_error(errmsg, terminate=.false.)
422  call store_error_filename(this%input_fname)
423  end if
424  deallocate (polyverts)
425  end subroutine validate_release_point
426 
427  !> Release a particle at the specified time.
428  !!
429  !! Releasing a particle entails validating the particle's
430  !! coordinates and settings, transforming its coordinates
431  !! if needed, initializing the particle's initial tracking
432  !! time to the given release time, storing the particle in
433  !! the particle store (from which the PRT model will later
434  !! retrieve it, apply the tracking method, and check it in
435  !! again), and accumulating the particle's mass (the total
436  !! mass released from each release point is calculated for
437  !! budget reporting).
438  !<
439  subroutine release(this, ip, trelease)
440  ! dummy
441  class(prtprptype), intent(inout) :: this !< this instance
442  integer(I4B), intent(in) :: ip !< particle index
443  real(DP), intent(in) :: trelease !< release time
444  ! local
445  integer(I4B) :: np
446  type(particletype), pointer :: particle
447 
448  call this%initialize_particle(particle, ip, trelease)
449  np = this%nparticles + 1
450  this%nparticles = np
451  call this%particles%put(particle, np)
452  deallocate (particle)
453  this%rptm(ip) = this%rptm(ip) + done ! TODO configurable mass
454 
455  end subroutine release
456 
457  subroutine initialize_particle(this, particle, ip, trelease)
459  class(prtprptype), intent(inout) :: this !< this instance
460  type(particletype), pointer, intent(inout) :: particle !< the particle
461  integer(I4B), intent(in) :: ip !< particle index
462  real(DP), intent(in) :: trelease !< release time
463  ! local
464  integer(I4B) :: irow, icol, ilay, icpl
465  integer(I4B) :: ic, icu, ic_old
466  real(DP) :: x, y, z
467  real(DP) :: top, bot, hds
468  ! formats
469  character(len=*), parameter :: fmticterr = &
470  "('Error in ',a,': Flow model interface does not contain ICELLTYPE. &
471  &ICELLTYPE is required for PRT to distinguish convertible cells &
472  &from confined cells if LOCAL_Z release coordinates are provided. &
473  &Make sure a GWFGRID entry is configured in the PRT FMI package.')"
474 
475  ic = this%rptnode(ip)
476  icu = this%dis%get_nodeuser(ic)
477 
478  call create_particle(particle)
479 
480  if (size(this%boundname) /= 0) then
481  particle%name = this%boundname(ip)
482  else
483  particle%name = ''
484  end if
485 
486  particle%irpt = ip
487  particle%istopweaksink = this%istopweaksink
488  particle%istopzone = this%istopzone
489  particle%idrymeth = this%idrymeth
490  particle%icu = icu
491 
492  select type (dis => this%dis)
493  type is (distype)
494  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
495  type is (disvtype)
496  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
497  end select
498  particle%ilay = ilay
499  particle%izone = this%rptzone(ic)
500  particle%istatus = 0 ! status 0 until tracking starts
501  ! If the cell is inactive, either drape the particle
502  ! to the top-most active cell beneath it if drape is
503  ! enabled, or else terminate permanently unreleased.
504  if (this%ibound(ic) == 0) then
505  ic_old = ic
506  if (this%drape) then
507  call this%dis%highest_active(ic, this%ibound)
508  if (ic == ic_old .or. this%ibound(ic) == 0) then
509  ! negative unreleased status signals to the
510  ! tracking method that we haven't yet saved
511  ! a termination record, it needs to do so.
512  particle%istatus = -1 * term_unreleased
513  end if
514  else
515  particle%istatus = -1 * term_unreleased
516  end if
517  end if
518 
519  ! load coordinates
520  x = this%rptx(ip)
521  y = this%rpty(ip)
522  if (this%localz) then
523  ! make sure FMI has cell type array. we need
524  ! it to distinguish convertible and confined
525  ! cells if release z coordinates are local
526  if (this%fmi%igwfceltyp /= 1) then
527  write (errmsg, fmticterr) trim(this%text)
528  call store_error(errmsg, terminate=.true.)
529  end if
530 
531  ! calculate model z coord from local z coord.
532  ! if cell is confined (icelltype == 0) use the
533  ! actual cell height (geometric top - bottom).
534  ! otherwise use head as cell top, clamping to
535  ! the cell bottom if head is below the bottom
536  top = this%fmi%dis%top(ic)
537  bot = this%fmi%dis%bot(ic)
538  hds = this%fmi%gwfhead(ic)
539  if (this%fmi%gwfceltyp(icu) /= 0) top = hds
540  if (top < bot) top = bot
541  z = bot + this%rptz(ip) * (top - bot)
542  else
543  z = this%rptz(ip)
544  end if
545 
546  if (this%ichkmeth > 0) &
547  call this%validate_release_point(ic, x, y, z)
548 
549  particle%x = x
550  particle%y = y
551  particle%z = z
552  particle%trelease = trelease
553 
554  ! Set stop time to earlier of STOPTIME and STOPTRAVELTIME
555  if (this%stoptraveltime == huge(1d0)) then
556  particle%tstop = this%stoptime
557  else
558  particle%tstop = particle%trelease + this%stoptraveltime
559  if (this%stoptime < particle%tstop) particle%tstop = this%stoptime
560  end if
561 
562  particle%ttrack = particle%trelease
563  particle%itrdomain(level_model) = 0
564  particle%iboundary(level_model) = 0
565  particle%itrdomain(level_feature) = ic
566  particle%iboundary(level_feature) = 0
567  particle%itrdomain(level_subfeature) = 0
568  particle%iboundary(level_subfeature) = 0
569  particle%frctrn = this%frctrn
570  particle%iexmeth = this%iexmeth
571  particle%extend = this%extend
572  particle%icycwin = this%icycwin
573  particle%extol = this%extol
574  end subroutine initialize_particle
575 
576  !> @ brief Read and prepare period data for particle input
577  subroutine prp_rp(this)
578  ! modules
579  use tdismodule, only: kper, nper
582  ! dummy variables
583  class(prtprptype), intent(inout) :: this
584  ! local variables
585  type(characterstringtype), dimension(:), contiguous, &
586  pointer :: settings
587  integer(I4B), pointer :: iper, ionper, nlist
588  character(len=LINELENGTH), allocatable :: lines(:)
589  integer(I4B) :: n
590 
591  ! set pointer to last and next period loaded
592  call mem_setptr(iper, 'IPER', this%input_mempath)
593  call mem_setptr(ionper, 'IONPER', this%input_mempath)
594 
595  if (kper == 1 .and. &
596  (iper == 0) .and. &
597  (ionper > nper) .and. &
598  size(this%schedule%time_select%times) == 0) then
599  ! If the user hasn't provided any release settings (neither
600  ! explicit release times, release time frequency, or period
601  ! block release settings), default to a single release at the
602  ! start of the first period's first time step.
603  allocate (lines(1))
604  lines(1) = "FIRST"
605  call this%schedule%advance(lines=lines)
606  deallocate (lines)
607  return
608  else if (iper /= kper) then
609  return
610  end if
611 
612  ! set input context pointers
613  call mem_setptr(nlist, 'NBOUND', this%input_mempath)
614  call mem_setptr(settings, 'SETTING', this%input_mempath)
615 
616  ! allocate and set input
617  allocate (lines(nlist))
618  do n = 1, nlist
619  lines(n) = settings(n)
620  end do
621 
622  ! update schedule
623  if (size(lines) > 0) &
624  call this%schedule%advance(lines=lines)
625 
626  ! cleanup
627  deallocate (lines)
628  end subroutine prp_rp
629 
630  !> @ brief Calculate flow between package and model.
631  subroutine prp_cq_simrate(this, hnew, flowja, imover)
632  ! modules
633  use tdismodule, only: delt
634  ! dummy variables
635  class(prtprptype) :: this
636  real(DP), dimension(:), intent(in) :: hnew
637  real(DP), dimension(:), intent(inout) :: flowja !< flow between package and model
638  integer(I4B), intent(in) :: imover !< flag indicating if the mover package is active
639  ! local variables
640  integer(I4B) :: i
641  integer(I4B) :: node
642  integer(I4B) :: idiag
643  real(DP) :: rrate
644 
645  ! If no boundaries, skip flow calculations.
646  if (this%nbound <= 0) return
647 
648  ! Loop through each boundary calculating flow.
649  do i = 1, this%nbound
650  node = this%nodelist(i)
651  rrate = dzero
652  ! If cell is no-flow or constant-head, then ignore it.
653  if (node > 0) then
654  ! Calculate the flow rate into the cell.
655  idiag = this%dis%con%ia(node)
656  rrate = this%rptm(i) * (done / delt) ! reciprocal of tstp length
657  flowja(idiag) = flowja(idiag) + rrate
658  end if
659 
660  ! Save simulated value to simvals array.
661  this%simvals(i) = rrate
662  end do
663  end subroutine prp_cq_simrate
664 
665  subroutine define_listlabel(this)
666  class(prtprptype), intent(inout) :: this
667  ! not implemented, not used
668  end subroutine define_listlabel
669 
670  !> @brief Indicates whether observations are supported.
671  logical function prp_obs_supported(this)
672  class(prtprptype) :: this
673  prp_obs_supported = .true.
674  end function prp_obs_supported
675 
676  !> @brief Store supported observations
677  subroutine prp_df_obs(this)
678  ! dummy
679  class(prtprptype) :: this
680  ! local
681  integer(I4B) :: indx
682  call this%obs%StoreObsType('prp', .true., indx)
683  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
684 
685  ! Store obs type and assign procedure pointer
686  ! for to-mvr observation type.
687  call this%obs%StoreObsType('to-mvr', .true., indx)
688  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
689  end subroutine prp_df_obs
690 
691  !> @ brief Set options specific to PrtPrpType
692  subroutine prp_options(this)
693  ! -- modules
696  use openspecmodule, only: access, form
699  ! -- dummy variables
700  class(prtprptype), intent(inout) :: this
701  ! -- local variables
702  character(len=LENVARNAME), dimension(3) :: drytrack_method = &
703  &[character(len=LENVARNAME) :: 'DROP', 'STOP', 'STAY']
704  character(len=lenvarname), dimension(2) :: coorcheck_method = &
705  &[character(len=LENVARNAME) :: 'NONE', 'EAGER']
706  character(len=LINELENGTH) :: trackfile, trackcsvfile, fname
707  type(prtprpparamfoundtype) :: found
708  character(len=*), parameter :: fmtextolwrn = &
709  "('WARNING: EXIT_SOLVE_TOLERANCE is set to ',g10.3,' &
710  &which is much greater than the default value of ',g10.3,'. &
711  &The tolerance that strikes the best balance between accuracy &
712  &and runtime is problem-dependent. Since the variable being &
713  &solved varies from 0 to 1, tolerance values much less than 1 &
714  &typically give the best results.')"
715 
716  ! -- source base class options
717  call this%BndExtType%source_options()
718 
719  ! -- update defaults from input context
720  call mem_set_value(this%stoptime, 'STOPTIME', this%input_mempath, &
721  found%stoptime)
722  call mem_set_value(this%stoptraveltime, 'STOPTRAVELTIME', &
723  this%input_mempath, found%stoptraveltime)
724  call mem_set_value(this%istopweaksink, 'ISTOPWEAKSINK', this%input_mempath, &
725  found%istopweaksink)
726  call mem_set_value(this%istopzone, 'ISTOPZONE', this%input_mempath, &
727  found%istopzone)
728  call mem_set_value(this%drape, 'DRAPE', this%input_mempath, &
729  found%drape)
730  call mem_set_value(this%idrymeth, 'IDRYMETH', this%input_mempath, &
731  drytrack_method, found%idrymeth)
732  call mem_set_value(trackfile, 'TRACKFILE', this%input_mempath, &
733  found%trackfile)
734  call mem_set_value(trackcsvfile, 'TRACKCSVFILE', this%input_mempath, &
735  found%trackcsvfile)
736  call mem_set_value(this%localz, 'LOCALZ', this%input_mempath, &
737  found%localz)
738  call mem_set_value(this%extend, 'EXTEND', this%input_mempath, &
739  found%extend)
740  call mem_set_value(this%extol, 'EXTOL', this%input_mempath, &
741  found%extol)
742  call mem_set_value(this%rttol, 'RTTOL', this%input_mempath, &
743  found%rttol)
744  call mem_set_value(this%rtfreq, 'RTFREQ', this%input_mempath, &
745  found%rtfreq)
746  call mem_set_value(this%frctrn, 'FRCTRN', this%input_mempath, &
747  found%frctrn)
748  call mem_set_value(this%iexmeth, 'IEXMETH', this%input_mempath, &
749  found%iexmeth)
750  call mem_set_value(this%ichkmeth, 'ICHKMETH', this%input_mempath, &
751  coorcheck_method, found%ichkmeth)
752  call mem_set_value(this%icycwin, 'ICYCWIN', this%input_mempath, found%icycwin)
753 
754  ! update internal state and validate input
755  if (found%idrymeth) then
756  if (this%idrymeth == 0) then
757  write (errmsg, '(a)') 'Unsupported dry tracking method. &
758  &DRY_TRACKING_METHOD must be "DROP", "STOP", or "STAY"'
759  call store_error(errmsg)
760  else
761  ! adjust for method zero indexing
762  this%idrymeth = this%idrymeth - 1
763  end if
764  end if
765 
766  if (found%extol) then
767  if (this%extol <= dzero) &
768  call store_error('EXIT_SOLVE_TOLERANCE MUST BE POSITIVE')
769  if (this%extol > dem2) then
770  write (warnmsg, fmt=fmtextolwrn) &
771  this%extol, default_exit_solve_tolerance
772  call store_warning(warnmsg)
773  end if
774  end if
775 
776  if (found%rttol) then
777  if (this%rttol <= dzero) &
778  call store_error('RELEASE_TIME_TOLERANCE MUST BE POSITIVE')
779  end if
780 
781  if (found%rtfreq) then
782  if (this%rtfreq <= dzero) &
783  call store_error('RELEASE_TIME_FREQUENCY MUST BE POSITIVE')
784  end if
785 
786  if (found%iexmeth) then
787  if (.not. (this%iexmeth /= 1 .or. this%iexmeth /= 2)) &
788  call store_error('DEV_EXIT_SOLVE_METHOD MUST BE &
789  &1 (BRENT) OR 2 (CHANDRUPATLA)')
790  end if
791 
792  if (found%ichkmeth) then
793  if (this%ichkmeth == 0) then
794  write (errmsg, '(a)') 'Unsupported coordinate check method. &
795  &COORDINATE_CHECK_METHOD must be "NONE" or "EAGER"'
796  call store_error(errmsg)
797  else
798  ! adjust for method zero based indexing
799  this%ichkmeth = this%ichkmeth - 1
800  end if
801  end if
802 
803  if (found%icycwin) then
804  if (this%icycwin < 0) &
805  call store_error('CYCLE_DETECTION_WINDOW MUST BE NON-NEGATIVE')
806  end if
807 
808  ! fileout options
809  if (found%trackfile) then
810  this%itrkout = getunit()
811  call openfile(this%itrkout, this%iout, trackfile, 'DATA(BINARY)', &
812  form, access, filstat_opt='REPLACE', &
813  mode_opt=mnormal)
814  ! open and write ascii header spec file
815  this%itrkhdr = getunit()
816  fname = trim(trackfile)//'.hdr'
817  call openfile(this%itrkhdr, this%iout, fname, 'CSV', &
818  filstat_opt='REPLACE', mode_opt=mnormal)
819  write (this%itrkhdr, '(a,/,a)') trackheader, trackdtypes
820  end if
821 
822  if (found%trackcsvfile) then
823  this%itrkcsv = getunit()
824  call openfile(this%itrkcsv, this%iout, trackcsvfile, 'CSV', &
825  filstat_opt='REPLACE')
826  write (this%itrkcsv, '(a)') trackheader
827  end if
828 
829  ! terminate if any errors were detected
830  if (count_errors() > 0) then
831  call store_error_filename(this%input_fname)
832  end if
833 
834  ! log found options
835  call this%prp_log_options(found, trackfile, trackcsvfile)
836 
837  ! Create release schedule now that we know
838  ! the coincident release time tolerance
839  this%schedule => create_release_schedule(tolerance=this%rttol)
840  end subroutine prp_options
841 
842  !> @ brief Log options specific to PrtPrpType
843  subroutine prp_log_options(this, found, trackfile, trackcsvfile)
844  ! -- modules
846  ! -- dummy variables
847  class(prtprptype), intent(inout) :: this
848  type(prtprpparamfoundtype), intent(in) :: found
849  character(len=*), intent(in) :: trackfile
850  character(len=*), intent(in) :: trackcsvfile
851  ! -- local variables
852  ! formats
853  character(len=*), parameter :: fmttrkbin = &
854  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO BINARY FILE: ', a, /4x, &
855  &'OPENED ON UNIT: ', I0)"
856  character(len=*), parameter :: fmttrkcsv = &
857  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO CSV FILE: ', a, /4x, &
858  &'OPENED ON UNIT: ', I0)"
859 
860  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
861 
862  if (found%frctrn) then
863  write (this%iout, '(4x,a)') &
864  'IF DISV, TRACKING WILL USE THE TERNARY METHOD REGARDLESS OF CELL TYPE'
865  end if
866 
867  if (found%trackfile) then
868  write (this%iout, fmttrkbin) trim(adjustl(trackfile)), this%itrkout
869  end if
870 
871  if (found%trackcsvfile) then
872  write (this%iout, fmttrkcsv) trim(adjustl(trackcsvfile)), this%itrkcsv
873  end if
874 
875  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
876  end subroutine prp_log_options
877 
878  !> @ brief Set dimensions specific to PrtPrpType
879  subroutine prp_dimensions(this)
880  ! -- modules
883  ! -- dummy variables
884  class(prtprptype), intent(inout) :: this
885  ! -- local variables
886  type(prtprpparamfoundtype) :: found
887 
888  call mem_set_value(this%nreleasepoints, 'NRELEASEPTS', this%input_mempath, &
889  found%nreleasepts)
890  call mem_set_value(this%nreleasetimes, 'NRELEASETIMES', this%input_mempath, &
891  found%nreleasetimes)
892 
893  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
894  write (this%iout, '(4x,a,i0)') 'NRELEASEPTS = ', this%nreleasepoints
895  write (this%iout, '(4x,a,i0)') 'NRELEASETIMES = ', this%nreleasetimes
896  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
897 
898  ! set maxbound and nbound to nreleasepts
899  this%maxbound = this%nreleasepoints
900  this%nbound = this%nreleasepoints
901 
902  ! allocate arrays for prp package
903  call this%prp_allocate_arrays()
904 
905  ! read packagedata and releasetimes blocks
906  call this%prp_packagedata()
907  call this%prp_releasetimes()
908  call this%prp_load_releasetimefrequency()
909  end subroutine prp_dimensions
910 
911  !> @brief Load package data (release points).
912  subroutine prp_packagedata(this)
914  use geomutilmodule, only: get_node
916  ! dummy
917  class(prtprptype), intent(inout) :: this
918  ! local
919  integer(I4B), dimension(:), pointer, contiguous :: irptno
920  integer(I4B), dimension(:, :), pointer, contiguous :: cellids
921  real(DP), dimension(:), pointer, contiguous :: xrpts, yrpts, zrpts
922  type(characterstringtype), dimension(:), pointer, &
923  contiguous :: boundnames
924  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
925  character(len=9) :: cno
926  character(len=20) :: cellidstr
927  integer(I4B), dimension(:), allocatable :: nboundchk
928  integer(I4B), dimension(:), pointer :: cellid
929  integer(I4B) :: n, noder, nodeu, rptno
930 
931  ! set input context pointers
932  call mem_setptr(irptno, 'IRPTNO', this%input_mempath)
933  call mem_setptr(cellids, 'CELLID', this%input_mempath)
934  call mem_setptr(xrpts, 'XRPT', this%input_mempath)
935  call mem_setptr(yrpts, 'YRPT', this%input_mempath)
936  call mem_setptr(zrpts, 'ZRPT', this%input_mempath)
937  call mem_setptr(boundnames, 'BOUNDNAME', this%input_mempath)
938 
939  ! allocate and initialize temporary variables
940  allocate (nboundchk(this%nreleasepoints))
941  do n = 1, this%nreleasepoints
942  nboundchk(n) = 0
943  end do
944 
945  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%packName)) &
946  //' PACKAGEDATA'
947 
948  do n = 1, size(irptno)
949 
950  rptno = irptno(n)
951 
952  if (rptno < 1 .or. rptno > this%nreleasepoints) then
953  write (errmsg, '(a,i0,a,i0,a)') &
954  'Expected ', this%nreleasepoints, ' release points. &
955  &Points must be numbered from 1 to ', this%nreleasepoints, '.'
956  call store_error(errmsg)
957  cycle
958  end if
959 
960  ! increment nboundchk
961  nboundchk(rptno) = nboundchk(rptno) + 1
962 
963  ! set cellid
964  cellid => cellids(:, n)
965 
966  ! set node user
967  if (this%dis%ndim == 1) then
968  nodeu = cellid(1)
969  elseif (this%dis%ndim == 2) then
970  nodeu = get_node(cellid(1), 1, cellid(2), &
971  this%dis%mshape(1), 1, &
972  this%dis%mshape(2))
973  else
974  nodeu = get_node(cellid(1), cellid(2), cellid(3), &
975  this%dis%mshape(1), &
976  this%dis%mshape(2), &
977  this%dis%mshape(3))
978  end if
979 
980  ! set noder
981  noder = this%dis%get_nodenumber(nodeu, 1)
982  if (noder <= 0) then
983  call this%dis%nodeu_to_string(nodeu, cellidstr)
984  write (errmsg, '(a)') &
985  'Particle release point configured for nonexistent cell: '// &
986  trim(adjustl(cellidstr))//'. This cell has IDOMAIN <= 0 and '&
987  &'therefore does not exist in the model grid.'
988  call store_error(errmsg)
989  cycle
990  else
991  this%rptnode(rptno) = noder
992  end if
993 
994  if (this%localz .and. (zrpts(n) < 0 .or. zrpts(n) > 1)) then
995  call store_error('Local z coordinate must fall in the interval [0, 1]')
996  cycle
997  end if
998 
999  ! set coordinates
1000  this%rptx(rptno) = xrpts(n)
1001  this%rpty(rptno) = yrpts(n)
1002  this%rptz(rptno) = zrpts(n)
1003 
1004  ! set default boundname
1005  write (cno, '(i9.9)') rptno
1006  bndname = 'PRP'//cno
1007 
1008  ! read boundnames from file, if provided
1009  if (this%inamedbound /= 0) then
1010  bndnametemp = boundnames(n)
1011  if (bndnametemp /= '') bndname = bndnametemp
1012  else
1013  bndname = ''
1014  end if
1015 
1016  ! set boundname
1017  this%rptname(rptno) = bndname
1018  end do
1019 
1020  write (this%iout, '(1x,a)') &
1021  'END OF '//trim(adjustl(this%packName))//' PACKAGEDATA'
1022 
1023  ! check for duplicate or missing particle release points
1024  do n = 1, this%nreleasepoints
1025  if (nboundchk(n) == 0) then
1026  write (errmsg, '(a,a,1x,i0,a)') 'No data specified for particle ', &
1027  'release point', n, '.'
1028  call store_error(errmsg)
1029  else if (nboundchk(n) > 1) then
1030  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
1031  'Data for particle release point', n, 'specified', nboundchk(n), &
1032  'times.'
1033  call store_error(errmsg)
1034  end if
1035  end do
1036 
1037  ! terminate if any errors were detected
1038  if (count_errors() > 0) then
1039  call store_error_filename(this%input_fname)
1040  end if
1041 
1042  ! cleanup
1043  deallocate (nboundchk)
1044  end subroutine prp_packagedata
1045 
1046  !> @brief Load explicitly specified release times.
1047  subroutine prp_releasetimes(this)
1049  ! dummy
1050  class(prtprptype), intent(inout) :: this
1051  ! local
1052  real(DP), dimension(:), pointer, contiguous :: time
1053  integer(I4B) :: n, isize
1054  real(DP), allocatable :: times(:)
1055 
1056  if (this%nreleasetimes <= 0) return
1057 
1058  ! allocate times array
1059  allocate (times(this%nreleasetimes))
1060 
1061  ! check if input array was read
1062  call get_isize('TIME', this%input_mempath, isize)
1063 
1064  if (isize <= 0) then
1065  errmsg = "RELEASTIMES block expected when &
1066  &NRELEASETIMES dimension is non-zero."
1067  call store_error(errmsg)
1068  call store_error_filename(this%input_fname)
1069  end if
1070 
1071  ! set input context pointer
1072  call mem_setptr(time, 'TIME', this%input_mempath)
1073 
1074  ! set input data
1075  do n = 1, size(time)
1076  times(n) = time(n)
1077  end do
1078 
1079  ! register times with the release schedule
1080  call this%schedule%time_select%extend(times)
1081 
1082  ! make sure times strictly increase
1083  if (.not. this%schedule%time_select%increasing()) then
1084  errmsg = "RELEASTIMES block entries must strictly increase."
1085  call store_error(errmsg)
1086  call store_error_filename(this%input_fname)
1087  end if
1088 
1089  ! deallocate
1090  deallocate (times)
1091  end subroutine prp_releasetimes
1092 
1093  !> @brief Load regularly spaced release times if configured.
1095  ! modules
1096  use tdismodule, only: totalsimtime
1097  ! dummy
1098  class(prtprptype), intent(inout) :: this
1099  ! local
1100  real(DP), allocatable :: times(:)
1101 
1102  ! check if a release time frequency is configured
1103  if (this%rtfreq <= dzero) return
1104 
1105  ! create array of regularly-spaced release times
1106  times = arange( &
1107  start=dzero, &
1108  stop=totalsimtime, &
1109  step=this%rtfreq)
1110 
1111  ! register times with release schedule
1112  call this%schedule%time_select%extend(times)
1113 
1114  ! make sure times strictly increase
1115  if (.not. this%schedule%time_select%increasing()) then
1116  errmsg = "Release times must strictly increase"
1117  call store_error(errmsg)
1118  call store_error_filename(this%input_fname)
1119  end if
1120 
1121  ! deallocate
1122  deallocate (times)
1123 
1124  end subroutine prp_load_releasetimefrequency
1125 
1126 end module prtprpmodule
This module contains the extended boundary package.
This module contains the base boundary package.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dsame
real constant for values that are considered the same based on machine precision
Definition: Constants.f90:122
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
@ tabcenter
centered table column
Definition: Constants.f90:172
@ tableft
left justified table column
Definition: Constants.f90:171
@ mnormal
normal output mode
Definition: Constants.f90:206
real(dp), parameter dep3
real constant 1000
Definition: Constants.f90:88
integer(i4b), parameter lenpakloc
maximum length of a package location
Definition: Constants.f90:50
real(dp), parameter dem1
real constant 1e-1
Definition: Constants.f90:103
real(dp), parameter dep9
real constant 1e9
Definition: Constants.f90:90
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
integer(i4b), parameter lenftype
maximum length of a package type (DIS, WEL, OC, etc.)
Definition: Constants.f90:39
integer(i4b), parameter lenboundname
maximum length of a bound name
Definition: Constants.f90:36
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dem5
real constant 1e-5
Definition: Constants.f90:108
real(dp), parameter dem2
real constant 1e-2
Definition: Constants.f90:105
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Definition: Dis.f90:1
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
integer(i4b) function, public get_node(ilay, irow, icol, nlay, nrow, ncol)
Get node number, given layer, row, and column indices for a structured grid. If any argument is inval...
Definition: GeomUtil.f90:83
logical function, public point_in_polygon(x, y, poly)
Check if a point is within a polygon.
Definition: GeomUtil.f90:27
subroutine, public get_ijk(nodenumber, nrow, ncol, nlay, irow, icol, ilay)
Get row, column and layer indices from node number and grid dimensions. If nodenumber is invalid,...
Definition: GeomUtil.f90:100
subroutine, public get_jk(nodenumber, ncpl, nlay, icpl, ilay)
Get layer index and within-layer node index from node number and grid dimensions. If nodenumber is in...
Definition: GeomUtil.f90:128
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module defines variable data types.
Definition: kind.f90:8
pure real(dp) function, dimension(:), allocatable, public arange(start, stop, step)
Return reals separated by the given step over the given interval.
Definition: MathUtil.f90:384
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Definition: MathUtil.f90:46
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
Particle tracking strategies.
Definition: Method.f90:2
@, public level_feature
Definition: Method.f90:41
@, public level_subfeature
Definition: Method.f90:42
@, public level_model
Definition: Method.f90:40
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public defaultobsidprocessor(obsrv, dis, inunitobs, iout)
@ brief Process IDstring provided for each observation
Definition: Obs.f90:246
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
subroutine create_particle_store(store, np, mempath)
Allocate particle store.
Definition: Particle.f90:152
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:36
subroutine create_particle(particle)
Create a new particle.
Definition: Particle.f90:145
Particle release scheduling.
type(particlereleasescheduletype) function, pointer, public create_release_schedule(tolerance)
Create a new release schedule.
Particle track output module.
character(len= *), parameter, public trackheader
character(len= *), parameter, public trackdtypes
subroutine, public prp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, input_mempath, fmi)
Create a new particle release point package.
Definition: prt-prp.f90:104
subroutine prp_allocate_arrays(this, nodelist, auxvar)
Allocate arrays.
Definition: prt-prp.f90:208
subroutine prp_rp(this)
@ brief Read and prepare period data for particle input
Definition: prt-prp.f90:578
subroutine prp_load_releasetimefrequency(this)
Load regularly spaced release times if configured.
Definition: prt-prp.f90:1095
subroutine prp_cq_simrate(this, hnew, flowja, imover)
@ brief Calculate flow between package and model.
Definition: prt-prp.f90:632
character(len=lenftype) ftype
Definition: prt-prp.f90:34
subroutine prp_df_obs(this)
Store supported observations.
Definition: prt-prp.f90:678
real(dp), parameter default_exit_solve_tolerance
Definition: prt-prp.f90:36
subroutine define_listlabel(this)
Definition: prt-prp.f90:666
subroutine log_release(this)
Log the release scheduled for this time step.
Definition: prt-prp.f90:382
subroutine prp_ad(this)
Advance a time step and release particles if scheduled.
Definition: prt-prp.f90:321
subroutine prp_allocate_scalars(this)
Allocate scalars.
Definition: prt-prp.f90:243
subroutine prp_dimensions(this)
@ brief Set dimensions specific to PrtPrpType
Definition: prt-prp.f90:880
subroutine prp_set_pointers(this, ibound, izone)
@ brief Set pointers to model variables
Definition: prt-prp.f90:198
subroutine initialize_particle(this, particle, ip, trelease)
Definition: prt-prp.f90:458
subroutine prp_da(this)
Deallocate memory.
Definition: prt-prp.f90:151
character(len=16) text
Definition: prt-prp.f90:35
subroutine prp_releasetimes(this)
Load explicitly specified release times.
Definition: prt-prp.f90:1048
subroutine prp_options(this)
@ brief Set options specific to PrtPrpType
Definition: prt-prp.f90:693
subroutine prp_log_options(this, found, trackfile, trackcsvfile)
@ brief Log options specific to PrtPrpType
Definition: prt-prp.f90:844
logical function prp_obs_supported(this)
Indicates whether observations are supported.
Definition: prt-prp.f90:672
subroutine prp_ar(this)
@ brief Allocate and read period data
Definition: prt-prp.f90:302
subroutine release(this, ip, trelease)
Release a particle at the specified time.
Definition: prt-prp.f90:440
subroutine validate_release_point(this, ic, x, y, z)
Verify that the release point is in the cell.
Definition: prt-prp.f90:397
subroutine prp_packagedata(this)
Load package data (release points).
Definition: prt-prp.f90:913
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
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=maxcharlen) errmsg
error message string
character(len=maxcharlen) warnmsg
warning message string
real(dp), pointer, public totalsimtime
time at end of simulation
Definition: tdis.f90:37
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
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21
@ brief BndType
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
Vertex grid discretization.
Definition: Disv.f90:24
Structure of arrays to store particles.
Definition: Particle.f90:104
Particle tracked by the PRT model.
Definition: Particle.f90:57
Manages particle track output (logging/writing).
Particle release point (PRP) package.
Definition: prt-prp.f90:39