MODFLOW 6  version 6.7.0.dev0
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, fmi)
 Create a new particle release point package. More...
 
subroutine prp_da (this)
 Deallocate memory. More...
 
subroutine prp_set_pointers (this, ibound, izone, trackctl)
 @ 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, option, found)
 Set options specific to PrtPrpType. More...
 
subroutine prp_read_dimensions (this)
 Read package dimensions. More...
 
subroutine prp_read_packagedata (this)
 Load package data (release points). More...
 
subroutine prp_read_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 668 of file prt-prp.f90.

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

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

372  class(PrtPrpType), intent(inout) :: this !< prp
373  if (this%iprpak > 0) then
374  write (this%iout, "(1x,/1x,a,1x,i0)") &
375  'PARTICLE RELEASE FOR PRP', this%ibcnum
376  call this%schedule%log(this%iout)
377  end if

◆ prp_ad()

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

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

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

◆ prp_allocate_scalars()

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

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

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 

◆ prp_ar()

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

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

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

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

635  ! modules
636  use tdismodule, only: delt
637  ! dummy variables
638  class(PrtPrpType) :: this
639  real(DP), dimension(:), intent(in) :: hnew
640  real(DP), dimension(:), intent(inout) :: flowja !< flow between package and model
641  integer(I4B), intent(in) :: imover !< flag indicating if the mover package is active
642  ! local variables
643  integer(I4B) :: i
644  integer(I4B) :: node
645  integer(I4B) :: idiag
646  real(DP) :: rrate
647 
648  ! If no boundaries, skip flow calculations.
649  if (this%nbound <= 0) return
650 
651  ! Loop through each boundary calculating flow.
652  do i = 1, this%nbound
653  node = this%nodelist(i)
654  rrate = dzero
655  ! If cell is no-flow or constant-head, then ignore it.
656  if (node > 0) then
657  ! Calculate the flow rate into the cell.
658  idiag = this%dis%con%ia(node)
659  rrate = this%rptm(i) * (done / delt) ! reciprocal of tstp length
660  flowja(idiag) = flowja(idiag) + rrate
661  end if
662 
663  ! Save simulated value to simvals array.
664  this%simvals(i) = rrate
665  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,
type(prtfmitype), pointer  fmi 
)

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

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
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%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)

◆ prp_df_obs()

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

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

681  ! dummy
682  class(PrtPrpType) :: this
683  ! local
684  integer(I4B) :: indx
685  call this%obs%StoreObsType('prp', .true., indx)
686  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
687 
688  ! Store obs type and assign procedure pointer
689  ! for to-mvr observation type.
690  call this%obs%StoreObsType('to-mvr', .true., indx)
691  this%obs%obsData(indx)%ProcessIdPtr => defaultobsidprocessor
Here is the call graph for this function:

◆ prp_load_releasetimefrequency()

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

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

1090  ! modules
1091  use tdismodule, only: totalsimtime
1092  ! dummy
1093  class(PrtPrpType), intent(inout) :: this
1094  ! local
1095  real(DP), allocatable :: times(:)
1096 
1097  ! check if a release time frequency is configured
1098  if (this%rtfreq <= dzero) return
1099 
1100  ! create array of regularly-spaced release times
1101  times = arange( &
1102  start=dzero, &
1103  stop=totalsimtime, &
1104  step=this%rtfreq)
1105 
1106  ! register times with release schedule
1107  call this%schedule%time_select%extend(times)
1108 
1109  ! make sure times strictly increase
1110  if (.not. this%schedule%time_select%increasing()) then
1111  errmsg = "Release times must strictly increase"
1112  call store_error(errmsg)
1113  call this%parser%StoreErrorUnit(terminate=.true.)
1114  end if
1115 
1116  ! deallocate
1117  deallocate (times)
1118 
Here is the call graph for this function:

◆ prp_obs_supported()

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

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

675  class(PrtPrpType) :: this
676  prp_obs_supported = .true.

◆ prp_options()

subroutine prtprpmodule::prp_options ( class(prtprptype), intent(inout)  this,
character(len=*), intent(inout)  option,
logical(lgp), intent(inout)  found 
)
private

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

696  use openspecmodule, only: access, form
697  use constantsmodule, only: maxcharlen, dzero
700  ! dummy
701  class(PrtPrpType), intent(inout) :: this
702  character(len=*), intent(inout) :: option
703  logical(LGP), intent(inout) :: found
704  ! locals
705  character(len=MAXCHARLEN) :: fname
706  character(len=MAXCHARLEN) :: keyword
707  ! formats
708  character(len=*), parameter :: fmttrkbin = &
709  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO BINARY FILE: ', a, /4x, &
710  &'OPENED ON UNIT: ', I0)"
711  character(len=*), parameter :: fmttrkcsv = &
712  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO CSV FILE: ', a, /4x, &
713  &'OPENED ON UNIT: ', I0)"
714 
715  select case (option)
716  case ('STOPTIME')
717  this%stoptime = this%parser%GetDouble()
718  found = .true.
719  case ('STOPTRAVELTIME')
720  this%stoptraveltime = this%parser%GetDouble()
721  found = .true.
722  case ('STOP_AT_WEAK_SINK')
723  this%istopweaksink = 1
724  found = .true.
725  case ('ISTOPZONE')
726  this%istopzone = this%parser%GetInteger()
727  found = .true.
728  case ('DRAPE')
729  this%idrape = 1
730  found = .true.
731  case ('DRY_TRACKING_METHOD')
732  call this%parser%GetStringCaps(keyword)
733  select case (keyword)
734  case ('DROP')
735  this%idrymeth = 0
736  case ('STOP')
737  this%idrymeth = 1
738  case ('STAY')
739  this%idrymeth = 2
740  case default
741  write (errmsg, '(a, a)') &
742  'Unknown dry tracking method: ', trim(keyword)
743  call store_error(errmsg)
744  write (errmsg, '(a, a)') &
745  'DRY must be "DROP", "STOP" or "STAY"'
746  call store_error(errmsg)
747  call this%parser%StoreErrorUnit()
748  end select
749  found = .true.
750  case ('TRACK')
751  call this%parser%GetStringCaps(keyword)
752  if (keyword == 'FILEOUT') then
753  ! parse filename
754  call this%parser%GetString(fname)
755  ! open binary output file
756  this%itrkout = getunit()
757  call openfile(this%itrkout, this%iout, fname, 'DATA(BINARY)', &
758  form, access, filstat_opt='REPLACE', &
759  mode_opt=mnormal)
760  write (this%iout, fmttrkbin) trim(adjustl(fname)), this%itrkout
761  ! open and write ascii header spec file
762  this%itrkhdr = getunit()
763  fname = trim(fname)//'.hdr'
764  call openfile(this%itrkhdr, this%iout, fname, 'CSV', &
765  filstat_opt='REPLACE', mode_opt=mnormal)
766  write (this%itrkhdr, '(a,/,a)') trackheader, trackdtypes
767  else
768  call store_error('OPTIONAL TRACK KEYWORD MUST BE '// &
769  'FOLLOWED BY FILEOUT')
770  end if
771  found = .true.
772  case ('TRACKCSV')
773  call this%parser%GetStringCaps(keyword)
774  if (keyword == 'FILEOUT') then
775  ! parse filename
776  call this%parser%GetString(fname)
777  ! open CSV output file and write headers
778  this%itrkcsv = getunit()
779  call openfile(this%itrkcsv, this%iout, fname, 'CSV', &
780  filstat_opt='REPLACE')
781  write (this%iout, fmttrkcsv) trim(adjustl(fname)), this%itrkcsv
782  write (this%itrkcsv, '(a)') trackheader
783  else
784  call store_error('OPTIONAL TRACKCSV KEYWORD MUST BE &
785  &FOLLOWED BY FILEOUT')
786  end if
787  found = .true.
788  case ('LOCAL_Z')
789  this%ilocalz = 1
790  found = .true.
791  case ('EXTEND_TRACKING')
792  this%iextend = 1
793  found = .true.
794  case ('EXIT_SOLVE_TOLERANCE')
795  this%extol = this%parser%GetDouble()
796  if (this%extol <= dzero) &
797  call store_error('EXIT_SOLVE_TOLERANCE MUST BE POSITIVE')
798  found = .true.
799  case ('RELEASE_TIME_TOLERANCE')
800  this%rttol = this%parser%GetDouble()
801  if (this%rttol <= dzero) &
802  call store_error('RELEASE_TIME_TOLERANCE MUST BE POSITIVE')
803  found = .true.
804  case ('RELEASE_TIME_FREQUENCY')
805  this%rtfreq = this%parser%GetDouble()
806  if (this%rtfreq <= dzero) &
807  call store_error('RELEASE_TIME_FREQUENCY MUST BE POSITIVE')
808  found = .true.
809  case ('DEV_FORCETERNARY')
810  call this%parser%DevOpt()
811  this%ifrctrn = 1
812  write (this%iout, '(4x,a)') &
813  'TRACKING WILL BE DONE USING THE TERNARY METHOD REGARDLESS OF CELL TYPE'
814  found = .true.
815  case ('DEV_EXIT_SOLVE_METHOD')
816  call this%parser%DevOpt()
817  this%iexmeth = this%parser%GetInteger()
818  if (.not. (this%iexmeth /= 1 .or. this%iexmeth /= 2)) &
819  call store_error('DEV_EXIT_SOLVE_METHOD MUST BE &
820  &1 (BRENT) OR 2 (CHANDRUPATLA)')
821  found = .true.
822  case default
823  found = .false.
824  end select
825 
826  ! Catch unrecognized options
827  if (.not. found) then
828  errmsg = "UNKNOWN PRP OPTION '"//trim(keyword)//"'."
829  call store_error(errmsg)
830  call this%parser%StoreErrorUnit()
831  end if
832 
833  ! Create release schedule now that we know
834  ! the coincident release time tolerance
835  this%schedule => create_release_schedule(tol=this%rttol)
836 
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:47
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.
character(len=20) access
Definition: OpenSpec.f90:7
character(len=20) form
Definition: OpenSpec.f90:7
character(len= *), parameter, public trackdtypes
Definition: TrackFile.f90:81
character(len= *), parameter, public trackheader
Definition: TrackFile.f90:77
Here is the call graph for this function:

◆ prp_read_dimensions()

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

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

841  ! dummy
842  class(PrtPrpType), intent(inout) :: this
843  ! local
844  character(len=LINELENGTH) :: errmsg, keyword
845  integer(I4B) :: ierr
846  logical :: isfound, endOfBlock
847 
848  ! get dimension block
849  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
850  supportopenclose=.true.)
851 
852  ! parse dimension block if detected
853  if (isfound) then
854  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
855  do
856  call this%parser%GetNextLine(endofblock)
857  if (endofblock) exit
858  call this%parser%GetStringCaps(keyword)
859  select case (keyword)
860  case ('NRELEASEPTS')
861  this%nreleasepoints = this%parser%GetInteger()
862  case ('NRELEASETIMES')
863  this%nreleasetimes = this%parser%GetInteger()
864  case default
865  write (errmsg, &
866  '(4x,a,a)') '****ERROR. UNKNOWN PARTICLE INPUT DIMENSION: ', &
867  trim(keyword)
868  call store_error(errmsg)
869  call this%parser%StoreErrorUnit()
870  end select
871  end do
872  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
873  else
874  call store_error('ERROR. REQUIRED DIMENSIONS BLOCK NOT FOUND.')
875  end if
876 
877  ! set maxbound and nbound to nreleasepts
878  this%maxbound = this%nreleasepoints
879  this%nbound = this%nreleasepoints
880 
881  ! allocate arrays for prp package
882  call this%prp_allocate_arrays()
883 
884  ! read packagedata and releasetimes blocks
885  call this%prp_read_packagedata()
886  call this%prp_read_releasetimes()
887  call this%prp_load_releasetimefrequency()
Here is the call graph for this function:

◆ prp_read_packagedata()

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

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

892  ! dummy
893  class(PrtPrpType), intent(inout) :: this
894  ! local
895  character(len=LINELENGTH) :: cellid
896  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
897  character(len=9) :: cno
898  logical :: isfound
899  logical :: endOfBlock
900  integer(I4B) :: ival
901  integer(I4B) :: n
902  integer(I4B) :: ierr
903  character(len=LENBOUNDNAME), dimension(:), allocatable :: nametxt
904  integer(I4B), dimension(:), allocatable :: nboundchk
905  integer(I4B), dimension(:), allocatable :: noder
906  real(DP), dimension(:), allocatable :: x
907  real(DP), dimension(:), allocatable :: y
908  real(DP), dimension(:), allocatable :: z
909  real(DP), dimension(:), allocatable :: tstop
910  ! format
911  character(len=*), parameter :: fmttend = &
912  "('end time (', G0, ') must be greater than or equal to the &
913  &begin time (', G0, ').')"
914 
915  ! allocate temporary variables
916  allocate (noder(this%nreleasepoints))
917  allocate (x(this%nreleasepoints))
918  allocate (y(this%nreleasepoints))
919  allocate (z(this%nreleasepoints))
920  allocate (tstop(this%nreleasepoints))
921  allocate (nametxt(this%nreleasepoints))
922  allocate (nboundchk(this%nreleasepoints))
923 
924  ! initialize temporary variables
925  do n = 1, this%nreleasepoints
926  nboundchk(n) = 0
927  end do
928 
929  ! read particle release point data
930  ! get particle release points block
931  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
932  supportopenclose=.true.)
933 
934  ! parse block if detected
935  if (isfound) then
936  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%packName)) &
937  //' PACKAGEDATA'
938  do
939  call this%parser%GetNextLine(endofblock)
940  if (endofblock) exit
941  ival = this%parser%GetInteger()
942  n = ival
943 
944  if (n < 1 .or. n > this%nreleasepoints) then
945  write (errmsg, '(a,i0,a,i0,a)') &
946  'Expected ', this%nreleasepoints, ' release points. &
947  &Points must be numbered from 1 to ', this%nreleasepoints, '.'
948  call store_error(errmsg)
949  cycle
950  end if
951 
952  ! increment nboundchk
953  nboundchk(n) = nboundchk(n) + 1
954 
955  ! node number
956  call this%parser%GetCellid(this%dis%ndim, cellid)
957  noder(n) = this%dis%noder_from_cellid(cellid, this%inunit, this%iout)
958 
959  ! x, y, z coordinates
960  x(n) = this%parser%GetDouble()
961  y(n) = this%parser%GetDouble()
962  z(n) = this%parser%GetDouble()
963 
964  if (this%ilocalz > 0 .and. (z(n) < 0 .or. z(n) > 1)) then
965  call store_error('Local z coordinate must fall in the interval [0, 1]')
966  cycle
967  end if
968 
969  ! set default boundname
970  write (cno, '(i9.9)') n
971  bndname = 'PRP'//cno
972 
973  ! read boundnames from file, if provided
974  if (this%inamedbound /= 0) then
975  call this%parser%GetStringCaps(bndnametemp)
976  if (bndnametemp /= '') &
977  bndname = bndnametemp
978  else
979  bndname = ''
980  end if
981 
982  ! store temp boundnames
983  nametxt(n) = bndname
984  end do
985 
986  write (this%iout, '(1x,a)') &
987  'END OF '//trim(adjustl(this%packName))//' PACKAGEDATA'
988 
989  ! check for duplicate or missing particle release points
990  do n = 1, this%nreleasepoints
991  if (nboundchk(n) == 0) then
992  write (errmsg, '(a,a,1x,i0,a)') 'No data specified for particle ', &
993  'release point', n, '.'
994  call store_error(errmsg)
995  else if (nboundchk(n) > 1) then
996  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
997  'Data for particle release point', n, 'specified', nboundchk(n), &
998  'times.'
999  call store_error(errmsg)
1000  end if
1001  end do
1002  else
1003  call store_error('Required packagedata block not found.')
1004  end if
1005 
1006  ! terminate if any errors were detected
1007  if (count_errors() > 0) then
1008  call this%parser%StoreErrorUnit()
1009  end if
1010 
1011  ! fill particle release point data with data stored in temporary local arrays
1012  do n = 1, this%nreleasepoints
1013  this%rptnode(n) = noder(n)
1014  this%rptx(n) = x(n)
1015  this%rpty(n) = y(n)
1016  this%rptz(n) = z(n)
1017  this%rptname(n) = nametxt(n)
1018  end do
1019 
1020  ! deallocate local storage
1021  deallocate (noder)
1022  deallocate (x)
1023  deallocate (y)
1024  deallocate (z)
1025  deallocate (tstop)
1026  deallocate (nametxt)
1027  deallocate (nboundchk)
Here is the call graph for this function:

◆ prp_read_releasetimes()

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

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

1032  ! dummy
1033  class(PrtPrpType), intent(inout) :: this
1034  ! local
1035  integer(I4B) :: i, ierr
1036  logical(LGP) :: eob, found, success
1037  real(DP) :: t
1038  real(DP), allocatable :: times(:)
1039 
1040  ! get releasetimes block
1041  call this%parser%GetBlock('RELEASETIMES', found, ierr, &
1042  supportopenclose=.true., &
1043  blockrequired=.false.)
1044 
1045  ! raise an error if releasetimes has a dimension
1046  ! but no block was found, otherwise return early
1047  if (.not. found) then
1048  if (this%nreleasetimes <= 0) return
1049  write (errmsg, '(a, i0)') &
1050  "Expected RELEASETIMES with length ", this%nreleasetimes
1051  call store_error(errmsg)
1052  call this%parser%StoreErrorUnit(terminate=.true.)
1053  end if
1054 
1055  ! allocate times array
1056  allocate (times(this%nreleasetimes))
1057 
1058  ! read times from the block
1059  write (this%iout, '(/1x,a)') &
1060  'PROCESSING '//trim(adjustl(this%text))//' RELEASETIMES'
1061  do i = 1, this%nreleasetimes
1062  call this%parser%GetNextLine(eob)
1063  if (eob) exit
1064  call this%parser%TryGetDouble(t, success)
1065  if (.not. success) then
1066  errmsg = "Failed to read double precision value"
1067  call store_error(errmsg)
1068  call this%parser%StoreErrorUnit(terminate=.true.)
1069  end if
1070  times(i) = t
1071  end do
1072 
1073  ! register times with the release schedule
1074  call this%schedule%time_select%extend(times)
1075 
1076  ! make sure times strictly increase
1077  if (.not. this%schedule%time_select%increasing()) then
1078  errmsg = "Release times must strictly increase"
1079  call store_error(errmsg)
1080  call this%parser%StoreErrorUnit(terminate=.true.)
1081  end if
1082 
1083  ! deallocate
1084  deallocate (times)
1085 
Here is the call graph for this function:

◆ prp_rp()

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

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

547  ! modules
548  use tdismodule, only: kper, nper
549  use inputoutputmodule, only: urword
550  ! dummy variables
551  class(PrtPrpType), intent(inout) :: this
552  ! local variables
553  integer(I4B) :: ierr
554  logical(LGP) :: is_found
555  logical(LGP) :: end_of_block
556  logical(LGP) :: no_blocks
557  character(len=LINELENGTH) :: line
558  character(len=LINELENGTH), allocatable :: lines(:)
559  ! formats
560  character(len=*), parameter :: fmtblkerr = &
561  "('Looking for BEGIN PERIOD iper. &
562  &Found ', a, ' instead.')"
563  character(len=*), parameter :: fmt_steps = &
564  "(6x,'TIME STEP(S) ',50(I0,' '))" ! 50 limit is similar to STEPS in OC
565  character(len=*), parameter :: fmt_freq = &
566  "(6x,'EVERY ',I0,' TIME STEP(S)')"
567  character(len=*), parameter :: fmt_fracs = &
568  "(6x,50(f10.3,' '))"
569 
570  ! Set ionper to the stress period number for which a new block of data
571  ! will be read.
572  if (this%inunit == 0) return
573 
574  ! get stress period data
575  no_blocks = .false.
576  if (this%ionper < kper) then
577  ! get period block
578  call this%parser%GetBlock('PERIOD', is_found, ierr, &
579  supportopenclose=.true., &
580  blockrequired=.false.)
581  if (is_found) then
582  ! read ionper and check for increasing period numbers
583  call this%read_check_ionper()
584  else
585  ! PERIOD block not found
586  if (ierr < 0) then
587  if (kper == 1) then
588  ! End of file found; no period data for the simulation.
589  no_blocks = .true.
590  else
591  ! End of file found; no more period data.
592  this%ionper = nper + 1
593  end if
594  else
595  ! Found invalid block
596  call this%parser%GetCurrentLine(line)
597  write (errmsg, fmtblkerr) adjustl(trim(line))
598  call store_error(errmsg, terminate=.true.)
599  end if
600  end if
601  end if
602 
603  ! If the user hasn't provided any release settings (neither
604  ! explicit release times, release time frequency, or period
605  ! block release settings), default to a single release at the
606  ! start of the first period's first time step.
607  if (no_blocks .and. &
608  kper == 1 .and. &
609  size(this%schedule%time_select%times) == 0) then
610  allocate (lines(1))
611  line = "FIRST"
612  lines(1) = line
613  call this%schedule%advance(lines=lines)
614  else if (this%ionper == kper) then
615  ! If the current stress period matches the
616  ! block we are reading, parse the setting
617  ! and register it with the schedule.
618  allocate (lines(0))
619  recordloop: do
620  call this%parser%GetNextLine(end_of_block)
621  if (end_of_block) exit recordloop
622  call this%parser%GetCurrentLine(line)
623  call expandarray(lines)
624  lines(size(lines)) = line
625  end do recordloop
626  if (size(lines) > 0) &
627  call this%schedule%advance(lines=lines)
628  deallocate (lines)
629  end if
630 
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
Here is the call graph for this function:

◆ prp_set_pointers()

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

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

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

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

430  ! dummy
431  class(PrtPrpType), intent(inout) :: this !< this instance
432  integer(I4B), intent(in) :: ip !< particle index
433  real(DP), intent(in) :: trelease !< release time
434  ! local
435  integer(I4B) :: np
436  type(ParticleType), pointer :: particle
437 
438  call this%initialize_particle(particle, ip, trelease)
439 
440  ! Increment cumulative particle count
441  np = this%nparticles + 1
442  this%nparticles = np
443 
444  ! Save the particle to the store
445  call this%particles%save_particle(particle, np)
446  deallocate (particle)
447 
448  ! Accumulate mass for release point
449  this%rptm(ip) = this%rptm(ip) + done
450 

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

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

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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