MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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, mempath, fmi)
 Create a new particle release point package. More...
 
subroutine prp_da (this)
 Deallocate memory. More...
 
subroutine prp_set_pointers (this, ibound, izone, trackfilectl)
 @ 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 appropriate. More...
 
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)
 @ brief Define list heading written with PRINT_INPUT option More...
 
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_packagedata (this)
 Read the packagedata for this package. More...
 
subroutine prp_read_dimensions (this)
 Read package dimensions. 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 651 of file prt-prp.f90.

652  class(PrtPrpType), intent(inout) :: this
653 
654  ! -- create the header list label
655  this%listlabel = trim(this%filtyp)//' NO.'
656  if (this%dis%ndim == 3) then
657  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
658  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'ROW'
659  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'COL'
660  elseif (this%dis%ndim == 2) then
661  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'LAYER'
662  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'CELL2D'
663  else
664  write (this%listlabel, '(a, a7)') trim(this%listlabel), 'NODE'
665  end if
666  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'STRESS RATE'
667  if (this%inamedbound == 1) then
668  write (this%listlabel, '(a, a16)') trim(this%listlabel), 'BOUNDARY NAME'
669  end if

◆ prp_ad()

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

Releases may be scheduled via a global RELEASETIME, or within a stress period via ALL, FIRST, FREQUENCY or STEPS (with optional FRACTION). If no release option is specified, a single release is conducted at the first moment of the first time step of the first stress period.

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

303  ! -- modules
304  use tdismodule, only: totimc, delt, kstp
305  use dismodule, only: distype
306  use disvmodule, only: disvtype
307  ! -- dummy
308  class(PrtPrpType) :: this
309  ! -- local
310  character(len=LINELENGTH) :: errmsg
311  integer(I4B) :: ic, icu, nps, nts, nrel, &
312  nreleasets, np, irow, icol, ilay, icpl
313  real(DP) :: x, y, z, trelease, tend, top, bot, hds
314  real(DP), allocatable :: polyverts(:, :)
315  type(ParticleType), pointer :: particle
316 
317  ! -- Check if there's a release to make
318  if (.not. ( &
319  ! all time steps?
320  this%rlsall .or. &
321  ! first time step?
322  (this%rlsfirst .and. kstp == 1) .or. &
323  ! specified time steps?
324  any(this%rlskstp == kstp) .or. &
325  ! specified release times?
326  this%rlstimelist)) return
327 
328  if (this%rlstimelist) then
329  nreleasets = size(this%releasetimes%times)
330  else
331  nreleasets = 1
332  end if
333  nrel = this%nreleasepts * nreleasets
334 
335  ! -- Reset mass release for time step
336  do nps = 1, this%nreleasepts
337  this%rptmass(nps) = dzero
338  end do
339 
340  ! -- Resize particle store if another set
341  ! of particles will exceed its capacity
342  if ((this%nparticles + nrel) > size(this%particles%irpt)) &
343  call this%particles%resize( &
344  size(this%particles%irpt) + nrel, &
345  this%memoryPath)
346 
347  ! -- Release a particle from each point...
348  do nps = 1, this%nreleasepts
349  ic = this%rptnode(nps)
350  icu = this%dis%get_nodeuser(ic)
351  ! -- ...for each release time in the current time step
352  tsloop: do nts = 1, nreleasets
353  if (this%rlstimelist) then
354  trelease = this%releasetimes%times(nts)
355  tend = totimc + delt
356  if (trelease < totimc .or. trelease >= tend) cycle tsloop
357  else
358  trelease = totimc + this%offset * delt
359  end if
360 
361  np = this%nparticles + 1
362  this%nparticles = np
363 
364  ! -- Check release point is within the specified cell
365  ! and not above/below grid top/bottom respectively
366  x = this%rptx(nps)
367  y = this%rpty(nps)
368  z = this%rptz(nps)
369  call this%dis%get_polyverts(ic, polyverts)
370  if (.not. point_in_polygon(x, y, polyverts)) then
371  write (errmsg, '(a,g0,a,g0,a,i0)') &
372  'Error: release point (x=', x, ', y=', y, ') is not in cell ', icu
373  call store_error(errmsg, terminate=.false.)
374  call store_error_unit(this%inunit, terminate=.true.)
375  end if
376  if (z > maxval(this%dis%top)) then
377  write (errmsg, '(a,g0,a,g0,a,i0)') &
378  'Error: release point (z=', z, ') is above grid top ', &
379  maxval(this%dis%top)
380  call store_error(errmsg, terminate=.false.)
381  call store_error_unit(this%inunit, terminate=.true.)
382  else if (z < minval(this%dis%bot)) then
383  write (errmsg, '(a,g0,a,g0,a,i0)') &
384  'Error: release point (z=', z, ') is below grid bottom ', &
385  minval(this%dis%bot)
386  call store_error(errmsg, terminate=.false.)
387  call store_error_unit(this%inunit, terminate=.true.)
388  end if
389 
390  ! -- Initialize particle and add it to particle store
391  ! -- Todo: branch depending on exchange PRP or a normal PRP.
392  ! if exchange PRP, particle identity properties should be
393  ! passed in (e.g. imdl, iprp, irpt, trelease, name).
394  ! if normal PRP, imdl and iprp should be set from pointers
395  ! provided to PRP by PRT model; irpt and trelease as below.
396  allocate (particle)
397  call create_particle(particle)
398  if (size(this%boundname) /= 0) then
399  particle%name = this%boundname(nps)
400  else
401  particle%name = ''
402  end if
403  particle%irpt = nps
404  particle%istopweaksink = this%istopweaksink
405  particle%istopzone = this%istopzone
406  particle%icu = icu
407  select type (dis => this%dis)
408  type is (distype)
409  call get_ijk(icu, dis%nrow, dis%ncol, dis%nlay, irow, icol, ilay)
410  type is (disvtype)
411  call get_jk(icu, dis%ncpl, dis%nlay, icpl, ilay)
412  end select
413  particle%ilay = ilay
414  particle%izone = this%rptzone(ic)
415  particle%istatus = 0
416  ! Handle inactive cells
417  if (this%ibound(ic) == 0) then
418  ! -- If drape option activated, release in highest active
419  ! cell vertically below release point.
420  if (this%idrape /= 0) &
421  call this%dis%highest_active(ic, this%ibound)
422  ! -- If returned cell is inactive, do not release particle
423  if (this%ibound(ic) == 0) &
424  particle%istatus = 8 ! permanently unreleased
425  end if
426  particle%x = x
427  particle%y = y
428  if (this%localz) then
429  top = this%fmi%dis%top(ic)
430  bot = this%fmi%dis%bot(ic)
431  hds = this%fmi%gwfhead(ic)
432  particle%z = bot + this%rptz(nps) * (hds - bot)
433  else
434  particle%z = this%rptz(nps)
435  end if
436  particle%trelease = trelease
437  ! Set stopping time to earlier of times specified by STOPTIME and STOPTRAVELTIME
438  if (this%stoptraveltime == huge(1d0)) then
439  particle%tstop = this%stoptime
440  else
441  particle%tstop = particle%trelease + this%stoptraveltime
442  if (this%stoptime < particle%tstop) particle%tstop = this%stoptime
443  end if
444  particle%ttrack = particle%trelease
445  particle%idomain(0) = 0
446  particle%iboundary(0) = 0
447  particle%idomain(1) = 0
448  particle%iboundary(1) = 0
449  particle%idomain(2) = ic
450  particle%iboundary(2) = 0
451  particle%idomain(3) = 0
452  particle%iboundary(3) = 0
453  call this%particles%load_from_particle(particle, np)
454 
455  ! -- Accumulate mass release from this point
456  this%rptmass(nps) = this%rptmass(nps) + done
457  end do tsloop
458  end do
Definition: Dis.f90:1
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
Structured grid discretization.
Definition: Dis.f90:23
Vertex grid discretization.
Definition: Disv.f90:24
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 192 of file prt-prp.f90.

193  ! -- dummy
194  class(PrtPrpType) :: this
195  integer(I4B), dimension(:), pointer, contiguous, optional :: nodelist
196  real(DP), dimension(:, :), pointer, contiguous, optional :: auxvar
197  ! -- local
198  integer(I4B) :: nps
199 
200  ! -- Allocate particle store, starting with the number
201  ! of release points (arrays resized if/when needed)
202  call create_particle_store(this%particles, this%nreleasepts, this%memoryPath)
203 
204  ! -- Allocate arrays
205  call mem_allocate(this%rptx, this%nreleasepts, 'RPTX', this%memoryPath)
206  call mem_allocate(this%rpty, this%nreleasepts, 'RPTY', this%memoryPath)
207  call mem_allocate(this%rptz, this%nreleasepts, 'RPTZ', this%memoryPath)
208  call mem_allocate(this%locz, this%nreleasepts, 'LOCZ', this%memoryPath)
209  call mem_allocate(this%rptmass, this%nreleasepts, 'RPTMASS', this%memoryPath)
210  call mem_allocate(this%rptnode, this%nreleasepts, 'RPTNODER', &
211  this%memoryPath)
212  call mem_allocate(this%rlskstp, 1, 'RLSKSTP', this%memoryPath)
213  call mem_allocate(this%rptname, lenboundname, this%nreleasepts, &
214  'RPTNAME', this%memoryPath)
215 
216  ! -- Initialize arrays
217  this%rlskstp(1) = 1 ! single release in first time step by default
218  do nps = 1, this%nreleasepts
219  this%rptmass(nps) = dzero
220  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 224 of file prt-prp.f90.

225  class(PrtPrpType) :: this
226 
227  ! -- Allocate release time selection
228  allocate (this%releasetimes)
229  call this%releasetimes%init()
230 
231  ! -- call standard BndType allocate scalars
232  call this%BndType%allocate_scalars()
233 
234  ! -- Allocate scalars for this type
235  call mem_allocate(this%rlsall, 'RLSALL', this%memoryPath)
236  call mem_allocate(this%rlsfirst, 'RLSFIRST', this%memoryPath)
237  call mem_allocate(this%rlstimelist, 'RELEASETIME', this%memoryPath)
238  call mem_allocate(this%localz, 'LOCALZ', this%memoryPath)
239  call mem_allocate(this%offset, 'OFFSET', this%memoryPath)
240  call mem_allocate(this%stoptime, 'STOPTIME', this%memoryPath)
241  call mem_allocate(this%stoptraveltime, 'STOPTRAVELTIME', this%memoryPath)
242  call mem_allocate(this%istopweaksink, 'ISTOPWEAKSINK', this%memoryPath)
243  call mem_allocate(this%istopzone, 'ISTOPZONE', this%memoryPath)
244  call mem_allocate(this%idrape, 'IDRAPE', this%memoryPath)
245  call mem_allocate(this%nreleasepts, 'NRELEASEPTS', this%memoryPath)
246  call mem_allocate(this%nparticles, 'NPART', this%memoryPath)
247  call mem_allocate(this%itrkout, 'ITRKOUT', this%memoryPath)
248  call mem_allocate(this%itrkhdr, 'ITRKHDR', this%memoryPath)
249  call mem_allocate(this%itrkcsv, 'ITRKCSV', this%memoryPath)
250  call mem_allocate(this%irlstls, 'IRLSTLS', this%memoryPath)
251 
252  ! -- Set values
253  this%rlsall = .false.
254  this%rlsfirst = .false.
255  this%rlstimelist = .false.
256  this%localz = .false.
257  this%offset = dzero
258  this%stoptime = huge(1d0)
259  this%stoptraveltime = huge(1d0)
260  this%istopweaksink = 0
261  this%istopzone = 0
262  this%idrape = 0
263  this%nreleasepts = 0
264  this%nparticles = 0
265  this%itrkout = 0
266  this%itrkhdr = 0
267  this%itrkcsv = 0
268  this%irlstls = 0

◆ prp_ar()

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

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

273  ! -- dummy variables
274  class(PrtPrpType), intent(inout) :: this
275  ! -- local variables
276  integer(I4B) :: n
277 
278  call this%obs%obs_ar()
279  call this%BndType%allocate_arrays()
280  if (this%inamedbound /= 0) then
281  do n = 1, this%nreleasepts
282  this%boundname(n) = this%rptname(n)
283  end do
284  end if
285  do n = 1, this%nreleasepts
286  this%nodelist(n) = this%rptnode(n)
287  end do
288  ! if (this%imover /= 0) then
289  ! allocate(this%pakmvrobj)
290  ! call this%pakmvrobj%ar(this%maxbound, this%maxbound, this%memoryPath)
291  ! endif

◆ 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]hnewtodo: mass concentration?
[in,out]flowjaflow between package and model
[in]imoverflag indicating if the mover package is active

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

615  ! -- modules
616  use tdismodule, only: delt
617  ! -- dummy variables
618  class(PrtPrpType) :: this
619  real(DP), dimension(:), intent(in) :: hnew !< todo: mass concentration?
620  real(DP), dimension(:), intent(inout) :: flowja !< flow between package and model
621  integer(I4B), intent(in) :: imover !< flag indicating if the mover package is active
622  ! -- local variables
623  integer(I4B) :: i
624  integer(I4B) :: node
625  integer(I4B) :: idiag
626  real(DP) :: rrate
627 
628  ! -- If no boundaries, skip flow calculations.
629  if (this%nbound <= 0) return
630 
631  ! -- Loop through each boundary calculating flow.
632  do i = 1, this%nbound
633  node = this%nodelist(i)
634  rrate = dzero
635  ! -- If cell is no-flow or constant-head, then ignore it.
636  ! todo: think about condition(s) under which to ignore cell
637  if (node > 0) then
638  idiag = this%dis%con%ia(node)
639  ! todo: think about condition(s) under which to ignore cell
640  ! -- Calculate the flow rate into the cell.
641  rrate = this%rptmass(i) * (done / delt) ! reciprocal of tstp length
642  flowja(idiag) = flowja(idiag) + rrate
643  end if
644 
645  ! -- Save simulated value to simvals array.
646  this%simvals(i) = rrate
647  end do

◆ 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)  mempath,
type(prtfmitype), pointer  fmi 
)

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

89  ! -- dummy
90  class(BndType), pointer :: packobj
91  integer(I4B), intent(in) :: id
92  integer(I4B), intent(in) :: ibcnum
93  integer(I4B), intent(in) :: inunit
94  integer(I4B), intent(in) :: iout
95  character(len=*), intent(in) :: namemodel
96  character(len=*), intent(in) :: pakname
97  character(len=*), intent(in) :: mempath
98  type(PrtFmiType), pointer :: fmi
99  ! -- local
100  type(PrtPrpType), pointer :: prpobj
101  ! -- formats
102  character(len=*), parameter :: fmtheader = &
103  "(1x, /1x, 'PRP -- PARTICLE RELEASE POINT PACKAGE', &
104  &' INPUT READ FROM MEMPATH: ', A, /)"
105 
106  ! -- allocate the object and assign values to object variables
107  allocate (prpobj)
108  packobj => prpobj
109 
110  ! -- create name and memory path
111  call packobj%set_names(ibcnum, namemodel, pakname, ftype, mempath)
112  prpobj%text = text
113 
114  ! -- allocate scalars
115  call prpobj%prp_allocate_scalars()
116 
117  ! -- initialize package
118  call packobj%pack_initialize()
119 
120  packobj%inunit = inunit
121  packobj%iout = iout
122  packobj%id = id
123  packobj%ibcnum = ibcnum
124  packobj%ncolbnd = 4
125  packobj%iscloc = 1
126 
127  ! -- store pointer to flow model interface
128  prpobj%fmi => fmi
129 
130  ! -- if prp is enabled, print a message identifying it
131  if (inunit > 0) write (iout, fmtheader) mempath
Here is the caller graph for this function:

◆ prp_da()

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

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

136  class(PrtPrpType) :: this
137 
138  ! -- deallocate parent
139  call this%BndType%bnd_da()
140 
141  ! -- deallocate scalars
142  call mem_deallocate(this%rlsall)
143  call mem_deallocate(this%rlsfirst)
144  call mem_deallocate(this%rlstimelist)
145  call mem_deallocate(this%localz)
146  call mem_deallocate(this%offset)
147  call mem_deallocate(this%stoptime)
148  call mem_deallocate(this%stoptraveltime)
149  call mem_deallocate(this%istopweaksink)
150  call mem_deallocate(this%istopzone)
151  call mem_deallocate(this%idrape)
152  call mem_deallocate(this%nreleasepts)
153  call mem_deallocate(this%nparticles)
154  call mem_deallocate(this%itrkout)
155  call mem_deallocate(this%itrkhdr)
156  call mem_deallocate(this%itrkcsv)
157  call mem_deallocate(this%irlstls)
158 
159  ! -- deallocate arrays
160  call mem_deallocate(this%rptx)
161  call mem_deallocate(this%rpty)
162  call mem_deallocate(this%rptz)
163  call mem_deallocate(this%locz)
164  call mem_deallocate(this%rptnode)
165  call mem_deallocate(this%rptmass)
166  call mem_deallocate(this%rlskstp)
167  call mem_deallocate(this%rptname, 'RPTNAME', this%memoryPath)
168 
169  ! -- deallocate particle store
170  call this%particles%destroy(this%memoryPath)
171  deallocate (this%particles)
172 
173  ! -- deallocate release time selection
174  call this%releasetimes%destroy()
175  deallocate (this%releasetimes)

◆ prp_df_obs()

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

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

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

◆ prp_obs_supported()

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

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

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

◆ prp_options()

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

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

695  use openspecmodule, only: access, form
696  use constantsmodule, only: maxcharlen, dzero
699  ! -- dummy
700  class(PrtPrpType), intent(inout) :: this
701  character(len=*), intent(inout) :: option
702  logical, intent(inout) :: found
703  ! -- locals
704  real(DP) :: dval
705  integer(I4B) :: i, ios, nlines
706  logical(LGP) :: success
707  character(len=MAXCHARLEN) :: fname
708  character(len=MAXCHARLEN) :: keyword
709  character(len=:), allocatable :: line
710  ! -- formats
711  character(len=*), parameter :: fmttrkbin = &
712  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO BINARY FILE: ', a, /4x, &
713  &'OPENED ON UNIT: ', I0)"
714  character(len=*), parameter :: fmttrkcsv = &
715  "(4x, 'PARTICLE TRACKS WILL BE SAVED TO CSV FILE: ', a, /4x, &
716  &'OPENED ON UNIT: ', I0)"
717 
718  select case (option)
719  case ('STOPTIME')
720  this%stoptime = this%parser%GetDouble()
721  found = .true.
722  case ('STOPTRAVELTIME')
723  this%stoptraveltime = this%parser%GetDouble()
724  found = .true.
725  case ('STOP_AT_WEAK_SINK')
726  this%istopweaksink = 1
727  found = .true.
728  case ('ISTOPZONE')
729  this%istopzone = this%parser%GetInteger()
730  found = .true.
731  case ('DRAPE')
732  this%idrape = 1
733  found = .true.
734  case ('RELEASE_TIMES')
735  rtloop: do
736  success = .false.
737  call this%parser%TryGetDouble(dval, success)
738  if (.not. success) exit rtloop
739  call this%releasetimes%expand()
740  this%releasetimes%times(size(this%releasetimes%times)) = dval
741  end do rtloop
742  if (.not. this%releasetimes%increasing()) then
743  errmsg = "RELEASE TIMES MUST STRICTLY INCREASE"
744  call store_error(errmsg)
745  call this%parser%StoreErrorUnit(terminate=.true.)
746  end if
747  this%rlstimelist = .true.
748  found = .true.
749  case ('RELEASE_TIMESFILE')
750  call this%parser%GetString(fname)
751  call openfile(this%irlstls, this%iout, fname, 'TLS')
752  nlines = 0
753  rtfloop: do
754  read (this%irlstls, '(A)', iostat=ios) line
755  if (ios /= 0) exit rtfloop
756  nlines = nlines + 1
757  end do rtfloop
758  call this%releasetimes%expand(nlines)
759  rewind(this%irlstls)
760  allocate (character(len=LINELENGTH) :: line)
761  do i = 1, nlines
762  read (this%irlstls, '(A)') line
763  read (line, '(f30.0)') dval
764  this%releasetimes%times(i) = dval
765  end do
766  if (.not. this%releasetimes%increasing()) then
767  errmsg = "RELEASE TIMES MUST STRICTLY INCREASE"
768  call store_error(errmsg)
769  call this%parser%StoreErrorUnit(terminate=.true.)
770  end if
771  this%rlstimelist = .true.
772  found = .true.
773  case ('TRACK')
774  call this%parser%GetStringCaps(keyword)
775  if (keyword == 'FILEOUT') then
776  ! parse filename
777  call this%parser%GetString(fname)
778  ! open binary output file
779  this%itrkout = getunit()
780  call openfile(this%itrkout, this%iout, fname, 'DATA(BINARY)', &
781  form, access, filstat_opt='REPLACE', &
782  mode_opt=mnormal)
783  write (this%iout, fmttrkbin) trim(adjustl(fname)), this%itrkout
784  ! open and write ascii header spec file
785  this%itrkhdr = getunit()
786  fname = trim(fname)//'.hdr'
787  call openfile(this%itrkhdr, this%iout, fname, 'CSV', &
788  filstat_opt='REPLACE', mode_opt=mnormal)
789  write (this%itrkhdr, '(a,/,a)') trackheader, trackdtypes
790  else
791  call store_error('OPTIONAL TRACK KEYWORD MUST BE '// &
792  'FOLLOWED BY FILEOUT')
793  end if
794  found = .true.
795  case ('TRACKCSV')
796  call this%parser%GetStringCaps(keyword)
797  if (keyword == 'FILEOUT') then
798  ! parse filename
799  call this%parser%GetString(fname)
800  ! open CSV output file and write headers
801  this%itrkcsv = getunit()
802  call openfile(this%itrkcsv, this%iout, fname, 'CSV', &
803  filstat_opt='REPLACE')
804  write (this%iout, fmttrkcsv) trim(adjustl(fname)), this%itrkcsv
805  write (this%itrkcsv, '(a)') trackheader
806  else
807  call store_error('OPTIONAL TRACKCSV KEYWORD MUST BE &
808  &FOLLOWED BY FILEOUT')
809  end if
810  found = .true.
811  case ('LOCAL_Z')
812  this%localz = .true.
813  found = .true.
814  case default
815  found = .false.
816  end select
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter maxcharlen
maximum length of char string
Definition: Constants.f90:46
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 trackheader
Definition: TrackData.f90:56
character(len= *), parameter, public trackdtypes
Definition: TrackData.f90:61
Here is the call graph for this function:

◆ prp_read_dimensions()

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

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

961  ! -- dummy
962  class(PrtPrpType), intent(inout) :: this
963  ! -- local
964  character(len=LINELENGTH) :: errmsg, keyword
965  integer(I4B) :: ierr
966  logical :: isfound, endOfBlock
967 
968  ! -- get dimension block
969  call this%parser%GetBlock('DIMENSIONS', isfound, ierr, &
970  supportopenclose=.true.)
971 
972  ! -- parse dimension block if detected
973  if (isfound) then
974  write (this%iout, '(1x,a)') 'PROCESSING PARTICLE INPUT DIMENSIONS'
975  do
976  call this%parser%GetNextLine(endofblock)
977  if (endofblock) exit
978  call this%parser%GetStringCaps(keyword)
979  select case (keyword)
980  case ('NRELEASEPTS')
981  this%nreleasepts = this%parser%GetInteger()
982  case default
983  write (errmsg, &
984  '(4x,a,a)') '****ERROR. UNKNOWN PARTICLE INPUT DIMENSION: ', &
985  trim(keyword)
986  call store_error(errmsg)
987  call this%parser%StoreErrorUnit()
988  end select
989  end do
990  write (this%iout, '(1x,a)') 'END OF PARTICLE INPUT DIMENSIONS'
991  else
992  call store_error('ERROR. REQUIRED DIMENSIONS BLOCK NOT FOUND.')
993  end if
994 
995  ! -- set maxbound and nbound to nreleasepts
996  this%maxbound = this%nreleasepts
997  this%nbound = this%nreleasepts
998 
999  ! -- allocate arrays for prp package
1000  call this%prp_allocate_arrays()
1001 
1002  ! -- read packagedata
1003  call this%prp_read_packagedata()
Here is the call graph for this function:

◆ prp_read_packagedata()

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

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

821  ! -- dummy
822  class(PrtPrpType), intent(inout) :: this
823  ! -- local
824  character(len=LINELENGTH) :: cellid
825  character(len=LENBOUNDNAME) :: bndName, bndNameTemp
826  character(len=9) :: cno
827  logical :: isfound
828  logical :: endOfBlock
829  integer(I4B) :: ival
830  integer(I4B) :: n
831  integer(I4B) :: ierr
832  character(len=LENBOUNDNAME), dimension(:), allocatable :: nametxt
833  integer(I4B), dimension(:), allocatable :: nboundchk
834  integer(I4B), dimension(:), allocatable :: noder
835  real(DP), dimension(:), allocatable :: x
836  real(DP), dimension(:), allocatable :: y
837  real(DP), dimension(:), allocatable :: z
838  real(DP), dimension(:), allocatable :: tstop
839  ! -- format
840  character(len=*), parameter :: fmttend = &
841  "('end time (', G0, ') must be greater than or equal to the &
842  &begin time (', G0, ').')"
843 
844  ! -- allocate and initialize temporary variables
845  allocate (noder(this%nreleasepts))
846  allocate (x(this%nreleasepts))
847  allocate (y(this%nreleasepts))
848  allocate (z(this%nreleasepts))
849  allocate (tstop(this%nreleasepts))
850  allocate (nametxt(this%nreleasepts))
851  allocate (nboundchk(this%nreleasepts))
852 
853  ! -- initialize temporary variables
854  do n = 1, this%nreleasepts
855  nboundchk(n) = 0
856  end do
857 
858  ! -- read particle release point data
859  ! -- get particle release points block
860  call this%parser%GetBlock('PACKAGEDATA', isfound, ierr, &
861  supportopenclose=.true.)
862 
863  ! -- parse block if detected
864  if (isfound) then
865  write (this%iout, '(/1x,a)') 'PROCESSING '//trim(adjustl(this%packName)) &
866  //' PACKAGEDATA'
867  do
868  call this%parser%GetNextLine(endofblock)
869  if (endofblock) exit
870  ival = this%parser%GetInteger()
871  n = ival
872 
873  if (n < 1 .or. n > this%nreleasepts) then
874  write (errmsg, '(a,1x,i0,a)') &
875  'Release point number must be greater than 0 and less than ', &
876  'or equal to', this%nreleasepts, '.'
877  call store_error(errmsg)
878  cycle
879  end if
880 
881  ! -- increment nboundchk
882  nboundchk(n) = nboundchk(n) + 1
883 
884  ! -- node number
885  call this%parser%GetCellid(this%dis%ndim, cellid)
886  noder(n) = this%dis%noder_from_cellid(cellid, this%inunit, this%iout)
887 
888  ! -- x, y, z coordinates
889  x(n) = this%parser%GetDouble()
890  y(n) = this%parser%GetDouble()
891  z(n) = this%parser%GetDouble()
892 
893  if (this%localz .and. (z(n) < 0 .or. z(n) > 1)) then
894  call store_error('Local z coordinate must fall in the interval [0, 1]')
895  cycle
896  end if
897 
898  ! -- set default boundname
899  write (cno, '(i9.9)') n
900  bndname = 'PRP'//cno
901 
902  ! -- read boundnames from file, if provided
903  if (this%inamedbound /= 0) then
904  call this%parser%GetStringCaps(bndnametemp)
905  if (bndnametemp /= '') &
906  bndname = bndnametemp
907  else
908  bndname = ''
909  end if
910 
911  ! -- store temp boundnames
912  nametxt(n) = bndname
913  end do
914 
915  write (this%iout, '(1x,a)') &
916  'END OF '//trim(adjustl(this%packName))//' PACKAGEDATA'
917 
918  ! -- check for duplicate or missing particle release points
919  do n = 1, this%nreleasepts
920  if (nboundchk(n) == 0) then
921  write (errmsg, '(a,a,1x,i0,a)') 'No data specified for particle ', &
922  'release point', n, '.'
923  call store_error(errmsg)
924  else if (nboundchk(n) > 1) then
925  write (errmsg, '(a,1x,i0,1x,a,1x,i0,1x,a)') &
926  'Data for particle release point', n, 'specified', nboundchk(n), &
927  'times.'
928  call store_error(errmsg)
929  end if
930  end do
931  else
932  call store_error('Required packagedata block not found.')
933  end if
934 
935  ! -- terminate if any errors were detected
936  if (count_errors() > 0) then
937  call this%parser%StoreErrorUnit()
938  end if
939 
940  ! -- fill particle release point data with data stored in temporary local arrays
941  do n = 1, this%nreleasepts
942  this%rptnode(n) = noder(n)
943  this%rptx(n) = x(n)
944  this%rpty(n) = y(n)
945  this%rptz(n) = z(n)
946  this%rptname(n) = nametxt(n)
947  end do
948 
949  ! -- deallocate local storage
950  deallocate (noder)
951  deallocate (x)
952  deallocate (y)
953  deallocate (z)
954  deallocate (tstop)
955  deallocate (nametxt)
956  deallocate (nboundchk)
Here is the call graph for this function:

◆ prp_rp()

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

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

463  ! -- modules
464  use tdismodule, only: kper, nper, nstp
465  use inputoutputmodule, only: urword
466  ! -- dummy variables
467  class(PrtPrpType), intent(inout) :: this
468  ! -- local variables
469  integer(I4B) :: ierr
470  integer(I4B) :: n, i
471  integer(I4B) :: lloc, istart, istop, ival
472  real(DP) :: dval
473  logical(LGP) :: isfound
474  logical(LGP) :: endOfBlock
475  logical(LGP) :: use_last
476  logical(LGP) :: noperiodblocks
477  character(len=LINELENGTH) :: keyword
478  character(len=:), allocatable :: line
479  ! -- formats
480  character(len=*), parameter :: fmtblkerr = &
481  "('Looking for BEGIN PERIOD iper. &
482  &Found ', a, ' instead.')"
483  character(len=*), parameter :: fmt_steps = &
484  "(6x,'TIME STEP(S) ',50(I0,' '))" ! kluge 50 (similar to STEPS in OC)?
485  character(len=*), parameter :: fmt_freq = &
486  "(6x,'EVERY ',I0,' TIME STEP(S)')"
487  character(len=*), parameter :: fmt_fracs = &
488  "(6x,50(f10.3,' '))"
489 
490  ! -- Set ionper to the stress period number for which a new block of data
491  ! will be read.
492  if (this%inunit == 0) return
493 
494  ! -- get stress period data
495  noperiodblocks = .false.
496  if (this%ionper < kper) then
497  ! -- get period block
498  call this%parser%GetBlock('PERIOD', isfound, ierr, &
499  supportopenclose=.true., &
500  blockrequired=.false.)
501  if (isfound) then
502  ! -- read ionper and check for increasing period numbers
503  call this%read_check_ionper()
504  else
505  ! -- PERIOD block not found
506  if (ierr < 0) then
507  if (kper == 1) then
508  ! -- End of file found; no period data for the simulation.
509  noperiodblocks = .true.
510  else
511  ! -- End of file found; no more period data.
512  this%ionper = nper + 1
513  end if
514  else
515  ! -- Found invalid block
516  call this%parser%GetCurrentLine(line)
517  write (errmsg, fmtblkerr) adjustl(trim(line))
518  call store_error(errmsg, terminate=.true.)
519  end if
520  end if
521  end if
522 
523  ! -- If no period data for the simulation default to single
524  ! release at beginning of first period's first time step.
525  ! Otherwise read release timing settings from the period
526  ! data block of the package input file.
527  if (noperiodblocks) then
528  if (kper == 1) then
529  call mem_reallocate(this%rlskstp, 1, &
530  "RLSKSTP", this%memoryPath)
531  this%rlsfirst = .true.
532  use_last = .false.
533  end if
534  ! -- If the current stress period matches the
535  ! block we are reading continue parsing it
536  else if (this%ionper == kper) then
537  use_last = .false.
538  recordloop: do
539  call this%parser%GetNextLine(endofblock)
540  if (endofblock) exit
541  call this%parser%GetStringCaps(keyword)
542  select case (keyword)
543  case ('ALL')
544  this%rlsall = .true.
545  case ('STEPS')
546  call mem_reallocate(this%rlskstp, 0, &
547  "RLSKSTP", this%memoryPath)
548  call this%parser%GetRemainingLine(line)
549  lloc = 1
550  stepslistsearch: do
551  call urword(line, lloc, istart, istop, 2, ival, dval, -1, 0)
552  if (ival > 0) then
553  n = size(this%rlskstp)
554  call mem_reallocate(this%rlskstp, n + 1, &
555  'RLSKSTP', this%memoryPath)
556  this%rlskstp(n + 1) = ival
557  cycle stepslistsearch
558  end if
559  exit stepslistsearch
560  end do stepslistsearch
561  case ('FIRST')
562  this%rlsfirst = .true.
563  case ('FREQUENCY')
564  ival = this%parser%GetInteger()
565  if (ival < 0) then
566  errmsg = "FREQUENCY must be non-negative"
567  call store_error(errmsg)
568  call this%parser%StoreErrorUnit(terminate=.true.)
569  end if
570  do i = 1, nstp(this%ionper)
571  if (mod(i, ival) == 0) then
572  n = size(this%rlskstp)
573  call mem_reallocate(this%rlskstp, n + 1, &
574  'RLSKSTP', this%memoryPath)
575  this%rlskstp(n + 1) = i
576  end if
577  end do
578  case ('FRACTION')
579  dval = this%parser%GetDouble()
580  this%offset = dval
581  case default
582  write (errmsg, '(2a)') &
583  'Looking for ALL, STEPS, FIRST, FREQUENCY, or FRACTION. Found: ', &
584  trim(adjustl(keyword))
585  call store_error(errmsg, terminate=.true.)
586  end select
587  end do recordloop
588  else
589  ! -- else repeat period settings
590  use_last = .true.
591  end if
592 
593  ! -- write settings to list file
594  if (.not. any(this%rlskstp > 0)) then
595  write (this%iout, "(1x,/1x,a)") 'NO PARTICLE RELEASES IN THIS STRESS '// &
596  'PERIOD'
597  else if (use_last) then
598  write (this%iout, "(1x,/1x,a)") 'REUSING PARTICLE RELEASE SETTINGS '// &
599  'FROM LAST STRESS PERIOD'
600  else
601  ! -- write particle release setting
602  write (this%iout, "(1x,/1x,a)", advance='no') 'PARTICLE RELEASE:'
603  if (any(this%rlskstp > 0)) then
604  n = size(this%rlskstp)
605  if (n > 0) write (this%iout, fmt_steps, advance='no') this%rlskstp
606  end if
607  write (this%iout, "(1x,a)", advance='no') 'AT OFFSET'
608  write (this%iout, fmt_fracs) (/this%offset/)
609  write (this%iout, '(A)')
610  end if
integer(i4b), dimension(:), pointer, public, contiguous nstp
number of time steps in each stress period
Definition: tdis.f90:39
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(trackfilecontroltype), pointer  trackfilectl 
)
private

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

180  ! -- dummy variables
181  class(PrtPrpType) :: this
182  integer(I4B), dimension(:), pointer, contiguous :: ibound
183  integer(I4B), dimension(:), pointer, contiguous :: izone
184  type(TrackFileControlType), pointer :: trackfilectl
185 
186  this%ibound => ibound
187  this%rptzone => izone
188  this%trackfilectl => trackfilectl

Variable Documentation

◆ ftype

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

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

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

◆ text

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

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

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