MODFLOW 6  version 6.7.0.dev0
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
8  use tablemodule, only: tabletype, table_cr
14  use prtfmimodule, only: prtfmitype
25  use dismodule, only: distype
26  use disvmodule, only: disvtype
27  use errorutilmodule, only: pstop
28  use mathutilmodule, only: arange, is_close
30 
31  implicit none
32 
33  private
34  public :: prtprptype
35  public :: prp_create
36 
37  character(len=LENFTYPE) :: ftype = 'PRP'
38  character(len=16) :: text = ' PRP'
39 
40  !> @brief Particle release point (PRP) package
41  type, extends(bndtype) :: prtprptype
42  type(prtfmitype), pointer :: fmi => null() !< flow model interface
43  type(particlestoretype), pointer :: particles => null() !< particle store
44  type(trackcontroltype), pointer :: trackctl => null() !< track control
45  type(releasescheduletype), pointer :: schedule !< particle release schedule
46  integer(I4B), pointer :: nreleasepoints => null() !< number of release points
47  integer(I4B), pointer :: nreleasetimes => null() !< number of user-specified particle release times
48  integer(I4B), pointer :: nparticles => null() !< number of particles released
49  integer(I4B), pointer :: istopweaksink => null() !< weak sink option: 0 = no stop, 1 = stop
50  integer(I4B), pointer :: istopzone => null() !< optional stop zone number: 0 = no stop zone
51  integer(I4B), pointer :: idrape => null() !< drape option: 0 = do not drape, 1 = drape to topmost active cell
52  integer(I4B), pointer :: idrymeth => null() !< dry tracking method: 0 = drop, 1 = stop, 2 = stay
53  integer(I4B), pointer :: itrkout => null() !< binary track file
54  integer(I4B), pointer :: itrkhdr => null() !< track header file
55  integer(I4B), pointer :: itrkcsv => null() !< CSV track file
56  integer(I4B), pointer :: irlstls => null() !< release time file
57  integer(I4B), pointer :: ilocalz => null() !< compute z coordinates local to the cell
58  integer(I4B), pointer :: iextend => null() !< extend tracking beyond simulation's end
59  integer(I4B), pointer :: ifrctrn => null() !< force ternary solution for quad grids
60  integer(I4B), pointer :: iexmeth => null() !< method for iterative solution of particle exit location and time in generalized Pollock's method
61  real(dp), pointer :: extol => null() !< tolerance for iterative solution of particle exit location and time in generalized Pollock's method
62  real(dp), pointer :: rttol => null() !< tolerance for coincident particle release times
63  real(dp), pointer :: rtfreq => null() !< frequency for regularly spaced release times
64  real(dp), pointer :: offset => null() !< release time offset
65  real(dp), pointer :: stoptime => null() !< stop time for all release points
66  real(dp), pointer :: stoptraveltime => null() !< stop travel time for all points
67  integer(I4B), pointer, contiguous :: rptnode(:) => null() !< release point reduced nns
68  integer(I4B), pointer, contiguous :: rptzone(:) => null() !< release point zone numbers
69  real(dp), pointer, contiguous :: rptx(:) => null() !< release point x coordinates
70  real(dp), pointer, contiguous :: rpty(:) => null() !< release point y coordinates
71  real(dp), pointer, contiguous :: rptz(:) => null() !< release point z coordinates
72  real(dp), pointer, contiguous :: rptm(:) => null() !< total mass released from point
73  character(len=LENBOUNDNAME), pointer, contiguous :: rptname(:) => null() !< release point names
74  contains
75  procedure :: prp_allocate_arrays
76  procedure :: prp_allocate_scalars
77  procedure :: bnd_ar => prp_ar
78  procedure :: bnd_ad => prp_ad
79  procedure :: bnd_rp => prp_rp
80  procedure :: bnd_cq_simrate => prp_cq_simrate
81  procedure :: bnd_da => prp_da
82  procedure :: define_listlabel
83  procedure :: prp_set_pointers
84  procedure :: bnd_options => prp_options
85  procedure :: read_dimensions => prp_read_dimensions
86  procedure :: prp_read_packagedata
87  procedure :: prp_read_releasetimes
89  procedure :: release
90  procedure :: log_release
92  procedure :: initialize_particle
93  procedure, public :: bnd_obs_supported => prp_obs_supported
94  procedure, public :: bnd_df_obs => prp_df_obs
95  end type prtprptype
96 
97 contains
98 
99  !> @brief Create a new particle release point package
100  subroutine prp_create(packobj, id, ibcnum, inunit, iout, namemodel, &
101  pakname, fmi)
102  ! dummy
103  class(bndtype), pointer :: packobj
104  integer(I4B), intent(in) :: id
105  integer(I4B), intent(in) :: ibcnum
106  integer(I4B), intent(in) :: inunit
107  integer(I4B), intent(in) :: iout
108  character(len=*), intent(in) :: namemodel
109  character(len=*), intent(in) :: pakname
110  type(prtfmitype), pointer :: fmi
111  ! local
112  type(prtprptype), pointer :: prpobj
113  ! formats
114  character(len=*), parameter :: fmtheader = &
115  "(1x, /1x, 'PRP PARTICLE RELEASE POINT PACKAGE', &
116  &' INPUT READ FROM UNIT ', i0, /)"
117 
118  ! allocate the object and assign values to object variables
119  allocate (prpobj)
120  packobj => prpobj
121 
122  ! create name and memory path
123  call packobj%set_names(ibcnum, namemodel, pakname, ftype)
124  prpobj%text = text
125 
126  ! allocate scalars
127  call prpobj%prp_allocate_scalars()
128 
129  ! initialize package
130  call packobj%pack_initialize()
131 
132  packobj%inunit = inunit
133  packobj%iout = iout
134  packobj%id = id
135  packobj%ibcnum = ibcnum
136  packobj%ncolbnd = 4
137  packobj%iscloc = 1
138 
139  ! store pointer to flow model interface
140  prpobj%fmi => fmi
141 
142  ! if prp is enabled, print a message identifying it
143  if (inunit > 0) write (iout, fmtheader) inunit
144  end subroutine prp_create
145 
146  !> @brief Deallocate memory
147  subroutine prp_da(this)
148  class(prtprptype) :: this
149 
150  ! Deallocate parent
151  call this%BndType%bnd_da()
152 
153  ! Deallocate scalars
154  call mem_deallocate(this%ilocalz)
155  call mem_deallocate(this%iextend)
156  call mem_deallocate(this%offset)
157  call mem_deallocate(this%stoptime)
158  call mem_deallocate(this%stoptraveltime)
159  call mem_deallocate(this%istopweaksink)
160  call mem_deallocate(this%istopzone)
161  call mem_deallocate(this%idrape)
162  call mem_deallocate(this%idrymeth)
163  call mem_deallocate(this%nreleasepoints)
164  call mem_deallocate(this%nreleasetimes)
165  call mem_deallocate(this%nparticles)
166  call mem_deallocate(this%itrkout)
167  call mem_deallocate(this%itrkhdr)
168  call mem_deallocate(this%itrkcsv)
169  call mem_deallocate(this%irlstls)
170  call mem_deallocate(this%ifrctrn)
171  call mem_deallocate(this%iexmeth)
172  call mem_deallocate(this%extol)
173  call mem_deallocate(this%rttol)
174  call mem_deallocate(this%rtfreq)
175 
176  ! Deallocate arrays
177  call mem_deallocate(this%rptx)
178  call mem_deallocate(this%rpty)
179  call mem_deallocate(this%rptz)
180  call mem_deallocate(this%rptnode)
181  call mem_deallocate(this%rptm)
182  call mem_deallocate(this%rptname, 'RPTNAME', this%memoryPath)
183 
184  ! Deallocate particle store, release time and release step selections
185  call this%particles%deallocate(this%memoryPath)
186  call this%schedule%deallocate()
187  deallocate (this%particles)
188  deallocate (this%schedule)
189  end subroutine prp_da
190 
191  !> @ brief Set pointers to model variables
192  subroutine prp_set_pointers(this, ibound, izone, trackctl)
193  class(prtprptype) :: this
194  integer(I4B), dimension(:), pointer, contiguous :: ibound
195  integer(I4B), dimension(:), pointer, contiguous :: izone
196  type(trackcontroltype), pointer :: trackctl
197 
198  this%ibound => ibound
199  this%rptzone => izone
200  this%trackctl => trackctl
201  end subroutine prp_set_pointers
202 
203  !> @brief Allocate arrays
204  subroutine prp_allocate_arrays(this, nodelist, auxvar)
205  ! dummy
206  class(prtprptype) :: this
207  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
208  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
209  ! local
210  integer(I4B) :: nps
211 
212  ! Allocate particle store, starting with the number
213  ! of release points (arrays resized if/when needed)
215  this%particles, &
216  this%nreleasepoints, &
217  this%memoryPath)
218 
219  ! Allocate arrays
220  call mem_allocate(this%rptx, this%nreleasepoints, 'RPTX', this%memoryPath)
221  call mem_allocate(this%rpty, this%nreleasepoints, 'RPTY', this%memoryPath)
222  call mem_allocate(this%rptz, this%nreleasepoints, 'RPTZ', this%memoryPath)
223  call mem_allocate(this%rptm, this%nreleasepoints, 'RPTMASS', this%memoryPath)
224  call mem_allocate(this%rptnode, this%nreleasepoints, 'RPTNODER', &
225  this%memoryPath)
226  call mem_allocate(this%rptname, lenboundname, this%nreleasepoints, &
227  'RPTNAME', this%memoryPath)
228 
229  ! Initialize arrays
230  do nps = 1, this%nreleasepoints
231  this%rptm(nps) = dzero
232  end do
233  end subroutine prp_allocate_arrays
234 
235  !> @brief Allocate scalars
236  subroutine prp_allocate_scalars(this)
237  class(prtprptype) :: this
238 
239  ! Allocate parent's scalars
240  call this%BndType%allocate_scalars()
241 
242  ! Allocate scalars for this type
243  call mem_allocate(this%ilocalz, 'ILOCALZ', this%memoryPath)
244  call mem_allocate(this%iextend, 'IEXTEND', this%memoryPath)
245  call mem_allocate(this%offset, 'OFFSET', this%memoryPath)
246  call mem_allocate(this%stoptime, 'STOPTIME', this%memoryPath)
247  call mem_allocate(this%stoptraveltime, 'STOPTRAVELTIME', this%memoryPath)
248  call mem_allocate(this%istopweaksink, 'ISTOPWEAKSINK', this%memoryPath)
249  call mem_allocate(this%istopzone, 'ISTOPZONE', this%memoryPath)
250  call mem_allocate(this%idrape, 'IDRAPE', this%memoryPath)
251  call mem_allocate(this%idrymeth, 'IDRYMETH', this%memoryPath)
252  call mem_allocate(this%nreleasepoints, 'NRELEASEPOINTS', this%memoryPath)
253  call mem_allocate(this%nreleasetimes, 'NRELEASETIMES', this%memoryPath)
254  call mem_allocate(this%nparticles, 'NPARTICLES', this%memoryPath)
255  call mem_allocate(this%itrkout, 'ITRKOUT', this%memoryPath)
256  call mem_allocate(this%itrkhdr, 'ITRKHDR', this%memoryPath)
257  call mem_allocate(this%itrkcsv, 'ITRKCSV', this%memoryPath)
258  call mem_allocate(this%irlstls, 'IRLSTLS', this%memoryPath)
259  call mem_allocate(this%ifrctrn, 'IFRCTRN', this%memoryPath)
260  call mem_allocate(this%iexmeth, 'IEXMETH', this%memoryPath)
261  call mem_allocate(this%extol, 'EXTOL', this%memoryPath)
262  call mem_allocate(this%rttol, 'RTTOL', this%memoryPath)
263  call mem_allocate(this%rtfreq, 'RTFREQ', this%memoryPath)
264 
265  ! Set values
266  this%ilocalz = 0
267  this%iextend = 0
268  this%offset = dzero
269  this%stoptime = huge(1d0)
270  this%stoptraveltime = huge(1d0)
271  this%istopweaksink = 0
272  this%istopzone = 0
273  this%idrape = 0
274  this%idrymeth = 0
275  this%nreleasepoints = 0
276  this%nreleasetimes = 0
277  this%nparticles = 0
278  this%itrkout = 0
279  this%itrkhdr = 0
280  this%itrkcsv = 0
281  this%irlstls = 0
282  this%ifrctrn = 0
283  this%iexmeth = 0
284  this%extol = dem5
285  this%rttol = dsame * dep9
286  this%rtfreq = dzero
287 
288  end subroutine prp_allocate_scalars
289 
290  !> @ brief Allocate and read period data
291  subroutine prp_ar(this)
292  ! dummy variables
293  class(prtprptype), intent(inout) :: this
294  ! local variables
295  integer(I4B) :: n
296 
297  call this%obs%obs_ar()
298  call this%BndType%allocate_arrays()
299  if (this%inamedbound /= 0) then
300  do n = 1, this%nreleasepoints
301  this%boundname(n) = this%rptname(n)
302  end do
303  end if
304  do n = 1, this%nreleasepoints
305  this%nodelist(n) = this%rptnode(n)
306  end do
307  end subroutine prp_ar
308 
309  !> @brief Advance a time step and release particles if scheduled.
310  subroutine prp_ad(this)
311  class(prtprptype) :: this
312  integer(I4B) :: ip, it
313 
314  ! Notes
315  ! -----
316  ! Each release point can be thought of as
317  ! a gumball machine with infinite supply:
318  ! a point can release an arbitrary number
319  ! of particles, but only one at any time.
320  ! Coincident release times are merged to
321  ! a single time by the release scheduler.
322 
323  ! Reset mass accumulators for this time step.
324  do ip = 1, this%nreleasepoints
325  this%rptm(ip) = dzero
326  end do
327 
328  ! Advance the release schedule and check if
329  ! any releases will be made this time step.
330  call this%schedule%advance()
331  if (.not. this%schedule%any()) return
332 
333  ! Log the schedule to the list file.
334  call this%log_release()
335 
336  ! Expand the particle store. We know from the
337  ! schedule how many particles will be released.
338  call this%particles%resize( &
339  this%particles%num_stored() + &
340  (this%nreleasepoints * this%schedule%count()), &
341  this%memoryPath)
342 
343  ! Release a particle from each point for
344  ! each release time in the current step.
345  do ip = 1, this%nreleasepoints
346  do it = 1, this%schedule%count()
347  call this%release(ip, this%schedule%times(it))
348  end do
349  end do
350  end subroutine prp_ad
351 
352  !> @brief Log the release scheduled for this time step.
353  subroutine log_release(this)
354  class(prtprptype), intent(inout) :: this !< prp
355  if (this%iprpak > 0) then
356  write (this%iout, "(1x,/1x,a,1x,i0)") &
357  'PARTICLE RELEASE FOR PRP', this%ibcnum
358  call this%schedule%log(this%iout)
359  end if
360  end subroutine log_release
361 
362  !> @brief Verify that the release point is in the cell.
363  !!
364  !! Terminate with an error if the release point lies outside the
365  !! given cell, or if the point is above or below the grid top or
366  !! bottom, respectively.
367  !<
368  subroutine validate_release_point(this, ic, x, y, z)
369  class(prtprptype), intent(inout) :: this !< this instance
370  integer(I4B), intent(in) :: ic !< cell index
371  real(DP), intent(in) :: x, y, z !< release point
372  ! local
373  real(DP), allocatable :: polyverts(:, :)
374 
375  call this%fmi%dis%get_polyverts(ic, polyverts)
376  if (.not. point_in_polygon(x, y, polyverts)) then
377  write (errmsg, '(a,g0,a,g0,a,i0)') &
378  'Error: release point (x=', x, ', y=', y, ') is not in cell ', &
379  this%dis%get_nodeuser(ic)
380  call store_error(errmsg, terminate=.false.)
381  call store_error_unit(this%inunit, terminate=.true.)
382  end if
383  if (z > maxval(this%dis%top)) then
384  write (errmsg, '(a,g0,a,g0,a,i0)') &
385  'Error: release point (z=', z, ') is above grid top ', &
386  maxval(this%dis%top)
387  call store_error(errmsg, terminate=.false.)
388  call store_error_unit(this%inunit, terminate=.true.)
389  else if (z < minval(this%dis%bot)) then
390  write (errmsg, '(a,g0,a,g0,a,i0)') &
391  'Error: release point (z=', z, ') is below grid bottom ', &
392  minval(this%dis%bot)
393  call store_error(errmsg, terminate=.false.)
394  call store_error_unit(this%inunit, terminate=.true.)
395  end if
396  deallocate (polyverts)
397  end subroutine validate_release_point
398 
399  !> Release a particle at the specified time.
400  !!
401  !! Releasing a particle entails validating the particle's
402  !! coordinates and settings, transforming its coordinates
403  !! if needed, initializing the particle's initial tracking
404  !! time to the given release time, storing the particle in
405  !! the particle store (from which the PRT model will later
406  !! retrieve it, apply the tracking method, and check it in
407  !! again), and accumulating the particle's mass (the total
408  !! mass released from each release point is calculated for
409  !! budget reporting).
410  !<
411  subroutine release(this, ip, trelease)
412  ! dummy
413  class(prtprptype), intent(inout) :: this !< this instance
414  integer(I4B), intent(in) :: ip !< particle index
415  real(DP), intent(in) :: trelease !< release time
416  ! local
417  integer(I4B) :: np
418  type(particletype), pointer :: particle
419 
420  call this%initialize_particle(particle, ip, trelease)
421 
422  ! Increment cumulative particle count
423  np = this%nparticles + 1
424  this%nparticles = np
425 
426  ! Save the particle to the store
427  call this%particles%save_particle(particle, np)
428  deallocate (particle)
429 
430  ! Accumulate mass for release point
431  this%rptm(ip) = this%rptm(ip) + done
432 
433  end subroutine release
434 
435  subroutine initialize_particle(this, particle, ip, trelease)
436  class(prtprptype), intent(inout) :: this !< this instance
437  type(particletype), pointer, intent(inout) :: particle !< the particle
438  integer(I4B), intent(in) :: ip !< particle index
439  real(DP), intent(in) :: trelease !< release time
440  ! local
441  integer(I4B) :: irow, icol, ilay, icpl
442  integer(I4B) :: ic, icu, ic_old
443  real(DP) :: x, y, z
444  real(DP) :: top, bot, hds
445 
446  ic = this%rptnode(ip)
447  icu = this%dis%get_nodeuser(ic)
448 
449  call create_particle(particle)
450 
451  if (size(this%boundname) /= 0) then
452  particle%name = this%boundname(ip)
453  else
454  particle%name = ''
455  end if
456 
457  particle%irpt = ip
458  particle%istopweaksink = this%istopweaksink
459  particle%istopzone = this%istopzone
460  particle%idrymeth = this%idrymeth
461  particle%icu = icu
462 
463  select type (dis => this%dis)
464  type is (distype)
465  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
466  type is (disvtype)
467  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
468  end select
469  particle%ilay = ilay
470  particle%izone = this%rptzone(ic)
471 
472  particle%istatus = 0
473  ! If the cell is inactive, either drape the particle
474  ! to the top-most active cell beneath it if drape is
475  ! enabled, or else terminate permanently unreleased.
476  if (this%ibound(ic) == 0) then
477  ic_old = ic
478  if (this%idrape > 0) then
479  call this%dis%highest_active(ic, this%ibound)
480  if (ic == ic_old .or. this%ibound(ic) == 0) &
481  particle%istatus = 8
482  else
483  particle%istatus = 8
484  end if
485  end if
486 
487  ! Load coordinates and transform if needed
488  x = this%rptx(ip)
489  y = this%rpty(ip)
490  if (this%ilocalz > 0) then
491  top = this%fmi%dis%top(ic)
492  bot = this%fmi%dis%bot(ic)
493  hds = this%fmi%gwfhead(ic)
494  z = bot + this%rptz(ip) * (hds - bot)
495  else
496  z = this%rptz(ip)
497  end if
498 
499  call this%validate_release_point(ic, x, y, z)
500 
501  particle%x = x
502  particle%y = y
503  particle%z = z
504  particle%trelease = trelease
505 
506  ! Set stop time to earlier of STOPTIME and STOPTRAVELTIME
507  if (this%stoptraveltime == huge(1d0)) then
508  particle%tstop = this%stoptime
509  else
510  particle%tstop = particle%trelease + this%stoptraveltime
511  if (this%stoptime < particle%tstop) particle%tstop = this%stoptime
512  end if
513 
514  particle%ttrack = particle%trelease
515  particle%idomain(1) = 0
516  particle%iboundary(1) = 0
517  particle%idomain(2) = ic
518  particle%iboundary(2) = 0
519  particle%idomain(3) = 0
520  particle%iboundary(3) = 0
521  particle%ifrctrn = this%ifrctrn
522  particle%iexmeth = this%iexmeth
523  particle%iextend = this%iextend
524  particle%extol = this%extol
525  end subroutine initialize_particle
526 
527  !> @ brief Read and prepare period data for particle input
528  subroutine prp_rp(this)
529  ! modules
530  use tdismodule, only: kper, nper
531  use inputoutputmodule, only: urword
532  ! dummy variables
533  class(prtprptype), intent(inout) :: this
534  ! local variables
535  integer(I4B) :: ierr
536  logical(LGP) :: is_found
537  logical(LGP) :: end_of_block
538  logical(LGP) :: no_blocks
539  character(len=LINELENGTH) :: line
540  character(len=LINELENGTH), allocatable :: lines(:)
541  ! formats
542  character(len=*), parameter :: fmtblkerr = &
543  "('Looking for BEGIN PERIOD iper. &
544  &Found ', a, ' instead.')"
545  character(len=*), parameter :: fmt_steps = &
546  "(6x,'TIME STEP(S) ',50(I0,' '))" ! 50 limit is similar to STEPS in OC
547  character(len=*), parameter :: fmt_freq = &
548  "(6x,'EVERY ',I0,' TIME STEP(S)')"
549  character(len=*), parameter :: fmt_fracs = &
550  "(6x,50(f10.3,' '))"
551 
552  ! Set ionper to the stress period number for which a new block of data
553  ! will be read.
554  if (this%inunit == 0) return
555 
556  ! get stress period data
557  no_blocks = .false.
558  if (this%ionper < kper) then
559  ! get period block
560  call this%parser%GetBlock('PERIOD', is_found, ierr, &
561  supportopenclose=.true., &
562  blockrequired=.false.)
563  if (is_found) then
564  ! read ionper and check for increasing period numbers
565  call this%read_check_ionper()
566  else
567  ! PERIOD block not found
568  if (ierr < 0) then
569  if (kper == 1) then
570  ! End of file found; no period data for the simulation.
571  no_blocks = .true.
572  else
573  ! End of file found; no more period data.
574  this%ionper = nper + 1
575  end if
576  else
577  ! Found invalid block
578  call this%parser%GetCurrentLine(line)
579  write (errmsg, fmtblkerr) adjustl(trim(line))
580  call store_error(errmsg, terminate=.true.)
581  end if
582  end if
583  end if
584 
585  ! If the user hasn't provided any release settings (neither
586  ! explicit release times, release time frequency, or period
587  ! block release settings), default to a single release at the
588  ! start of the first period's first time step.
589  if (no_blocks .and. &
590  kper == 1 .and. &
591  size(this%schedule%time_select%times) == 0) then
592  allocate (lines(1))
593  line = "FIRST"
594  lines(1) = line
595  call this%schedule%advance(lines=lines)
596  else if (this%ionper == kper) then
597  ! If the current stress period matches the
598  ! block we are reading, parse the setting
599  ! and register it with the schedule.
600  allocate (lines(0))
601  recordloop: do
602  call this%parser%GetNextLine(end_of_block)
603  if (end_of_block) exit recordloop
604  call this%parser%GetCurrentLine(line)
605  call expandarray(lines)
606  lines(size(lines)) = line
607  end do recordloop
608  if (size(lines) > 0) &
609  call this%schedule%advance(lines=lines)
610  deallocate (lines)
611  end if
612 
613  end subroutine prp_rp
614 
615  !> @ brief Calculate flow between package and model.
616  subroutine prp_cq_simrate(this, hnew, flowja, imover)
617  ! modules
618  use tdismodule, only: delt
619  ! dummy variables
620  class(prtprptype) :: this
621  real(DP), dimension(:), intent(in) :: hnew
622  real(DP), dimension(:), intent(inout) :: flowja !< flow between package and model
623  integer(I4B), intent(in) :: imover !< flag indicating if the mover package is active
624  ! local variables
625  integer(I4B) :: i
626  integer(I4B) :: node
627  integer(I4B) :: idiag
628  real(DP) :: rrate
629 
630  ! If no boundaries, skip flow calculations.
631  if (this%nbound <= 0) return
632 
633  ! Loop through each boundary calculating flow.
634  do i = 1, this%nbound
635  node = this%nodelist(i)
636  rrate = dzero
637  ! If cell is no-flow or constant-head, then ignore it.
638  if (node > 0) then
639  ! Calculate the flow rate into the cell.
640  idiag = this%dis%con%ia(node)
641  rrate = this%rptm(i) * (done / delt) ! reciprocal of tstp length
642  flowja(idiag) = flowja(idiag) + rrate
643  end if
644 
645  ! Save simulated value to simvals array.
646  this%simvals(i) = rrate
647  end do
648  end subroutine prp_cq_simrate
649 
650  subroutine define_listlabel(this)
651  class(prtprptype), intent(inout) :: this
652  ! not implemented, not used
653  end subroutine define_listlabel
654 
655  !> @brief Indicates whether observations are supported.
656  logical function prp_obs_supported(this)
657  class(prtprptype) :: this
658  prp_obs_supported = .true.
659  end function prp_obs_supported
660 
661  !> @brief Store supported observations
662  subroutine prp_df_obs(this)
663  ! dummy
664  class(prtprptype) :: this
665  ! local
666  integer(I4B) :: indx
667  call this%obs%StoreObsType('prp', .true., indx)
668  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
669 
670  ! Store obs type and assign procedure pointer
671  ! for to-mvr observation type.
672  call this%obs%StoreObsType('to-mvr', .true., indx)
673  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
674  end subroutine prp_df_obs
675 
676  !> @brief Set options specific to PrtPrpType
677  subroutine prp_options(this, option, found)
678  use openspecmodule, only: access, form
679  use constantsmodule, only: maxcharlen, dzero
682  ! dummy
683  class(prtprptype), intent(inout) :: this
684  character(len=*), intent(inout) :: option
685  logical(LGP), intent(inout) :: found
686  ! locals
687  character(len=MAXCHARLEN) :: fname
688  character(len=MAXCHARLEN) :: keyword
689  ! formats
690  character(len=*), parameter :: fmttrkbin = &
691  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO BINARY FILE: ', a, /4x, &
692  &'OPENED ON UNIT: ', I0)"
693  character(len=*), parameter :: fmttrkcsv = &
694  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO CSV FILE: ', a, /4x, &
695  &'OPENED ON UNIT: ', I0)"
696 
697  select case (option)
698  case ('STOPTIME')
699  this%stoptime = this%parser%GetDouble()
700  found = .true.
701  case ('STOPTRAVELTIME')
702  this%stoptraveltime = this%parser%GetDouble()
703  found = .true.
704  case ('STOP_AT_WEAK_SINK')
705  this%istopweaksink = 1
706  found = .true.
707  case ('ISTOPZONE')
708  this%istopzone = this%parser%GetInteger()
709  found = .true.
710  case ('DRAPE')
711  this%idrape = 1
712  found = .true.
713  case ('DRY_TRACKING_METHOD')
714  call this%parser%GetStringCaps(keyword)
715  select case (keyword)
716  case ('DROP')
717  this%idrymeth = 0
718  case ('STOP')
719  this%idrymeth = 1
720  case ('STAY')
721  this%idrymeth = 2
722  case default
723  write (errmsg, '(a, a)') &
724  'Unknown dry tracking method: ', trim(keyword)
725  call store_error(errmsg)
726  write (errmsg, '(a, a)') &
727  'DRY must be "DROP", "STOP" or "STAY"'
728  call store_error(errmsg)
729  call this%parser%StoreErrorUnit()
730  end select
731  found = .true.
732  case ('TRACK')
733  call this%parser%GetStringCaps(keyword)
734  if (keyword == 'FILEOUT') then
735  ! parse filename
736  call this%parser%GetString(fname)
737  ! open binary output file
738  this%itrkout = getunit()
739  call openfile(this%itrkout, this%iout, fname, 'DATA(BINARY)', &
740  form, access, filstat_opt='REPLACE', &
741  mode_opt=mnormal)
742  write (this%iout, fmttrkbin) trim(adjustl(fname)), this%itrkout
743  ! open and write ascii header spec file
744  this%itrkhdr = getunit()
745  fname = trim(fname)//'.hdr'
746  call openfile(this%itrkhdr, this%iout, fname, 'CSV', &
747  filstat_opt='REPLACE', mode_opt=mnormal)
748  write (this%itrkhdr, '(a,/,a)') trackheader, trackdtypes
749  else
750  call store_error('OPTIONAL TRACK KEYWORD MUST BE '// &
751  'FOLLOWED BY FILEOUT')
752  end if
753  found = .true.
754  case ('TRACKCSV')
755  call this%parser%GetStringCaps(keyword)
756  if (keyword == 'FILEOUT') then
757  ! parse filename
758  call this%parser%GetString(fname)
759  ! open CSV output file and write headers
760  this%itrkcsv = getunit()
761  call openfile(this%itrkcsv, this%iout, fname, 'CSV', &
762  filstat_opt='REPLACE')
763  write (this%iout, fmttrkcsv) trim(adjustl(fname)), this%itrkcsv
764  write (this%itrkcsv, '(a)') trackheader
765  else
766  call store_error('OPTIONAL TRACKCSV KEYWORD MUST BE &
767  &FOLLOWED BY FILEOUT')
768  end if
769  found = .true.
770  case ('LOCAL_Z')
771  this%ilocalz = 1
772  found = .true.
773  case ('EXTEND_TRACKING')
774  this%iextend = 1
775  found = .true.
776  case ('EXIT_SOLVE_TOLERANCE')
777  this%extol = this%parser%GetDouble()
778  if (this%extol <= dzero) &
779  call store_error('EXIT_SOLVE_TOLERANCE MUST BE POSITIVE')
780  found = .true.
781  case ('RELEASE_TIME_TOLERANCE')
782  this%rttol = this%parser%GetDouble()
783  if (this%rttol <= dzero) &
784  call store_error('RELEASE_TIME_TOLERANCE MUST BE POSITIVE')
785  found = .true.
786  case ('RELEASE_TIME_FREQUENCY')
787  this%rtfreq = this%parser%GetDouble()
788  if (this%rtfreq <= dzero) &
789  call store_error('RELEASE_TIME_FREQUENCY MUST BE POSITIVE')
790  found = .true.
791  case ('DEV_FORCETERNARY')
792  call this%parser%DevOpt()
793  this%ifrctrn = 1
794  write (this%iout, '(4x,a)') &
795  'TRACKING WILL BE DONE USING THE TERNARY METHOD REGARDLESS OF CELL TYPE'
796  found = .true.
797  case ('DEV_EXIT_SOLVE_METHOD')
798  call this%parser%DevOpt()
799  this%iexmeth = this%parser%GetInteger()
800  if (.not. (this%iexmeth /= 1 .or. this%iexmeth /= 2)) &
801  call store_error('DEV_EXIT_SOLVE_METHOD MUST BE &
802  &1 (BRENT) OR 2 (CHANDRUPATLA)')
803  found = .true.
804  case default
805  found = .false.
806  end select
807 
808  ! Catch unrecognized options
809  if (.not. found) then
810  errmsg = "UNKNOWN PRP OPTION '"//trim(keyword)//"'."
811  call store_error(errmsg)
812  call this%parser%StoreErrorUnit()
813  end if
814 
815  ! Create release schedule now that we know
816  ! the coincident release time tolerance
817  this%schedule => create_release_schedule(tol=this%rttol)
818 
819  end subroutine prp_options
820 
821  !> @brief Read package dimensions
822  subroutine prp_read_dimensions(this)
823  ! dummy
824  class(prtprptype), intent(inout) :: this
825  ! local
826  character(len=LINELENGTH) :: errmsg, keyword
827  integer(I4B) :: ierr
828  logical :: isfound, endOfBlock
829 
830  ! get dimension block
831  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
832  supportopenclose=.true.)
833 
834  ! parse dimension block if detected
835  if (isfound) then
836  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
837  do
838  call this%parser%GetNextLine(endofblock)
839  if (endofblock) exit
840  call this%parser%GetStringCaps(keyword)
841  select case (keyword)
842  case ('NRELEASEPTS')
843  this%nreleasepoints = this%parser%GetInteger()
844  case ('NRELEASETIMES')
845  this%nreleasetimes = this%parser%GetInteger()
846  case default
847  write (errmsg, &
848  '(4x,a,a)') '****ERROR. UNKNOWN PARTICLE INPUT DIMENSION: ', &
849  trim(keyword)
850  call store_error(errmsg)
851  call this%parser%StoreErrorUnit()
852  end select
853  end do
854  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
855  else
856  call store_error('ERROR. REQUIRED DIMENSIONS BLOCK NOT FOUND.')
857  end if
858 
859  ! set maxbound and nbound to nreleasepts
860  this%maxbound = this%nreleasepoints
861  this%nbound = this%nreleasepoints
862 
863  ! allocate arrays for prp package
864  call this%prp_allocate_arrays()
865 
866  ! read packagedata and releasetimes blocks
867  call this%prp_read_packagedata()
868  call this%prp_read_releasetimes()
869  call this%prp_load_releasetimefrequency()
870  end subroutine prp_read_dimensions
871 
872  !> @brief Load package data (release points).
873  subroutine prp_read_packagedata(this)
874  ! dummy
875  class(prtprptype), intent(inout) :: this
876  ! local
877  character(len=LINELENGTH) :: cellid
878  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
879  character(len=9) :: cno
880  logical :: isfound
881  logical :: endOfBlock
882  integer(I4B) :: ival
883  integer(I4B) :: n
884  integer(I4B) :: ierr
885  character(len=LENBOUNDNAME), dimension(:), allocatable :: nametxt
886  integer(I4B), dimension(:), allocatable :: nboundchk
887  integer(I4B), dimension(:), allocatable :: noder
888  real(DP), dimension(:), allocatable :: x
889  real(DP), dimension(:), allocatable :: y
890  real(DP), dimension(:), allocatable :: z
891  real(DP), dimension(:), allocatable :: tstop
892  ! format
893  character(len=*), parameter :: fmttend = &
894  "('end time (', G0, ') must be greater than or equal to the &
895  &begin time (', G0, ').')"
896 
897  ! allocate temporary variables
898  allocate (noder(this%nreleasepoints))
899  allocate (x(this%nreleasepoints))
900  allocate (y(this%nreleasepoints))
901  allocate (z(this%nreleasepoints))
902  allocate (tstop(this%nreleasepoints))
903  allocate (nametxt(this%nreleasepoints))
904  allocate (nboundchk(this%nreleasepoints))
905 
906  ! initialize temporary variables
907  do n = 1, this%nreleasepoints
908  nboundchk(n) = 0
909  end do
910 
911  ! read particle release point data
912  ! get particle release points block
913  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
914  supportopenclose=.true.)
915 
916  ! parse block if detected
917  if (isfound) then
918  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%packName)) &
919  //' PACKAGEDATA'
920  do
921  call this%parser%GetNextLine(endofblock)
922  if (endofblock) exit
923  ival = this%parser%GetInteger()
924  n = ival
925 
926  if (n < 1 .or. n > this%nreleasepoints) then
927  write (errmsg, '(a,1x,i0,a)') &
928  'Release point number must be greater than 0 and less than', &
929  'or equal to', this%nreleasepoints, '.'
930  call store_error(errmsg)
931  cycle
932  end if
933 
934  ! increment nboundchk
935  nboundchk(n) = nboundchk(n) + 1
936 
937  ! node number
938  call this%parser%GetCellid(this%dis%ndim, cellid)
939  noder(n) = this%dis%noder_from_cellid(cellid, this%inunit, this%iout)
940 
941  ! x, y, z coordinates
942  x(n) = this%parser%GetDouble()
943  y(n) = this%parser%GetDouble()
944  z(n) = this%parser%GetDouble()
945 
946  if (this%ilocalz > 0 .and. (z(n) < 0 .or. z(n) > 1)) then
947  call store_error('Local z coordinate must fall in the interval [0, 1]')
948  cycle
949  end if
950 
951  ! set default boundname
952  write (cno, '(i9.9)') n
953  bndname = 'PRP'//cno
954 
955  ! read boundnames from file, if provided
956  if (this%inamedbound /= 0) then
957  call this%parser%GetStringCaps(bndnametemp)
958  if (bndnametemp /= '') &
959  bndname = bndnametemp
960  else
961  bndname = ''
962  end if
963 
964  ! store temp boundnames
965  nametxt(n) = bndname
966  end do
967 
968  write (this%iout, '(1x,a)') &
969  'END OF '//trim(adjustl(this%packName))//' PACKAGEDATA'
970 
971  ! check for duplicate or missing particle release points
972  do n = 1, this%nreleasepoints
973  if (nboundchk(n) == 0) then
974  write (errmsg, '(a,a,1x,i0,a)') 'No data specified for particle ', &
975  'release point', n, '.'
976  call store_error(errmsg)
977  else if (nboundchk(n) > 1) then
978  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
979  'Data for particle release point', n, 'specified', nboundchk(n), &
980  'times.'
981  call store_error(errmsg)
982  end if
983  end do
984  else
985  call store_error('Required packagedata block not found.')
986  end if
987 
988  ! terminate if any errors were detected
989  if (count_errors() > 0) then
990  call this%parser%StoreErrorUnit()
991  end if
992 
993  ! fill particle release point data with data stored in temporary local arrays
994  do n = 1, this%nreleasepoints
995  this%rptnode(n) = noder(n)
996  this%rptx(n) = x(n)
997  this%rpty(n) = y(n)
998  this%rptz(n) = z(n)
999  this%rptname(n) = nametxt(n)
1000  end do
1001 
1002  ! deallocate local storage
1003  deallocate (noder)
1004  deallocate (x)
1005  deallocate (y)
1006  deallocate (z)
1007  deallocate (tstop)
1008  deallocate (nametxt)
1009  deallocate (nboundchk)
1010  end subroutine prp_read_packagedata
1011 
1012  !> @brief Load explicitly specified release times.
1013  subroutine prp_read_releasetimes(this)
1014  ! dummy
1015  class(prtprptype), intent(inout) :: this
1016  ! local
1017  integer(I4B) :: i, ierr
1018  logical(LGP) :: eob, found, success
1019  real(DP) :: t
1020  real(DP), allocatable :: times(:)
1021 
1022  ! get releasetimes block
1023  call this%parser%GetBlock('RELEASETIMES', found, ierr, &
1024  supportopenclose=.true., &
1025  blockrequired=.false.)
1026 
1027  ! raise an error if releasetimes has a dimension
1028  ! but no block was found, otherwise return early
1029  if (.not. found) then
1030  if (this%nreleasetimes <= 0) return
1031  write (errmsg, '(a, i0)') &
1032  "Expected RELEASETIMES with length ", this%nreleasetimes
1033  call store_error(errmsg)
1034  call this%parser%StoreErrorUnit(terminate=.true.)
1035  end if
1036 
1037  ! allocate times array
1038  allocate (times(this%nreleasetimes))
1039 
1040  ! read times from the block
1041  write (this%iout, '(/1x,a)') &
1042  'PROCESSING '//trim(adjustl(this%text))//' RELEASETIMES'
1043  do i = 1, this%nreleasetimes
1044  call this%parser%GetNextLine(eob)
1045  if (eob) exit
1046  call this%parser%TryGetDouble(t, success)
1047  if (.not. success) then
1048  errmsg = "Failed to read double precision value"
1049  call store_error(errmsg)
1050  call this%parser%StoreErrorUnit(terminate=.true.)
1051  end if
1052  times(i) = t
1053  end do
1054 
1055  ! register times with the release schedule
1056  call this%schedule%time_select%extend(times)
1057 
1058  ! make sure times strictly increase
1059  if (.not. this%schedule%time_select%increasing()) then
1060  errmsg = "Release times must strictly increase"
1061  call store_error(errmsg)
1062  call this%parser%StoreErrorUnit(terminate=.true.)
1063  end if
1064 
1065  ! deallocate
1066  deallocate (times)
1067 
1068  end subroutine prp_read_releasetimes
1069 
1070  !> @brief Load regularly spaced release times if configured.
1072  ! modules
1073  use tdismodule, only: totalsimtime
1074  ! dummy
1075  class(prtprptype), intent(inout) :: this
1076  ! local
1077  real(DP), allocatable :: times(:)
1078 
1079  ! check if a release time frequency is configured
1080  if (this%rtfreq <= dzero) return
1081 
1082  ! create array of regularly-spaced release times
1083  times = arange( &
1084  start=dzero, &
1085  stop=totalsimtime, &
1086  step=this%rtfreq)
1087 
1088  ! register times with release schedule
1089  call this%schedule%time_select%extend(times)
1090 
1091  ! make sure times strictly increase
1092  if (.not. this%schedule%time_select%increasing()) then
1093  errmsg = "Release times must strictly increase"
1094  call store_error(errmsg)
1095  call this%parser%StoreErrorUnit(terminate=.true.)
1096  end if
1097 
1098  ! deallocate
1099  deallocate (times)
1100 
1101  end subroutine prp_load_releasetimefrequency
1102 
1103 end module prtprpmodule
This module contains block parser methods.
Definition: BlockParser.f90:7
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 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
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
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
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
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
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
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, public allocate_particle_store(this, np, mempath)
Create a new particle store.
Definition: Particle.f90:127
subroutine, public create_particle(particle)
Create a new particle.
Definition: Particle.f90:119
subroutine prp_set_pointers(this, ibound, izone, trackctl)
@ brief Set pointers to model variables
Definition: prt-prp.f90:193
subroutine prp_allocate_arrays(this, nodelist, auxvar)
Allocate arrays.
Definition: prt-prp.f90:205
subroutine prp_rp(this)
@ brief Read and prepare period data for particle input
Definition: prt-prp.f90:529
subroutine prp_load_releasetimefrequency(this)
Load regularly spaced release times if configured.
Definition: prt-prp.f90:1072
subroutine prp_cq_simrate(this, hnew, flowja, imover)
@ brief Calculate flow between package and model.
Definition: prt-prp.f90:617
subroutine prp_read_dimensions(this)
Read package dimensions.
Definition: prt-prp.f90:823
character(len=lenftype) ftype
Definition: prt-prp.f90:37
subroutine prp_read_releasetimes(this)
Load explicitly specified release times.
Definition: prt-prp.f90:1014
subroutine prp_df_obs(this)
Store supported observations.
Definition: prt-prp.f90:663
subroutine, public prp_create(packobj, id, ibcnum, inunit, iout, namemodel, pakname, fmi)
Create a new particle release point package.
Definition: prt-prp.f90:102
subroutine define_listlabel(this)
Definition: prt-prp.f90:651
subroutine log_release(this)
Log the release scheduled for this time step.
Definition: prt-prp.f90:354
subroutine prp_ad(this)
Advance a time step and release particles if scheduled.
Definition: prt-prp.f90:311
subroutine prp_allocate_scalars(this)
Allocate scalars.
Definition: prt-prp.f90:237
subroutine initialize_particle(this, particle, ip, trelease)
Definition: prt-prp.f90:436
subroutine prp_da(this)
Deallocate memory.
Definition: prt-prp.f90:148
character(len=16) text
Definition: prt-prp.f90:38
logical function prp_obs_supported(this)
Indicates whether observations are supported.
Definition: prt-prp.f90:657
subroutine prp_ar(this)
@ brief Allocate and read period data
Definition: prt-prp.f90:292
subroutine release(this, ip, trelease)
Release a particle at the specified time.
Definition: prt-prp.f90:412
subroutine prp_read_packagedata(this)
Load package data (release points).
Definition: prt-prp.f90:874
subroutine validate_release_point(this, ic, x, y, z)
Verify that the release point is in the cell.
Definition: prt-prp.f90:369
subroutine prp_options(this, option, found)
Set options specific to PrtPrpType.
Definition: prt-prp.f90:678
Particle release scheduling.
type(releasescheduletype) function, pointer, public create_release_schedule(tol)
Create a new release schedule object.
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_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
character(len=maxcharlen) warnmsg
warning message string
subroutine, public table_cr(this, name, title)
Definition: Table.f90:87
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
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
character(len= *), parameter, public trackdtypes
Definition: TrackFile.f90:81
character(len= *), parameter, public trackheader
Definition: TrackFile.f90:77
@ brief BndType
Structured grid discretization.
Definition: Dis.f90:23
Vertex grid discretization.
Definition: Disv.f90:24
Structure of arrays to store particles.
Definition: Particle.f90:78
Particle tracked by the PRT model.
Definition: Particle.f90:32
Particle release point (PRP) package.
Definition: prt-prp.f90:41
Particle release scheduling utility.
Manages particle track (i.e. pathline) files.