MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
prtprpmodule Module Reference

Data Types

type  prtprptype
 Particle release point (PRP) package. More...
 

Functions/Subroutines

subroutine, public prp_create (packobj, id, ibcnum, inunit, iout, namemodel, pakname, input_mempath, fmi)
 Create a new particle release point package. More...
 
subroutine prp_da (this)
 Deallocate memory. More...
 
subroutine prp_set_pointers (this, ibound, izone)
 @ brief Set pointers to model variables More...
 
subroutine prp_allocate_arrays (this, nodelist, auxvar)
 Allocate arrays. More...
 
subroutine prp_allocate_scalars (this)
 Allocate scalars. More...
 
subroutine prp_ar (this)
 @ brief Allocate and read period data More...
 
subroutine prp_ad (this)
 Advance a time step and release particles if scheduled. More...
 
subroutine log_release (this)
 Log the release scheduled for this time step. More...
 
subroutine validate_release_point (this, ic, x, y, z)
 Verify that the release point is in the cell. More...
 
subroutine release (this, ip, trelease)
 Release a particle at the specified time. More...
 
subroutine initialize_particle (this, particle, ip, trelease)
 
subroutine prp_rp (this)
 @ brief Read and prepare period data for particle input More...
 
subroutine prp_cq_simrate (this, hnew, flowja, imover)
 @ brief Calculate flow between package and model. More...
 
subroutine define_listlabel (this)
 
logical function prp_obs_supported (this)
 Indicates whether observations are supported. More...
 
subroutine prp_df_obs (this)
 Store supported observations. More...
 
subroutine prp_options (this)
 @ brief Set options specific to PrtPrpType More...
 
subroutine prp_log_options (this, found, trackfile, trackcsvfile)
 @ brief Log options specific to PrtPrpType More...
 
subroutine prp_dimensions (this)
 @ brief Set dimensions specific to PrtPrpType More...
 
subroutine prp_packagedata (this)
 Load package data (release points). More...
 
subroutine prp_releasetimes (this)
 Load explicitly specified release times. More...
 
subroutine prp_load_releasetimefrequency (this)
 Load regularly spaced release times if configured. More...
 

Variables

character(len=lenftype) ftype = 'PRP'
 
character(len=16) text = ' PRP'
 
real(dp), parameter default_exit_solve_tolerance = DEM5
 

Function/Subroutine Documentation

◆ define_listlabel()

subroutine prtprpmodule::define_listlabel ( class(prtprptype), intent(inout)  this)

Definition at line 679 of file prt-prp.f90.

680  class(PrtPrpType), intent(inout) :: this
681  ! not implemented, not used

◆ initialize_particle()

subroutine prtprpmodule::initialize_particle ( class(prtprptype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b), intent(in)  ip,
real(dp), intent(in)  trelease 
)
private
Parameters
[in,out]thisthis instance
[in,out]particlethe particle
[in]ipparticle index
[in]treleaserelease time

Definition at line 475 of file prt-prp.f90.

477  class(PrtPrpType), intent(inout) :: this !< this instance
478  type(ParticleType), pointer, intent(inout) :: particle !< the particle
479  integer(I4B), intent(in) :: ip !< particle index
480  real(DP), intent(in) :: trelease !< release time
481  ! local
482  integer(I4B) :: irow, icol, ilay, icpl
483  integer(I4B) :: ic, icu, ic_old
484  real(DP) :: x, y, z
485  real(DP) :: top, bot, hds
486  ! formats
487  character(len=*), parameter :: fmticterr = &
488  "('Error in ',a,': Flow model interface does not contain ICELLTYPE. &
489  &ICELLTYPE is required for PRT to distinguish convertible cells &
490  &from confined cells if LOCAL_Z release coordinates are provided. &
491  &Make sure a GWFGRID entry is configured in the PRT FMI package.')"
492 
493  ic = this%rptnode(ip)
494  icu = this%dis%get_nodeuser(ic)
495 
496  call create_particle(particle)
497 
498  if (size(this%boundname) /= 0) then
499  particle%name = this%boundname(ip)
500  else
501  particle%name = ''
502  end if
503 
504  particle%irpt = ip
505  particle%istopweaksink = this%istopweaksink
506  particle%istopzone = this%istopzone
507  particle%idrymeth = this%idrymeth
508  particle%icu = icu
509 
510  select type (dis => this%dis)
511  type is (distype)
512  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
513  type is (disvtype)
514  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
515  end select
516  particle%ilay = ilay
517  particle%izone = this%rptzone(ic)
518  particle%istatus = 0 ! status 0 until tracking starts
519  ! If the cell is inactive, either drape the particle
520  ! to the top-most active cell beneath it if drape is
521  ! enabled, or else terminate permanently unreleased.
522  if (this%ibound(ic) == 0) then
523  ic_old = ic
524  if (this%drape) then
525  call this%dis%highest_active(ic, this%ibound)
526  if (ic == ic_old .or. this%ibound(ic) == 0) then
527  ! negative unreleased status signals to the
528  ! tracking method that we haven't yet saved
529  ! a termination record, it needs to do so.
530  particle%istatus = -1 * term_unreleased
531  end if
532  else
533  particle%istatus = -1 * term_unreleased
534  end if
535  end if
536 
537  ! load coordinates
538  x = this%rptx(ip)
539  y = this%rpty(ip)
540  if (this%localz) then
541  ! make sure FMI has cell type array. we need
542  ! it to distinguish convertible and confined
543  ! cells if release z coordinates are local
544  if (this%fmi%igwfceltyp /= 1) then
545  write (errmsg, fmticterr) trim(this%text)
546  call store_error(errmsg, terminate=.true.)
547  end if
548 
549  ! calculate model z coord from local z coord.
550  ! if cell is confined (icelltype == 0) use the
551  ! actual cell height (geometric top - bottom).
552  ! otherwise use head as cell top, clamping to
553  ! the cell bottom if head is below the bottom
554  ! and to geometric cell top if head is above top
555  top = this%fmi%dis%top(ic)
556  bot = this%fmi%dis%bot(ic)
557  if (this%fmi%gwfceltyp(icu) /= 0) then
558  hds = this%fmi%gwfhead(ic)
559  top = min(top, hds)
560  top = max(top, bot)
561  end if
562  z = bot + this%rptz(ip) * (top - bot)
563  else
564  z = this%rptz(ip)
565  end if
566 
567  if (this%ichkmeth > 0) &
568  call this%validate_release_point(ic, x, y, z)
569 
570  particle%x = x
571  particle%y = y
572  particle%z = z
573  particle%trelease = trelease
574 
575  ! Set stop time to earlier of STOPTIME and STOPTRAVELTIME
576  if (this%stoptraveltime == huge(1d0)) then
577  particle%tstop = this%stoptime
578  else
579  particle%tstop = particle%trelease + this%stoptraveltime
580  if (this%stoptime < particle%tstop) particle%tstop = this%stoptime
581  end if
582 
583  particle%ttrack = particle%trelease
584  particle%itrdomain(level_model) = 0
585  particle%iboundary(level_model) = 0
586  particle%itrdomain(level_feature) = ic
587  particle%iboundary(level_feature) = 0
588  particle%itrdomain(level_subfeature) = 0
589  particle%iboundary(level_subfeature) = 0
590  particle%frctrn = this%frctrn
591  particle%iexmeth = this%iexmeth
592  particle%extend = this%extend
593  particle%icycwin = this%icycwin
594  particle%extol = this%extol
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:36
Here is the call graph for this function:

◆ log_release()

subroutine prtprpmodule::log_release ( class(prtprptype), intent(inout)  this)
Parameters
[in,out]thisprp

Definition at line 399 of file prt-prp.f90.

400  class(PrtPrpType), intent(inout) :: this !< prp
401  if (this%iprpak > 0) then
402  write (this%iout, "(1x,/1x,a,1x,i0)") &
403  'PARTICLE RELEASE FOR PRP', this%ibcnum
404  call this%schedule%log(this%iout)
405  end if

◆ prp_ad()

subroutine prtprpmodule::prp_ad ( class(prtprptype this)
private

Definition at line 326 of file prt-prp.f90.

327  use tdismodule, only: totalsimtime, kstp, kper
328  class(PrtPrpType) :: this
329  integer(I4B) :: ip, it
330  real(DP) :: t
331 
332  ! Notes
333  ! -----
334  ! Each release point can be thought of as
335  ! a gumball machine with infinite supply:
336  ! a point can release an arbitrary number
337  ! of particles, but only one at any time.
338  ! Coincident release times are merged to
339  ! a single time by the release scheduler.
340 
341  ! Reset mass accumulators for this time step.
342  do ip = 1, this%nreleasepoints
343  this%rptm(ip) = dzero
344  end do
345 
346  ! Advance the release schedule. At the start of each period (kstp==1),
347  ! apply period block configuration if available and not yet applied.
348  ! This handles both new configuration and fill-forward periods.
349  ! For subsequent time steps, just advance without arguments to
350  ! advance the time selection object to the current time step.
351  if (kstp == 1 .and. &
352  kper /= this%applied_kper .and. &
353  allocated(this%period_block_lines)) then
354  call this%schedule%advance(lines=this%period_block_lines)
355  this%applied_kper = kper
356  else
357  call this%schedule%advance()
358  end if
359 
360  ! Check if any releases will be made this time step.
361  if (.not. this%schedule%any()) return
362 
363  ! Log the schedule to the list file.
364  call this%log_release()
365 
366  ! Expand the particle store. We know from the
367  ! schedule how many particles will be released.
368  call this%particles%resize( &
369  this%particles%num_stored() + &
370  (this%nreleasepoints * this%schedule%count()), &
371  this%memoryPath)
372 
373  ! Release a particle from each point for
374  ! each release time in the current step.
375  do ip = 1, this%nreleasepoints
376  do it = 1, this%schedule%count()
377  t = this%schedule%times(it)
378  ! Skip the release time if it's before the simulation
379  ! starts, or if no `extend_tracking`, after it ends.
380  if (t < dzero) then
381  write (warnmsg, '(a,g0,a)') &
382  'Skipping negative release time (t=', t, ').'
383  call store_warning(warnmsg)
384  cycle
385  else if (t > totalsimtime .and. .not. this%extend) then
386  write (warnmsg, '(a,g0,a)') &
387  'Skipping release time falling after the end of the &
388  &simulation (t=', t, '). Enable EXTEND_TRACKING to &
389  &release particles after the simulation end time.'
390  call store_warning(warnmsg)
391  cycle
392  end if
393  call this%release(ip, t)
394  end do
395  end do
real(dp), pointer, public totalsimtime
time at end of simulation
Definition: tdis.f90:37
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ prp_allocate_arrays()

subroutine prtprpmodule::prp_allocate_arrays ( class(prtprptype this,
integer(i4b), dimension(:), optional, pointer, contiguous  nodelist,
real(dp), dimension(:, :), optional, pointer, contiguous  auxvar 
)
private

Definition at line 212 of file prt-prp.f90.

213  ! dummy
214  class(PrtPrpType) :: this
215  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
216  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
217  ! local
218  integer(I4B) :: nps
219 
220  call this%BndExtType%allocate_arrays()
221 
222  ! Allocate particle store, starting with the number
223  ! of release points (arrays resized if/when needed)
224  call create_particle_store( &
225  this%particles, &
226  this%nreleasepoints, &
227  this%memoryPath)
228 
229  ! Allocate arrays
230  call mem_allocate(this%rptx, this%nreleasepoints, 'RPTX', this%memoryPath)
231  call mem_allocate(this%rpty, this%nreleasepoints, 'RPTY', this%memoryPath)
232  call mem_allocate(this%rptz, this%nreleasepoints, 'RPTZ', this%memoryPath)
233  call mem_allocate(this%rptm, this%nreleasepoints, 'RPTMASS', &
234  this%memoryPath)
235  call mem_allocate(this%rptnode, this%nreleasepoints, 'RPTNODER', &
236  this%memoryPath)
237  call mem_allocate(this%rptname, lenboundname, this%nreleasepoints, &
238  'RPTNAME', this%memoryPath)
239 
240  ! Initialize arrays
241  do nps = 1, this%nreleasepoints
242  this%rptm(nps) = dzero
243  end do
Here is the call graph for this function:

◆ prp_allocate_scalars()

subroutine prtprpmodule::prp_allocate_scalars ( class(prtprptype this)
private

Definition at line 247 of file prt-prp.f90.

248  class(PrtPrpType) :: this
249 
250  ! Allocate parent's scalars
251  call this%BndExtType%allocate_scalars()
252 
253  ! Allocate scalars for this type
254  call mem_allocate(this%localz, 'LOCALZ', this%memoryPath)
255  call mem_allocate(this%extend, 'EXTEND', this%memoryPath)
256  call mem_allocate(this%offset, 'OFFSET', this%memoryPath)
257  call mem_allocate(this%stoptime, 'STOPTIME', this%memoryPath)
258  call mem_allocate(this%stoptraveltime, 'STOPTRAVELTIME', this%memoryPath)
259  call mem_allocate(this%istopweaksink, 'ISTOPWEAKSINK', this%memoryPath)
260  call mem_allocate(this%istopzone, 'ISTOPZONE', this%memoryPath)
261  call mem_allocate(this%drape, 'DRAPE', this%memoryPath)
262  call mem_allocate(this%idrymeth, 'IDRYMETH', this%memoryPath)
263  call mem_allocate(this%nreleasepoints, 'NRELEASEPOINTS', this%memoryPath)
264  call mem_allocate(this%nreleasetimes, 'NRELEASETIMES', this%memoryPath)
265  call mem_allocate(this%nparticles, 'NPARTICLES', this%memoryPath)
266  call mem_allocate(this%itrkout, 'ITRKOUT', this%memoryPath)
267  call mem_allocate(this%itrkhdr, 'ITRKHDR', this%memoryPath)
268  call mem_allocate(this%itrkcsv, 'ITRKCSV', this%memoryPath)
269  call mem_allocate(this%irlstls, 'IRLSTLS', this%memoryPath)
270  call mem_allocate(this%frctrn, 'FRCTRN', this%memoryPath)
271  call mem_allocate(this%iexmeth, 'IEXMETH', this%memoryPath)
272  call mem_allocate(this%ichkmeth, 'ICHKMETH', this%memoryPath)
273  call mem_allocate(this%icycwin, 'ICYCWIN', this%memoryPath)
274  call mem_allocate(this%extol, 'EXTOL', this%memoryPath)
275  call mem_allocate(this%rttol, 'RTTOL', this%memoryPath)
276  call mem_allocate(this%rtfreq, 'RTFREQ', this%memoryPath)
277 
278  ! Set values
279  this%localz = .false.
280  this%extend = .false.
281  this%offset = dzero
282  this%stoptime = huge(1d0)
283  this%stoptraveltime = huge(1d0)
284  this%istopweaksink = 0
285  this%istopzone = 0
286  this%drape = .false.
287  this%idrymeth = 0
288  this%nreleasepoints = 0
289  this%nreleasetimes = 0
290  this%nparticles = 0
291  this%itrkout = 0
292  this%itrkhdr = 0
293  this%itrkcsv = 0
294  this%irlstls = 0
295  this%frctrn = .false.
296  this%iexmeth = 0
297  this%ichkmeth = 1
298  this%icycwin = 0
299  this%extol = default_exit_solve_tolerance
300  this%rttol = dsame * dep9
301  this%rtfreq = dzero
302  this%applied_kper = 0
303 

◆ prp_ar()

subroutine prtprpmodule::prp_ar ( class(prtprptype), intent(inout)  this)
private

Definition at line 307 of file prt-prp.f90.

308  ! dummy variables
309  class(PrtPrpType), intent(inout) :: this
310  ! local variables
311  integer(I4B) :: n
312 
313  call this%obs%obs_ar()
314 
315  if (this%inamedbound /= 0) then
316  do n = 1, this%nreleasepoints
317  this%boundname(n) = this%rptname(n)
318  end do
319  end if
320  do n = 1, this%nreleasepoints
321  this%nodelist(n) = this%rptnode(n)
322  end do

◆ prp_cq_simrate()

subroutine prtprpmodule::prp_cq_simrate ( class(prtprptype this,
real(dp), dimension(:), intent(in)  hnew,
real(dp), dimension(:), intent(inout)  flowja,
integer(i4b), intent(in)  imover 
)
Parameters
[in,out]flowjaflow between package and model
[in]imoverflag indicating if the mover package is active

Definition at line 645 of file prt-prp.f90.

646  ! modules
647  use tdismodule, only: delt
648  ! dummy variables
649  class(PrtPrpType) :: this
650  real(DP), dimension(:), intent(in) :: hnew
651  real(DP), dimension(:), intent(inout) :: flowja !< flow between package and model
652  integer(I4B), intent(in) :: imover !< flag indicating if the mover package is active
653  ! local variables
654  integer(I4B) :: i
655  integer(I4B) :: node
656  integer(I4B) :: idiag
657  real(DP) :: rrate
658 
659  ! If no boundaries, skip flow calculations.
660  if (this%nbound <= 0) return
661 
662  ! Loop through each boundary calculating flow.
663  do i = 1, this%nbound
664  node = this%nodelist(i)
665  rrate = dzero
666  ! If cell is no-flow or constant-head, then ignore it.
667  if (node > 0) then
668  ! Calculate the flow rate into the cell.
669  idiag = this%dis%con%ia(node)
670  rrate = this%rptm(i) * (done / delt) ! reciprocal of tstp length
671  flowja(idiag) = flowja(idiag) + rrate
672  end if
673 
674  ! Save simulated value to simvals array.
675  this%simvals(i) = rrate
676  end do
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29

◆ prp_create()

subroutine, public prtprpmodule::prp_create ( class(bndtype), pointer  packobj,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  ibcnum,
integer(i4b), intent(in)  inunit,
integer(i4b), intent(in)  iout,
character(len=*), intent(in)  namemodel,
character(len=*), intent(in)  pakname,
character(len=*), intent(in)  input_mempath,
type(prtfmitype), pointer  fmi 
)

Definition at line 104 of file prt-prp.f90.

106  ! dummy
107  class(BndType), pointer :: packobj
108  integer(I4B), intent(in) :: id
109  integer(I4B), intent(in) :: ibcnum
110  integer(I4B), intent(in) :: inunit
111  integer(I4B), intent(in) :: iout
112  character(len=*), intent(in) :: namemodel
113  character(len=*), intent(in) :: pakname
114  character(len=*), intent(in) :: input_mempath
115  type(PrtFmiType), pointer :: fmi
116  ! local
117  type(PrtPrpType), pointer :: prpobj
118  ! formats
119  character(len=*), parameter :: fmtheader = &
120  "(1x, /1x, 'PRP PARTICLE RELEASE POINT PACKAGE', &
121  &' INPUT READ FROM MEMPATH: ', a, /)"
122 
123  ! allocate the object and assign values to object variables
124  allocate (prpobj)
125  packobj => prpobj
126 
127  ! create name and memory path
128  call packobj%set_names(ibcnum, namemodel, pakname, ftype, input_mempath)
129  prpobj%text = text
130 
131  ! allocate scalars
132  call prpobj%prp_allocate_scalars()
133 
134  ! initialize package
135  call packobj%pack_initialize()
136 
137  packobj%inunit = inunit
138  packobj%iout = iout
139  packobj%id = id
140  packobj%ibcnum = ibcnum
141  packobj%ncolbnd = 4
142  packobj%iscloc = 1
143 
144  ! store pointer to flow model interface
145  prpobj%fmi => fmi
146 
147  ! if prp is enabled, print a message identifying it
148  if (inunit > 0) write (iout, fmtheader) input_mempath
Here is the caller graph for this function:

◆ prp_da()

subroutine prtprpmodule::prp_da ( class(prtprptype this)
private

Definition at line 152 of file prt-prp.f90.

153  class(PrtPrpType) :: this
154 
155  ! Deallocate parent
156  call this%BndExtType%bnd_da()
157 
158  ! Deallocate scalars
159  call mem_deallocate(this%localz)
160  call mem_deallocate(this%extend)
161  call mem_deallocate(this%offset)
162  call mem_deallocate(this%stoptime)
163  call mem_deallocate(this%stoptraveltime)
164  call mem_deallocate(this%istopweaksink)
165  call mem_deallocate(this%istopzone)
166  call mem_deallocate(this%drape)
167  call mem_deallocate(this%idrymeth)
168  call mem_deallocate(this%nreleasepoints)
169  call mem_deallocate(this%nreleasetimes)
170  call mem_deallocate(this%nparticles)
171  call mem_deallocate(this%itrkout)
172  call mem_deallocate(this%itrkhdr)
173  call mem_deallocate(this%itrkcsv)
174  call mem_deallocate(this%irlstls)
175  call mem_deallocate(this%frctrn)
176  call mem_deallocate(this%iexmeth)
177  call mem_deallocate(this%ichkmeth)
178  call mem_deallocate(this%icycwin)
179  call mem_deallocate(this%extol)
180  call mem_deallocate(this%rttol)
181  call mem_deallocate(this%rtfreq)
182 
183  ! Deallocate arrays
184  call mem_deallocate(this%rptx)
185  call mem_deallocate(this%rpty)
186  call mem_deallocate(this%rptz)
187  call mem_deallocate(this%rptnode)
188  call mem_deallocate(this%rptm)
189  call mem_deallocate(this%rptname, 'RPTNAME', this%memoryPath)
190 
191  ! Deallocate period block storage
192  if (allocated(this%period_block_lines)) deallocate (this%period_block_lines)
193 
194  ! Deallocate objects
195  call this%particles%destroy(this%memoryPath)
196  call this%schedule%destroy()
197  deallocate (this%particles)
198  deallocate (this%schedule)

◆ prp_df_obs()

subroutine prtprpmodule::prp_df_obs ( class(prtprptype this)
private

Definition at line 691 of file prt-prp.f90.

692  ! dummy
693  class(PrtPrpType) :: this
694  ! local
695  integer(I4B) :: indx
696  call this%obs%StoreObsType('prp', .true., indx)
697  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
698 
699  ! Store obs type and assign procedure pointer
700  ! for to-mvr observation type.
701  call this%obs%StoreObsType('to-mvr', .true., indx)
702  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ prp_dimensions()

subroutine prtprpmodule::prp_dimensions ( class(prtprptype), intent(inout)  this)

Definition at line 893 of file prt-prp.f90.

894  ! -- modules
897  ! -- dummy variables
898  class(PrtPrpType), intent(inout) :: this
899  ! -- local variables
900  type(PrtPrpParamFoundType) :: found
901 
902  call mem_set_value(this%nreleasepoints, 'NRELEASEPTS', this%input_mempath, &
903  found%nreleasepts)
904  call mem_set_value(this%nreleasetimes, 'NRELEASETIMES', this%input_mempath, &
905  found%nreleasetimes)
906 
907  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
908  write (this%iout, '(4x,a,i0)') 'NRELEASEPTS = ', this%nreleasepoints
909  write (this%iout, '(4x,a,i0)') 'NRELEASETIMES = ', this%nreleasetimes
910  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
911 
912  ! set maxbound and nbound to nreleasepts
913  this%maxbound = this%nreleasepoints
914  this%nbound = this%nreleasepoints
915 
916  ! allocate arrays for prp package
917  call this%prp_allocate_arrays()
918 
919  ! read packagedata and releasetimes blocks
920  call this%prp_packagedata()
921  call this%prp_releasetimes()
922  call this%prp_load_releasetimefrequency()

◆ prp_load_releasetimefrequency()

subroutine prtprpmodule::prp_load_releasetimefrequency ( class(prtprptype), intent(inout)  this)

Definition at line 1108 of file prt-prp.f90.

1109  ! modules
1110  use tdismodule, only: totalsimtime
1111  ! dummy
1112  class(PrtPrpType), intent(inout) :: this
1113  ! local
1114  real(DP), allocatable :: times(:)
1115 
1116  ! check if a release time frequency is configured
1117  if (this%rtfreq <= dzero) return
1118 
1119  ! create array of regularly-spaced release times
1120  times = arange( &
1121  start=dzero, &
1122  stop=totalsimtime, &
1123  step=this%rtfreq)
1124 
1125  ! register times with release schedule
1126  call this%schedule%time_select%extend(times)
1127 
1128  ! make sure times strictly increase
1129  if (.not. this%schedule%time_select%increasing()) then
1130  errmsg = "Release times must strictly increase"
1131  call store_error(errmsg)
1132  call store_error_filename(this%input_fname)
1133  end if
1134 
1135  ! deallocate
1136  deallocate (times)
1137 
Here is the call graph for this function:

◆ prp_log_options()

subroutine prtprpmodule::prp_log_options ( class(prtprptype), intent(inout)  this,
type(prtprpparamfoundtype), intent(in)  found,
character(len=*), intent(in)  trackfile,
character(len=*), intent(in)  trackcsvfile 
)

Definition at line 857 of file prt-prp.f90.

858  ! -- modules
860  ! -- dummy variables
861  class(PrtPrpType), intent(inout) :: this
862  type(PrtPrpParamFoundType), intent(in) :: found
863  character(len=*), intent(in) :: trackfile
864  character(len=*), intent(in) :: trackcsvfile
865  ! -- local variables
866  ! formats
867  character(len=*), parameter :: fmttrkbin = &
868  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO BINARY FILE: ', a, /4x, &
869  &'OPENED ON UNIT: ', I0)"
870  character(len=*), parameter :: fmttrkcsv = &
871  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO CSV FILE: ', a, /4x, &
872  &'OPENED ON UNIT: ', I0)"
873 
874  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
875 
876  if (found%frctrn) then
877  write (this%iout, '(4x,a)') &
878  'IF DISV, TRACKING WILL USE THE TERNARY METHOD REGARDLESS OF CELL TYPE'
879  end if
880 
881  if (found%trackfile) then
882  write (this%iout, fmttrkbin) trim(adjustl(trackfile)), this%itrkout
883  end if
884 
885  if (found%trackcsvfile) then
886  write (this%iout, fmttrkcsv) trim(adjustl(trackcsvfile)), this%itrkcsv
887  end if
888 
889  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'

◆ prp_obs_supported()

logical function prtprpmodule::prp_obs_supported ( class(prtprptype this)
private

Definition at line 685 of file prt-prp.f90.

686  class(PrtPrpType) :: this
687  prp_obs_supported = .true.

◆ prp_options()

subroutine prtprpmodule::prp_options ( class(prtprptype), intent(inout)  this)
private

Definition at line 706 of file prt-prp.f90.

707  ! -- modules
710  use openspecmodule, only: access, form
713  ! -- dummy variables
714  class(PrtPrpType), intent(inout) :: this
715  ! -- local variables
716  character(len=LENVARNAME), dimension(3) :: drytrack_method = &
717  &[character(len=LENVARNAME) :: 'DROP', 'STOP', 'STAY']
718  character(len=lenvarname), dimension(2) :: coorcheck_method = &
719  &[character(len=LENVARNAME) :: 'NONE', 'EAGER']
720  character(len=LINELENGTH) :: trackfile, trackcsvfile, fname
721  type(PrtPrpParamFoundType) :: found
722  character(len=*), parameter :: fmtextolwrn = &
723  "('WARNING: EXIT_SOLVE_TOLERANCE is set to ',g10.3,' &
724  &which is much greater than the default value of ',g10.3,'. &
725  &The tolerance that strikes the best balance between accuracy &
726  &and runtime is problem-dependent. Since the variable being &
727  &solved varies from 0 to 1, tolerance values much less than 1 &
728  &typically give the best results.')"
729 
730  ! -- source base class options
731  call this%BndExtType%source_options()
732 
733  ! -- update defaults from input context
734  call mem_set_value(this%stoptime, 'STOPTIME', this%input_mempath, &
735  found%stoptime)
736  call mem_set_value(this%stoptraveltime, 'STOPTRAVELTIME', &
737  this%input_mempath, found%stoptraveltime)
738  call mem_set_value(this%istopweaksink, 'ISTOPWEAKSINK', this%input_mempath, &
739  found%istopweaksink)
740  call mem_set_value(this%istopzone, 'ISTOPZONE', this%input_mempath, &
741  found%istopzone)
742  call mem_set_value(this%drape, 'DRAPE', this%input_mempath, &
743  found%drape)
744  call mem_set_value(this%idrymeth, 'IDRYMETH', this%input_mempath, &
745  drytrack_method, found%idrymeth)
746  call mem_set_value(trackfile, 'TRACKFILE', this%input_mempath, &
747  found%trackfile)
748  call mem_set_value(trackcsvfile, 'TRACKCSVFILE', this%input_mempath, &
749  found%trackcsvfile)
750  call mem_set_value(this%localz, 'LOCALZ', this%input_mempath, &
751  found%localz)
752  call mem_set_value(this%extend, 'EXTEND', this%input_mempath, &
753  found%extend)
754  call mem_set_value(this%extol, 'EXTOL', this%input_mempath, &
755  found%extol)
756  call mem_set_value(this%rttol, 'RTTOL', this%input_mempath, &
757  found%rttol)
758  call mem_set_value(this%rtfreq, 'RTFREQ', this%input_mempath, &
759  found%rtfreq)
760  call mem_set_value(this%frctrn, 'FRCTRN', this%input_mempath, &
761  found%frctrn)
762  call mem_set_value(this%iexmeth, 'IEXMETH', this%input_mempath, &
763  found%iexmeth)
764  call mem_set_value(this%ichkmeth, 'ICHKMETH', this%input_mempath, &
765  coorcheck_method, found%ichkmeth)
766  call mem_set_value(this%icycwin, 'ICYCWIN', this%input_mempath, found%icycwin)
767 
768  ! update internal state and validate input
769  if (found%idrymeth) then
770  if (this%idrymeth == 0) then
771  write (errmsg, '(a)') 'Unsupported dry tracking method. &
772  &DRY_TRACKING_METHOD must be "DROP", "STOP", or "STAY"'
773  call store_error(errmsg)
774  else
775  ! adjust for method zero indexing
776  this%idrymeth = this%idrymeth - 1
777  end if
778  end if
779 
780  if (found%extol) then
781  if (this%extol <= dzero) &
782  call store_error('EXIT_SOLVE_TOLERANCE MUST BE POSITIVE')
783  if (this%extol > dem2) then
784  write (warnmsg, fmt=fmtextolwrn) &
785  this%extol, default_exit_solve_tolerance
786  call store_warning(warnmsg)
787  end if
788  end if
789 
790  if (found%rttol) then
791  if (this%rttol <= dzero) &
792  call store_error('RELEASE_TIME_TOLERANCE MUST BE POSITIVE')
793  end if
794 
795  if (found%rtfreq) then
796  if (this%rtfreq <= dzero) &
797  call store_error('RELEASE_TIME_FREQUENCY MUST BE POSITIVE')
798  end if
799 
800  if (found%iexmeth) then
801  if (.not. (this%iexmeth /= 1 .or. this%iexmeth /= 2)) &
802  call store_error('DEV_EXIT_SOLVE_METHOD MUST BE &
803  &1 (BRENT) OR 2 (CHANDRUPATLA)')
804  end if
805 
806  if (found%ichkmeth) then
807  if (this%ichkmeth == 0) then
808  write (errmsg, '(a)') 'Unsupported coordinate check method. &
809  &COORDINATE_CHECK_METHOD must be "NONE" or "EAGER"'
810  call store_error(errmsg)
811  else
812  ! adjust for method zero based indexing
813  this%ichkmeth = this%ichkmeth - 1
814  end if
815  end if
816 
817  if (found%icycwin) then
818  if (this%icycwin < 0) &
819  call store_error('CYCLE_DETECTION_WINDOW MUST BE NON-NEGATIVE')
820  end if
821 
822  ! fileout options
823  if (found%trackfile) then
824  this%itrkout = getunit()
825  call openfile(this%itrkout, this%iout, trackfile, 'DATA(BINARY)', &
826  form, access, filstat_opt='REPLACE', &
827  mode_opt=mnormal)
828  ! open and write ascii header spec file
829  this%itrkhdr = getunit()
830  fname = trim(trackfile)//'.hdr'
831  call openfile(this%itrkhdr, this%iout, fname, 'CSV', &
832  filstat_opt='REPLACE', mode_opt=mnormal)
833  write (this%itrkhdr, '(a,/,a)') trackheader, trackdtypes
834  end if
835 
836  if (found%trackcsvfile) then
837  this%itrkcsv = getunit()
838  call openfile(this%itrkcsv, this%iout, trackcsvfile, 'CSV', &
839  filstat_opt='REPLACE')
840  write (this%itrkcsv, '(a)') trackheader
841  end if
842 
843  ! terminate if any errors were detected
844  if (count_errors() > 0) then
845  call store_error_filename(this%input_fname)
846  end if
847 
848  ! log found options
849  call this%prp_log_options(found, trackfile, trackcsvfile)
850 
851  ! Create release schedule now that we know
852  ! the coincident release time tolerance
853  this%schedule => create_release_schedule(tolerance=this%rttol)
This module contains simulation constants.
Definition: Constants.f90:9
@ mnormal
normal output mode
Definition: Constants.f90:206
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
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
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
Here is the call graph for this function:

◆ prp_packagedata()

subroutine prtprpmodule::prp_packagedata ( class(prtprptype), intent(inout)  this)

Definition at line 926 of file prt-prp.f90.

928  use geomutilmodule, only: get_node
930  ! dummy
931  class(PrtPrpType), intent(inout) :: this
932  ! local
933  integer(I4B), dimension(:), pointer, contiguous :: irptno
934  integer(I4B), dimension(:, :), pointer, contiguous :: cellids
935  real(DP), dimension(:), pointer, contiguous :: xrpts, yrpts, zrpts
936  type(CharacterStringType), dimension(:), pointer, &
937  contiguous :: boundnames
938  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
939  character(len=9) :: cno
940  character(len=20) :: cellidstr
941  integer(I4B), dimension(:), allocatable :: nboundchk
942  integer(I4B), dimension(:), pointer :: cellid
943  integer(I4B) :: n, noder, nodeu, rptno
944 
945  ! set input context pointers
946  call mem_setptr(irptno, 'IRPTNO', this%input_mempath)
947  call mem_setptr(cellids, 'CELLID', this%input_mempath)
948  call mem_setptr(xrpts, 'XRPT', this%input_mempath)
949  call mem_setptr(yrpts, 'YRPT', this%input_mempath)
950  call mem_setptr(zrpts, 'ZRPT', this%input_mempath)
951  call mem_setptr(boundnames, 'BOUNDNAME', this%input_mempath)
952 
953  ! allocate and initialize temporary variables
954  allocate (nboundchk(this%nreleasepoints))
955  do n = 1, this%nreleasepoints
956  nboundchk(n) = 0
957  end do
958 
959  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%packName)) &
960  //' PACKAGEDATA'
961 
962  do n = 1, size(irptno)
963 
964  rptno = irptno(n)
965 
966  if (rptno < 1 .or. rptno > this%nreleasepoints) then
967  write (errmsg, '(a,i0,a,i0,a)') &
968  'Expected ', this%nreleasepoints, ' release points. &
969  &Points must be numbered from 1 to ', this%nreleasepoints, '.'
970  call store_error(errmsg)
971  cycle
972  end if
973 
974  ! increment nboundchk
975  nboundchk(rptno) = nboundchk(rptno) + 1
976 
977  ! set cellid
978  cellid => cellids(:, n)
979 
980  ! set node user
981  if (this%dis%ndim == 1) then
982  nodeu = cellid(1)
983  elseif (this%dis%ndim == 2) then
984  nodeu = get_node(cellid(1), 1, cellid(2), &
985  this%dis%mshape(1), 1, &
986  this%dis%mshape(2))
987  else
988  nodeu = get_node(cellid(1), cellid(2), cellid(3), &
989  this%dis%mshape(1), &
990  this%dis%mshape(2), &
991  this%dis%mshape(3))
992  end if
993 
994  ! set noder
995  noder = this%dis%get_nodenumber(nodeu, 1)
996  if (noder <= 0) then
997  call this%dis%nodeu_to_string(nodeu, cellidstr)
998  write (errmsg, '(a)') &
999  'Particle release point configured for nonexistent cell: '// &
1000  trim(adjustl(cellidstr))//'. This cell has IDOMAIN <= 0 and '&
1001  &'therefore does not exist in the model grid.'
1002  call store_error(errmsg)
1003  cycle
1004  else
1005  this%rptnode(rptno) = noder
1006  end if
1007 
1008  if (this%localz .and. (zrpts(n) < 0 .or. zrpts(n) > 1)) then
1009  call store_error('Local z coordinate must fall in the interval [0, 1]')
1010  cycle
1011  end if
1012 
1013  ! set coordinates
1014  this%rptx(rptno) = xrpts(n)
1015  this%rpty(rptno) = yrpts(n)
1016  this%rptz(rptno) = zrpts(n)
1017 
1018  ! set default boundname
1019  write (cno, '(i9.9)') rptno
1020  bndname = 'PRP'//cno
1021 
1022  ! read boundnames from file, if provided
1023  if (this%inamedbound /= 0) then
1024  bndnametemp = boundnames(n)
1025  if (bndnametemp /= '') bndname = bndnametemp
1026  else
1027  bndname = ''
1028  end if
1029 
1030  ! set boundname
1031  this%rptname(rptno) = bndname
1032  end do
1033 
1034  write (this%iout, '(1x,a)') &
1035  'END OF '//trim(adjustl(this%packName))//' PACKAGEDATA'
1036 
1037  ! check for duplicate or missing particle release points
1038  do n = 1, this%nreleasepoints
1039  if (nboundchk(n) == 0) then
1040  write (errmsg, '(a,a,1x,i0,a)') 'No data specified for particle ', &
1041  'release point', n, '.'
1042  call store_error(errmsg)
1043  else if (nboundchk(n) > 1) then
1044  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
1045  'Data for particle release point', n, 'specified', nboundchk(n), &
1046  'times.'
1047  call store_error(errmsg)
1048  end if
1049  end do
1050 
1051  ! terminate if any errors were detected
1052  if (count_errors() > 0) then
1053  call store_error_filename(this%input_fname)
1054  end if
1055 
1056  ! cleanup
1057  deallocate (nboundchk)
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
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ prp_releasetimes()

subroutine prtprpmodule::prp_releasetimes ( class(prtprptype), intent(inout)  this)

Definition at line 1061 of file prt-prp.f90.

1063  ! dummy
1064  class(PrtPrpType), intent(inout) :: this
1065  ! local
1066  real(DP), dimension(:), pointer, contiguous :: time
1067  integer(I4B) :: n, isize
1068  real(DP), allocatable :: times(:)
1069 
1070  if (this%nreleasetimes <= 0) return
1071 
1072  ! allocate times array
1073  allocate (times(this%nreleasetimes))
1074 
1075  ! check if input array was read
1076  call get_isize('TIME', this%input_mempath, isize)
1077 
1078  if (isize <= 0) then
1079  errmsg = "RELEASTIMES block expected when &
1080  &NRELEASETIMES dimension is non-zero."
1081  call store_error(errmsg)
1082  call store_error_filename(this%input_fname)
1083  end if
1084 
1085  ! set input context pointer
1086  call mem_setptr(time, 'TIME', this%input_mempath)
1087 
1088  ! set input data
1089  do n = 1, size(time)
1090  times(n) = time(n)
1091  end do
1092 
1093  ! register times with the release schedule
1094  call this%schedule%time_select%extend(times)
1095 
1096  ! make sure times strictly increase
1097  if (.not. this%schedule%time_select%increasing()) then
1098  errmsg = "RELEASTIMES block entries must strictly increase."
1099  call store_error(errmsg)
1100  call store_error_filename(this%input_fname)
1101  end if
1102 
1103  ! deallocate
1104  deallocate (times)
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
Here is the call graph for this function:

◆ prp_rp()

subroutine prtprpmodule::prp_rp ( class(prtprptype), intent(inout)  this)

Definition at line 598 of file prt-prp.f90.

599  ! modules
600  use tdismodule, only: kper, nper
603  ! dummy variables
604  class(PrtPrpType), intent(inout) :: this
605  ! local variables
606  type(CharacterStringType), dimension(:), contiguous, &
607  pointer :: settings
608  integer(I4B), pointer :: iper, ionper, nlist
609  integer(I4B) :: n
610 
611  ! set pointer to last and next period loaded
612  call mem_setptr(iper, 'IPER', this%input_mempath)
613  call mem_setptr(ionper, 'IONPER', this%input_mempath)
614 
615  if (kper == 1 .and. &
616  (iper == 0) .and. &
617  (ionper > nper) .and. &
618  size(this%schedule%time_select%times) == 0) then
619  ! If the user hasn't provided any release settings (neither
620  ! explicit release times, release time frequency, or period
621  ! block release settings), default to a single release at the
622  ! start of the first period's first time step.
623  ! Store default configuration; advance() will be called in prp_ad().
624  if (allocated(this%period_block_lines)) deallocate (this%period_block_lines)
625  allocate (this%period_block_lines(1))
626  this%period_block_lines(1) = "FIRST"
627  return
628  else if (iper /= kper) then
629  return
630  end if
631 
632  ! set input context pointers
633  call mem_setptr(nlist, 'NBOUND', this%input_mempath)
634  call mem_setptr(settings, 'SETTING', this%input_mempath)
635 
636  ! Store period block configuration for fill-forward.
637  if (allocated(this%period_block_lines)) deallocate (this%period_block_lines)
638  allocate (this%period_block_lines(nlist))
639  do n = 1, nlist
640  this%period_block_lines(n) = settings(n)
641  end do
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21

◆ prp_set_pointers()

subroutine prtprpmodule::prp_set_pointers ( class(prtprptype this,
integer(i4b), dimension(:), pointer, contiguous  ibound,
integer(i4b), dimension(:), pointer, contiguous  izone 
)
private

Definition at line 202 of file prt-prp.f90.

203  class(PrtPrpType) :: this
204  integer(I4B), dimension(:), pointer, contiguous :: ibound
205  integer(I4B), dimension(:), pointer, contiguous :: izone
206 
207  this%ibound => ibound
208  this%rptzone => izone

◆ release()

subroutine prtprpmodule::release ( class(prtprptype), intent(inout)  this,
integer(i4b), intent(in)  ip,
real(dp), intent(in)  trelease 
)
private

Releasing a particle entails validating the particle's coordinates and settings, transforming its coordinates if needed, initializing the particle's initial tracking time to the given release time, storing the particle in the particle store (from which the PRT model will later retrieve it, apply the tracking method, and check it in again), and accumulating the particle's mass (the total mass released from each release point is calculated for budget reporting).

Parameters
[in,out]thisthis instance
[in]ipparticle index
[in]treleaserelease time

Definition at line 457 of file prt-prp.f90.

458  ! dummy
459  class(PrtPrpType), intent(inout) :: this !< this instance
460  integer(I4B), intent(in) :: ip !< particle index
461  real(DP), intent(in) :: trelease !< release time
462  ! local
463  integer(I4B) :: np
464  type(ParticleType), pointer :: particle
465 
466  call this%initialize_particle(particle, ip, trelease)
467  np = this%nparticles + 1
468  this%nparticles = np
469  call this%particles%put(particle, np)
470  deallocate (particle)
471  this%rptm(ip) = this%rptm(ip) + done ! TODO configurable mass
472 

◆ validate_release_point()

subroutine prtprpmodule::validate_release_point ( class(prtprptype), intent(inout)  this,
integer(i4b), intent(in)  ic,
real(dp), intent(in)  x,
real(dp), intent(in)  y,
real(dp), intent(in)  z 
)
private

Terminate with an error if the release point lies outside the given cell, or if the point is above or below the grid top or bottom, respectively.

Parameters
[in,out]thisthis instance
[in]iccell index
[in]zrelease point

Definition at line 414 of file prt-prp.f90.

415  class(PrtPrpType), intent(inout) :: this !< this instance
416  integer(I4B), intent(in) :: ic !< cell index
417  real(DP), intent(in) :: x, y, z !< release point
418  ! local
419  real(DP), allocatable :: polyverts(:, :)
420 
421  call this%fmi%dis%get_polyverts(ic, polyverts)
422  if (.not. point_in_polygon(x, y, polyverts)) then
423  write (errmsg, '(a,g0,a,g0,a,i0)') &
424  'Error: release point (x=', x, ', y=', y, ') is not in cell ', &
425  this%dis%get_nodeuser(ic)
426  call store_error(errmsg, terminate=.false.)
427  call store_error_filename(this%input_fname)
428  end if
429  if (z > maxval(this%dis%top)) then
430  write (errmsg, '(a,g0,a,g0,a,i0)') &
431  'Error: release point (z=', z, ') is above grid top ', &
432  maxval(this%dis%top)
433  call store_error(errmsg, terminate=.false.)
434  call store_error_filename(this%input_fname)
435  else if (z < minval(this%dis%bot)) then
436  write (errmsg, '(a,g0,a,g0,a,i0)') &
437  'Error: release point (z=', z, ') is below grid bottom ', &
438  minval(this%dis%bot)
439  call store_error(errmsg, terminate=.false.)
440  call store_error_filename(this%input_fname)
441  end if
442  deallocate (polyverts)
Here is the call graph for this function:

Variable Documentation

◆ default_exit_solve_tolerance

real(dp), parameter prtprpmodule::default_exit_solve_tolerance = DEM5
private

Definition at line 36 of file prt-prp.f90.

36  real(DP), parameter :: DEFAULT_EXIT_SOLVE_TOLERANCE = dem5

◆ ftype

character(len=lenftype) prtprpmodule::ftype = 'PRP'
private

Definition at line 34 of file prt-prp.f90.

34  character(len=LENFTYPE) :: ftype = 'PRP'

◆ text

character(len=16) prtprpmodule::text = ' PRP'
private

Definition at line 35 of file prt-prp.f90.

35  character(len=16) :: text = ' PRP'