MODFLOW 6  version 6.7.0.dev2
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'
 

Function/Subroutine Documentation

◆ define_listlabel()

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

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

642  class(PrtPrpType), intent(inout) :: this
643  ! 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 454 of file prt-prp.f90.

456  class(PrtPrpType), intent(inout) :: this !< this instance
457  type(ParticleType), pointer, intent(inout) :: particle !< the particle
458  integer(I4B), intent(in) :: ip !< particle index
459  real(DP), intent(in) :: trelease !< release time
460  ! local
461  integer(I4B) :: irow, icol, ilay, icpl
462  integer(I4B) :: ic, icu, ic_old
463  real(DP) :: x, y, z
464  real(DP) :: top, bot, hds
465 
466  ic = this%rptnode(ip)
467  icu = this%dis%get_nodeuser(ic)
468 
469  call create_particle(particle)
470 
471  if (size(this%boundname) /= 0) then
472  particle%name = this%boundname(ip)
473  else
474  particle%name = ''
475  end if
476 
477  particle%irpt = ip
478  particle%istopweaksink = this%istopweaksink
479  particle%istopzone = this%istopzone
480  particle%idrymeth = this%idrymeth
481  particle%icu = icu
482 
483  select type (dis => this%dis)
484  type is (distype)
485  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
486  type is (disvtype)
487  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
488  end select
489  particle%ilay = ilay
490  particle%izone = this%rptzone(ic)
491  particle%istatus = 0 ! status 0 until tracking starts
492  ! If the cell is inactive, either drape the particle
493  ! to the top-most active cell beneath it if drape is
494  ! enabled, or else terminate permanently unreleased.
495  if (this%ibound(ic) == 0) then
496  ic_old = ic
497  if (this%idrape > 0) then
498  call this%dis%highest_active(ic, this%ibound)
499  if (ic == ic_old .or. this%ibound(ic) == 0) then
500  ! negative unreleased status signals to the
501  ! tracking method that we haven't yet saved
502  ! a termination record, it needs to do so.
503  particle%istatus = -1 * term_unreleased
504  end if
505  else
506  particle%istatus = -1 * term_unreleased
507  end if
508  end if
509 
510  ! Load coordinates and transform if needed
511  x = this%rptx(ip)
512  y = this%rpty(ip)
513  if (this%ilocalz > 0) then
514  top = this%fmi%dis%top(ic)
515  bot = this%fmi%dis%bot(ic)
516  hds = this%fmi%gwfhead(ic)
517  z = bot + this%rptz(ip) * (hds - bot)
518  else
519  z = this%rptz(ip)
520  end if
521 
522  if (this%ichkmeth > 0) &
523  call this%validate_release_point(ic, x, y, z)
524 
525  particle%x = x
526  particle%y = y
527  particle%z = z
528  particle%trelease = trelease
529 
530  ! Set stop time to earlier of STOPTIME and STOPTRAVELTIME
531  if (this%stoptraveltime == huge(1d0)) then
532  particle%tstop = this%stoptime
533  else
534  particle%tstop = particle%trelease + this%stoptraveltime
535  if (this%stoptime < particle%tstop) particle%tstop = this%stoptime
536  end if
537 
538  particle%ttrack = particle%trelease
539  particle%itrdomain(level_model) = 0
540  particle%iboundary(level_model) = 0
541  particle%itrdomain(level_feature) = ic
542  particle%iboundary(level_feature) = 0
543  particle%itrdomain(level_subfeature) = 0
544  particle%iboundary(level_subfeature) = 0
545  particle%ifrctrn = this%ifrctrn
546  particle%iexmeth = this%iexmeth
547  particle%iextend = this%iextend
548  particle%icycwin = this%icycwin
549  particle%extol = this%extol
@ term_unreleased
terminated permanently unreleased
Definition: Particle.f90:35
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 378 of file prt-prp.f90.

379  class(PrtPrpType), intent(inout) :: this !< prp
380  if (this%iprpak > 0) then
381  write (this%iout, "(1x,/1x,a,1x,i0)") &
382  'PARTICLE RELEASE FOR PRP', this%ibcnum
383  call this%schedule%log(this%iout)
384  end if

◆ prp_ad()

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

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

318  use tdismodule, only: totalsimtime
319  class(PrtPrpType) :: this
320  integer(I4B) :: ip, it
321  real(DP) :: t
322 
323  ! Notes
324  ! -----
325  ! Each release point can be thought of as
326  ! a gumball machine with infinite supply:
327  ! a point can release an arbitrary number
328  ! of particles, but only one at any time.
329  ! Coincident release times are merged to
330  ! a single time by the release scheduler.
331 
332  ! Reset mass accumulators for this time step.
333  do ip = 1, this%nreleasepoints
334  this%rptm(ip) = dzero
335  end do
336 
337  ! Advance the release schedule and check if
338  ! any releases will be made this time step.
339  call this%schedule%advance()
340  if (.not. this%schedule%any()) return
341 
342  ! Log the schedule to the list file.
343  call this%log_release()
344 
345  ! Expand the particle store. We know from the
346  ! schedule how many particles will be released.
347  call this%particles%resize( &
348  this%particles%num_stored() + &
349  (this%nreleasepoints * this%schedule%count()), &
350  this%memoryPath)
351 
352  ! Release a particle from each point for
353  ! each release time in the current step.
354  do ip = 1, this%nreleasepoints
355  do it = 1, this%schedule%count()
356  t = this%schedule%times(it)
357  ! Skip the release time if it's before the simulation
358  ! starts, or if no `extend_tracking`, after it ends.
359  if (t < dzero) then
360  write (warnmsg, '(a,g0,a)') &
361  'Skipping negative release time (t=', t, ').'
362  call store_warning(warnmsg)
363  cycle
364  else if (t > totalsimtime .and. this%iextend == 0) then
365  write (warnmsg, '(a,g0,a)') &
366  'Skipping release time falling after the end of the &
367  &simulation (t=', t, '). Enable EXTEND_TRACKING to &
368  &release particles after the simulation end time.'
369  call store_warning(warnmsg)
370  cycle
371  end if
372  call this%release(ip, t)
373  end do
374  end do
real(dp), pointer, public totalsimtime
time at end of simulation
Definition: tdis.f90:37
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 204 of file prt-prp.f90.

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  call this%BndExtType%allocate_arrays()
213 
214  ! Allocate particle store, starting with the number
215  ! of release points (arrays resized if/when needed)
216  call create_particle_store( &
217  this%particles, &
218  this%nreleasepoints, &
219  this%memoryPath)
220 
221  ! Allocate arrays
222  call mem_allocate(this%rptx, this%nreleasepoints, 'RPTX', this%memoryPath)
223  call mem_allocate(this%rpty, this%nreleasepoints, 'RPTY', this%memoryPath)
224  call mem_allocate(this%rptz, this%nreleasepoints, 'RPTZ', this%memoryPath)
225  call mem_allocate(this%rptm, this%nreleasepoints, 'RPTMASS', &
226  this%memoryPath)
227  call mem_allocate(this%rptnode, this%nreleasepoints, 'RPTNODER', &
228  this%memoryPath)
229  call mem_allocate(this%rptname, lenboundname, this%nreleasepoints, &
230  'RPTNAME', this%memoryPath)
231 
232  ! Initialize arrays
233  do nps = 1, this%nreleasepoints
234  this%rptm(nps) = dzero
235  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 239 of file prt-prp.f90.

240  class(PrtPrpType) :: this
241 
242  ! Allocate parent's scalars
243  call this%BndExtType%allocate_scalars()
244 
245  ! Allocate scalars for this type
246  call mem_allocate(this%ilocalz, 'ILOCALZ', this%memoryPath)
247  call mem_allocate(this%iextend, 'IEXTEND', this%memoryPath)
248  call mem_allocate(this%offset, 'OFFSET', this%memoryPath)
249  call mem_allocate(this%stoptime, 'STOPTIME', this%memoryPath)
250  call mem_allocate(this%stoptraveltime, 'STOPTRAVELTIME', this%memoryPath)
251  call mem_allocate(this%istopweaksink, 'ISTOPWEAKSINK', this%memoryPath)
252  call mem_allocate(this%istopzone, 'ISTOPZONE', this%memoryPath)
253  call mem_allocate(this%idrape, 'IDRAPE', this%memoryPath)
254  call mem_allocate(this%idrymeth, 'IDRYMETH', this%memoryPath)
255  call mem_allocate(this%nreleasepoints, 'NRELEASEPOINTS', this%memoryPath)
256  call mem_allocate(this%nreleasetimes, 'NRELEASETIMES', this%memoryPath)
257  call mem_allocate(this%nparticles, 'NPARTICLES', this%memoryPath)
258  call mem_allocate(this%itrkout, 'ITRKOUT', this%memoryPath)
259  call mem_allocate(this%itrkhdr, 'ITRKHDR', this%memoryPath)
260  call mem_allocate(this%itrkcsv, 'ITRKCSV', this%memoryPath)
261  call mem_allocate(this%irlstls, 'IRLSTLS', this%memoryPath)
262  call mem_allocate(this%ifrctrn, 'IFRCTRN', this%memoryPath)
263  call mem_allocate(this%iexmeth, 'IEXMETH', this%memoryPath)
264  call mem_allocate(this%ichkmeth, 'ICHKMETH', this%memoryPath)
265  call mem_allocate(this%icycwin, 'ICYCWIN', this%memoryPath)
266  call mem_allocate(this%extol, 'EXTOL', this%memoryPath)
267  call mem_allocate(this%rttol, 'RTTOL', this%memoryPath)
268  call mem_allocate(this%rtfreq, 'RTFREQ', this%memoryPath)
269 
270  ! Set values
271  this%ilocalz = 0
272  this%iextend = 0
273  this%offset = dzero
274  this%stoptime = huge(1d0)
275  this%stoptraveltime = huge(1d0)
276  this%istopweaksink = 0
277  this%istopzone = 0
278  this%idrape = 0
279  this%idrymeth = 0
280  this%nreleasepoints = 0
281  this%nreleasetimes = 0
282  this%nparticles = 0
283  this%itrkout = 0
284  this%itrkhdr = 0
285  this%itrkcsv = 0
286  this%irlstls = 0
287  this%ifrctrn = 0
288  this%iexmeth = 0
289  this%ichkmeth = 1
290  this%icycwin = 0
291  this%extol = dem5
292  this%rttol = dsame * dep9
293  this%rtfreq = dzero
294 

◆ prp_ar()

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

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

299  ! dummy variables
300  class(PrtPrpType), intent(inout) :: this
301  ! local variables
302  integer(I4B) :: n
303 
304  call this%obs%obs_ar()
305 
306  if (this%inamedbound /= 0) then
307  do n = 1, this%nreleasepoints
308  this%boundname(n) = this%rptname(n)
309  end do
310  end if
311  do n = 1, this%nreleasepoints
312  this%nodelist(n) = this%rptnode(n)
313  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 607 of file prt-prp.f90.

608  ! modules
609  use tdismodule, only: delt
610  ! dummy variables
611  class(PrtPrpType) :: this
612  real(DP), dimension(:), intent(in) :: hnew
613  real(DP), dimension(:), intent(inout) :: flowja !< flow between package and model
614  integer(I4B), intent(in) :: imover !< flag indicating if the mover package is active
615  ! local variables
616  integer(I4B) :: i
617  integer(I4B) :: node
618  integer(I4B) :: idiag
619  real(DP) :: rrate
620 
621  ! If no boundaries, skip flow calculations.
622  if (this%nbound <= 0) return
623 
624  ! Loop through each boundary calculating flow.
625  do i = 1, this%nbound
626  node = this%nodelist(i)
627  rrate = dzero
628  ! If cell is no-flow or constant-head, then ignore it.
629  if (node > 0) then
630  ! Calculate the flow rate into the cell.
631  idiag = this%dis%con%ia(node)
632  rrate = this%rptm(i) * (done / delt) ! reciprocal of tstp length
633  flowja(idiag) = flowja(idiag) + rrate
634  end if
635 
636  ! Save simulated value to simvals array.
637  this%simvals(i) = rrate
638  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 99 of file prt-prp.f90.

101  ! dummy
102  class(BndType), pointer :: packobj
103  integer(I4B), intent(in) :: id
104  integer(I4B), intent(in) :: ibcnum
105  integer(I4B), intent(in) :: inunit
106  integer(I4B), intent(in) :: iout
107  character(len=*), intent(in) :: namemodel
108  character(len=*), intent(in) :: pakname
109  character(len=*), intent(in) :: input_mempath
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 MEMPATH: ', a, /)"
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, input_mempath)
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) input_mempath
Here is the caller graph for this function:

◆ prp_da()

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

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

148  class(PrtPrpType) :: this
149 
150  ! Deallocate parent
151  call this%BndExtType%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%ichkmeth)
173  call mem_deallocate(this%icycwin)
174  call mem_deallocate(this%extol)
175  call mem_deallocate(this%rttol)
176  call mem_deallocate(this%rtfreq)
177 
178  ! Deallocate arrays
179  call mem_deallocate(this%rptx)
180  call mem_deallocate(this%rpty)
181  call mem_deallocate(this%rptz)
182  call mem_deallocate(this%rptnode)
183  call mem_deallocate(this%rptm)
184  call mem_deallocate(this%rptname, 'RPTNAME', this%memoryPath)
185 
186  ! Deallocate objects
187  call this%particles%destroy(this%memoryPath)
188  call this%schedule%destroy()
189  deallocate (this%particles)
190  deallocate (this%schedule)

◆ prp_df_obs()

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

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

654  ! dummy
655  class(PrtPrpType) :: this
656  ! local
657  integer(I4B) :: indx
658  call this%obs%StoreObsType('prp', .true., indx)
659  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
660 
661  ! Store obs type and assign procedure pointer
662  ! for to-mvr observation type.
663  call this%obs%StoreObsType('to-mvr', .true., indx)
664  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 843 of file prt-prp.f90.

844  ! -- modules
847  ! -- dummy variables
848  class(PrtPrpType), intent(inout) :: this
849  ! -- local variables
850  type(PrtPrpParamFoundType) :: found
851 
852  call mem_set_value(this%nreleasepoints, 'NRELEASEPTS', this%input_mempath, &
853  found%nreleasepts)
854  call mem_set_value(this%nreleasetimes, 'NRELEASETIMES', this%input_mempath, &
855  found%nreleasetimes)
856 
857  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
858  write (this%iout, '(4x,a,i0)') 'NRELEASEPTS = ', this%nreleasepoints
859  write (this%iout, '(4x,a,i0)') 'NRELEASETIMES = ', this%nreleasetimes
860  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
861 
862  ! set maxbound and nbound to nreleasepts
863  this%maxbound = this%nreleasepoints
864  this%nbound = this%nreleasepoints
865 
866  ! allocate arrays for prp package
867  call this%prp_allocate_arrays()
868 
869  ! read packagedata and releasetimes blocks
870  call this%prp_packagedata()
871  call this%prp_releasetimes()
872  call this%prp_load_releasetimefrequency()

◆ prp_load_releasetimefrequency()

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

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

1052  ! modules
1053  use tdismodule, only: totalsimtime
1054  ! dummy
1055  class(PrtPrpType), intent(inout) :: this
1056  ! local
1057  real(DP), allocatable :: times(:)
1058 
1059  ! check if a release time frequency is configured
1060  if (this%rtfreq <= dzero) return
1061 
1062  ! create array of regularly-spaced release times
1063  times = arange( &
1064  start=dzero, &
1065  stop=totalsimtime, &
1066  step=this%rtfreq)
1067 
1068  ! register times with release schedule
1069  call this%schedule%time_select%extend(times)
1070 
1071  ! make sure times strictly increase
1072  if (.not. this%schedule%time_select%increasing()) then
1073  errmsg = "Release times must strictly increase"
1074  call store_error(errmsg)
1075  call store_error_filename(this%input_fname)
1076  end if
1077 
1078  ! deallocate
1079  deallocate (times)
1080 
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 807 of file prt-prp.f90.

808  ! -- modules
810  ! -- dummy variables
811  class(PrtPrpType), intent(inout) :: this
812  type(PrtPrpParamFoundType), intent(in) :: found
813  character(len=*), intent(in) :: trackfile
814  character(len=*), intent(in) :: trackcsvfile
815  ! -- local variables
816  ! formats
817  character(len=*), parameter :: fmttrkbin = &
818  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO BINARY FILE: ', a, /4x, &
819  &'OPENED ON UNIT: ', I0)"
820  character(len=*), parameter :: fmttrkcsv = &
821  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO CSV FILE: ', a, /4x, &
822  &'OPENED ON UNIT: ', I0)"
823 
824  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
825 
826  if (found%ifrctrn) then
827  write (this%iout, '(4x,a)') &
828  'IF DISV, TRACKING WILL USE THE TERNARY METHOD REGARDLESS OF CELL TYPE'
829  end if
830 
831  if (found%trackfile) then
832  write (this%iout, fmttrkbin) trim(adjustl(trackfile)), this%itrkout
833  end if
834 
835  if (found%trackcsvfile) then
836  write (this%iout, fmttrkcsv) trim(adjustl(trackcsvfile)), this%itrkcsv
837  end if
838 
839  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 647 of file prt-prp.f90.

648  class(PrtPrpType) :: this
649  prp_obs_supported = .true.

◆ prp_options()

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

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

669  ! -- modules
672  use openspecmodule, only: access, form
675  ! -- dummy variables
676  class(PrtPrpType), intent(inout) :: this
677  ! -- local variables
678  character(len=LENVARNAME), dimension(3) :: drytrack_method = &
679  &[character(len=LENVARNAME) :: 'DROP', 'STOP', 'STAY']
680  character(len=lenvarname), dimension(2) :: coorcheck_method = &
681  &[character(len=LENVARNAME) :: 'NONE', 'EAGER']
682  character(len=LINELENGTH) :: trackfile, trackcsvfile, fname
683  type(PrtPrpParamFoundType) :: found
684 
685  ! -- source base class options
686  call this%BndExtType%source_options()
687 
688  ! -- update defaults from input context
689  call mem_set_value(this%stoptime, 'STOPTIME', this%input_mempath, &
690  found%stoptime)
691  call mem_set_value(this%stoptraveltime, 'STOPTRAVELTIME', &
692  this%input_mempath, found%stoptraveltime)
693  call mem_set_value(this%istopweaksink, 'ISTOPWEAKSINK', this%input_mempath, &
694  found%istopweaksink)
695  call mem_set_value(this%istopzone, 'ISTOPZONE', this%input_mempath, &
696  found%istopzone)
697  call mem_set_value(this%idrape, 'DRAPE', this%input_mempath, &
698  found%drape)
699  call mem_set_value(this%idrymeth, 'IDRYMETH', this%input_mempath, &
700  drytrack_method, found%idrymeth)
701  call mem_set_value(trackfile, 'TRACKFILE', this%input_mempath, &
702  found%trackfile)
703  call mem_set_value(trackcsvfile, 'TRACKCSVFILE', this%input_mempath, &
704  found%trackcsvfile)
705  call mem_set_value(this%ilocalz, 'LOCAL_Z', this%input_mempath, &
706  found%local_z)
707  call mem_set_value(this%iextend, 'EXTEND_TRACKING', this%input_mempath, &
708  found%extend_tracking)
709  call mem_set_value(this%extol, 'EXTOL', this%input_mempath, &
710  found%extol)
711  call mem_set_value(this%rttol, 'RTTOL', this%input_mempath, &
712  found%rttol)
713  call mem_set_value(this%rtfreq, 'RTFREQ', this%input_mempath, &
714  found%rtfreq)
715  call mem_set_value(this%ifrctrn, 'IFRCTRN', this%input_mempath, &
716  found%ifrctrn)
717  call mem_set_value(this%iexmeth, 'IEXMETH', this%input_mempath, &
718  found%iexmeth)
719  call mem_set_value(this%ichkmeth, 'ICHKMETH', this%input_mempath, &
720  coorcheck_method, found%ichkmeth)
721  call mem_set_value(this%icycwin, 'ICYCWIN', this%input_mempath, found%icycwin)
722 
723  ! update internal state and validate input
724  if (found%idrymeth) then
725  if (this%idrymeth == 0) then
726  write (errmsg, '(a)') 'Unsupported dry tracking method. &
727  &DRY_TRACKING_METHOD must be "DROP", "STOP", or "STAY"'
728  call store_error(errmsg)
729  else
730  ! adjust for method zero indexing
731  this%idrymeth = this%idrymeth - 1
732  end if
733  end if
734 
735  if (found%extol) then
736  if (this%extol <= dzero) &
737  call store_error('EXIT_SOLVE_TOLERANCE MUST BE POSITIVE')
738  end if
739 
740  if (found%rttol) then
741  if (this%rttol <= dzero) &
742  call store_error('RELEASE_TIME_TOLERANCE MUST BE POSITIVE')
743  end if
744 
745  if (found%rtfreq) then
746  if (this%rtfreq <= dzero) &
747  call store_error('RELEASE_TIME_FREQUENCY MUST BE POSITIVE')
748  end if
749 
750  if (found%iexmeth) then
751  if (.not. (this%iexmeth /= 1 .or. this%iexmeth /= 2)) &
752  call store_error('DEV_EXIT_SOLVE_METHOD MUST BE &
753  &1 (BRENT) OR 2 (CHANDRUPATLA)')
754  end if
755 
756  if (found%ichkmeth) then
757  if (this%ichkmeth == 0) then
758  write (errmsg, '(a)') 'Unsupported coordinate check method. &
759  &COORDINATE_CHECK_METHOD must be "NONE" or "EAGER"'
760  call store_error(errmsg)
761  else
762  ! adjust for method zero based indexing
763  this%ichkmeth = this%ichkmeth - 1
764  end if
765  end if
766 
767  if (found%icycwin) then
768  if (this%icycwin < 0) &
769  call store_error('CYCLE_DETECTION_WINDOW MUST BE NON-NEGATIVE')
770  end if
771 
772  ! fileout options
773  if (found%trackfile) then
774  this%itrkout = getunit()
775  call openfile(this%itrkout, this%iout, trackfile, 'DATA(BINARY)', &
776  form, access, filstat_opt='REPLACE', &
777  mode_opt=mnormal)
778  ! open and write ascii header spec file
779  this%itrkhdr = getunit()
780  fname = trim(trackfile)//'.hdr'
781  call openfile(this%itrkhdr, this%iout, fname, 'CSV', &
782  filstat_opt='REPLACE', mode_opt=mnormal)
783  write (this%itrkhdr, '(a,/,a)') trackheader, trackdtypes
784  end if
785 
786  if (found%trackcsvfile) then
787  this%itrkcsv = getunit()
788  call openfile(this%itrkcsv, this%iout, trackcsvfile, 'CSV', &
789  filstat_opt='REPLACE')
790  write (this%itrkcsv, '(a)') trackheader
791  end if
792 
793  ! terminate if any errors were detected
794  if (count_errors() > 0) then
795  call store_error_filename(this%input_fname)
796  end if
797 
798  ! log found options
799  call this%prp_log_options(found, trackfile, trackcsvfile)
800 
801  ! Create release schedule now that we know
802  ! the coincident release time tolerance
803  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 876 of file prt-prp.f90.

878  use geomutilmodule, only: get_node
880  ! dummy
881  class(PrtPrpType), intent(inout) :: this
882  ! local
883  integer(I4B), dimension(:), pointer, contiguous :: irptno
884  integer(I4B), dimension(:, :), pointer, contiguous :: cellids
885  real(DP), dimension(:), pointer, contiguous :: xrpts, yrpts, zrpts
886  type(CharacterStringType), dimension(:), pointer, &
887  contiguous :: boundnames
888  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
889  character(len=9) :: cno
890  integer(I4B), dimension(:), allocatable :: nboundchk
891  integer(I4B), dimension(:), pointer :: cellid
892  integer(I4B) :: n, noder, nodeu, rptno
893 
894  ! set input context pointers
895  call mem_setptr(irptno, 'IRPTNO', this%input_mempath)
896  call mem_setptr(cellids, 'CELLID', this%input_mempath)
897  call mem_setptr(xrpts, 'XRPT', this%input_mempath)
898  call mem_setptr(yrpts, 'YRPT', this%input_mempath)
899  call mem_setptr(zrpts, 'ZRPT', this%input_mempath)
900  call mem_setptr(boundnames, 'BOUNDNAME', this%input_mempath)
901 
902  ! allocate and initialize temporary variables
903  allocate (nboundchk(this%nreleasepoints))
904  do n = 1, this%nreleasepoints
905  nboundchk(n) = 0
906  end do
907 
908  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%packName)) &
909  //' PACKAGEDATA'
910 
911  do n = 1, size(irptno)
912 
913  rptno = irptno(n)
914 
915  if (rptno < 1 .or. rptno > this%nreleasepoints) then
916  write (errmsg, '(a,i0,a,i0,a)') &
917  'Expected ', this%nreleasepoints, ' release points. &
918  &Points must be numbered from 1 to ', this%nreleasepoints, '.'
919  call store_error(errmsg)
920  cycle
921  end if
922 
923  ! increment nboundchk
924  nboundchk(rptno) = nboundchk(rptno) + 1
925 
926  ! set cellid
927  cellid => cellids(:, n)
928 
929  ! set node user
930  if (this%dis%ndim == 1) then
931  nodeu = cellid(1)
932  elseif (this%dis%ndim == 2) then
933  nodeu = get_node(cellid(1), 1, cellid(2), &
934  this%dis%mshape(1), 1, &
935  this%dis%mshape(2))
936  else
937  nodeu = get_node(cellid(1), cellid(2), cellid(3), &
938  this%dis%mshape(1), &
939  this%dis%mshape(2), &
940  this%dis%mshape(3))
941  end if
942 
943  ! set noder
944  noder = this%dis%get_nodenumber(nodeu, 1)
945  if (noder <= 0) then
946  cycle
947  else
948  this%rptnode(rptno) = noder
949  end if
950 
951  if (this%ilocalz > 0 .and. (zrpts(n) < 0 .or. zrpts(n) > 1)) then
952  call store_error('Local z coordinate must fall in the interval [0, 1]')
953  cycle
954  end if
955 
956  ! set coordinates
957  this%rptx(rptno) = xrpts(n)
958  this%rpty(rptno) = yrpts(n)
959  this%rptz(rptno) = zrpts(n)
960 
961  ! set default boundname
962  write (cno, '(i9.9)') rptno
963  bndname = 'PRP'//cno
964 
965  ! read boundnames from file, if provided
966  if (this%inamedbound /= 0) then
967  bndnametemp = boundnames(n)
968  if (bndnametemp /= '') bndname = bndnametemp
969  else
970  bndname = ''
971  end if
972 
973  ! set boundname
974  this%rptname(rptno) = bndname
975  end do
976 
977  write (this%iout, '(1x,a)') &
978  'END OF '//trim(adjustl(this%packName))//' PACKAGEDATA'
979 
980  ! check for duplicate or missing particle release points
981  do n = 1, this%nreleasepoints
982  if (nboundchk(n) == 0) then
983  write (errmsg, '(a,a,1x,i0,a)') 'No data specified for particle ', &
984  'release point', n, '.'
985  call store_error(errmsg)
986  else if (nboundchk(n) > 1) then
987  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
988  'Data for particle release point', n, 'specified', nboundchk(n), &
989  'times.'
990  call store_error(errmsg)
991  end if
992  end do
993 
994  ! terminate if any errors were detected
995  if (count_errors() > 0) then
996  call store_error_filename(this%input_fname)
997  end if
998 
999  ! cleanup
1000  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 1004 of file prt-prp.f90.

1006  ! dummy
1007  class(PrtPrpType), intent(inout) :: this
1008  ! local
1009  real(DP), dimension(:), pointer, contiguous :: time
1010  integer(I4B) :: n, isize
1011  real(DP), allocatable :: times(:)
1012 
1013  if (this%nreleasetimes <= 0) return
1014 
1015  ! allocate times array
1016  allocate (times(this%nreleasetimes))
1017 
1018  ! check if input array was read
1019  call get_isize('TIME', this%input_mempath, isize)
1020 
1021  if (isize <= 0) then
1022  errmsg = "RELEASTIMES block expected when &
1023  &NRELEASETIMES dimension is non-zero."
1024  call store_error(errmsg)
1025  call store_error_filename(this%input_fname)
1026  end if
1027 
1028  ! set input context pointer
1029  call mem_setptr(time, 'TIME', this%input_mempath)
1030 
1031  ! set input data
1032  do n = 1, size(time)
1033  times(n) = time(n)
1034  end do
1035 
1036  ! register times with the release schedule
1037  call this%schedule%time_select%extend(times)
1038 
1039  ! make sure times strictly increase
1040  if (.not. this%schedule%time_select%increasing()) then
1041  errmsg = "RELEASTIMES block entries must strictly increase."
1042  call store_error(errmsg)
1043  call store_error_filename(this%input_fname)
1044  end if
1045 
1046  ! deallocate
1047  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 553 of file prt-prp.f90.

554  ! modules
555  use tdismodule, only: kper, nper
558  ! dummy variables
559  class(PrtPrpType), intent(inout) :: this
560  ! local variables
561  type(CharacterStringType), dimension(:), contiguous, &
562  pointer :: settings
563  integer(I4B), pointer :: iper, ionper, nlist
564  character(len=LINELENGTH), allocatable :: lines(:)
565  integer(I4B) :: n
566 
567  ! set pointer to last and next period loaded
568  call mem_setptr(iper, 'IPER', this%input_mempath)
569  call mem_setptr(ionper, 'IONPER', this%input_mempath)
570 
571  if (kper == 1 .and. &
572  (iper == 0) .and. &
573  (ionper > nper) .and. &
574  size(this%schedule%time_select%times) == 0) then
575  ! If the user hasn't provided any release settings (neither
576  ! explicit release times, release time frequency, or period
577  ! block release settings), default to a single release at the
578  ! start of the first period's first time step.
579  allocate (lines(1))
580  lines(1) = "FIRST"
581  call this%schedule%advance(lines=lines)
582  deallocate (lines)
583  return
584  else if (iper /= kper) then
585  return
586  end if
587 
588  ! set input context pointers
589  call mem_setptr(nlist, 'NBOUND', this%input_mempath)
590  call mem_setptr(settings, 'SETTING', this%input_mempath)
591 
592  ! allocate and set input
593  allocate (lines(nlist))
594  do n = 1, nlist
595  lines(n) = settings(n)
596  end do
597 
598  ! update schedule
599  if (size(lines) > 0) &
600  call this%schedule%advance(lines=lines)
601 
602  ! cleanup
603  deallocate (lines)
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
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 194 of file prt-prp.f90.

195  class(PrtPrpType) :: this
196  integer(I4B), dimension(:), pointer, contiguous :: ibound
197  integer(I4B), dimension(:), pointer, contiguous :: izone
198 
199  this%ibound => ibound
200  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 436 of file prt-prp.f90.

437  ! dummy
438  class(PrtPrpType), intent(inout) :: this !< this instance
439  integer(I4B), intent(in) :: ip !< particle index
440  real(DP), intent(in) :: trelease !< release time
441  ! local
442  integer(I4B) :: np
443  type(ParticleType), pointer :: particle
444 
445  call this%initialize_particle(particle, ip, trelease)
446  np = this%nparticles + 1
447  this%nparticles = np
448  call this%particles%put(particle, np)
449  deallocate (particle)
450  this%rptm(ip) = this%rptm(ip) + done ! TODO configurable mass
451 

◆ 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 393 of file prt-prp.f90.

394  class(PrtPrpType), intent(inout) :: this !< this instance
395  integer(I4B), intent(in) :: ic !< cell index
396  real(DP), intent(in) :: x, y, z !< release point
397  ! local
398  real(DP), allocatable :: polyverts(:, :)
399 
400  call this%fmi%dis%get_polyverts(ic, polyverts)
401  if (.not. point_in_polygon(x, y, polyverts)) then
402  write (errmsg, '(a,g0,a,g0,a,i0)') &
403  'Error: release point (x=', x, ', y=', y, ') is not in cell ', &
404  this%dis%get_nodeuser(ic)
405  call store_error(errmsg, terminate=.false.)
406  call store_error_filename(this%input_fname)
407  end if
408  if (z > maxval(this%dis%top)) then
409  write (errmsg, '(a,g0,a,g0,a,i0)') &
410  'Error: release point (z=', z, ') is above grid top ', &
411  maxval(this%dis%top)
412  call store_error(errmsg, terminate=.false.)
413  call store_error_filename(this%input_fname)
414  else if (z < minval(this%dis%bot)) then
415  write (errmsg, '(a,g0,a,g0,a,i0)') &
416  'Error: release point (z=', z, ') is below grid bottom ', &
417  minval(this%dis%bot)
418  call store_error(errmsg, terminate=.false.)
419  call store_error_filename(this%input_fname)
420  end if
421  deallocate (polyverts)
Here is the call graph for this function:

Variable Documentation

◆ 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'