32 character(len=40),
pointer,
public :: name
33 logical(LGP),
public :: delegates
35 class(
celltype),
pointer,
public :: cell => null()
39 integer(I4B),
dimension(:),
pointer,
contiguous,
public :: izone => null()
40 real(dp),
dimension(:),
pointer,
contiguous,
public :: flowja => null()
41 real(dp),
dimension(:),
pointer,
contiguous,
public :: porosity => null()
42 real(dp),
dimension(:),
pointer,
contiguous,
public :: retfactor => null()
58 subroutine apply(this, particle, tmax)
64 real(DP),
intent(in) :: tmax
74 subroutine init(this, fmi, cell, subcell, trackctl, tracktimes, &
75 izone, flowja, porosity, retfactor)
77 type(
prtfmitype),
intent(in),
pointer,
optional :: fmi
78 class(
celltype),
intent(in),
pointer,
optional :: cell
79 class(
subcelltype),
intent(in),
pointer,
optional :: subcell
82 integer(I4B),
intent(in),
pointer,
optional :: izone(:)
83 real(DP),
intent(in),
pointer,
optional :: flowja(:)
84 real(DP),
intent(in),
pointer,
optional :: porosity(:)
85 real(DP),
intent(in),
pointer,
optional :: retfactor(:)
87 if (
present(fmi)) this%fmi => fmi
88 if (
present(cell)) this%cell => cell
89 if (
present(subcell)) this%subcell => subcell
90 if (
present(trackctl)) this%trackctl => trackctl
91 if (
present(tracktimes)) this%tracktimes => tracktimes
92 if (
present(izone)) this%izone => izone
93 if (
present(flowja)) this%flowja => flowja
94 if (
present(porosity)) this%porosity => porosity
95 if (
present(retfactor)) this%retfactor => retfactor
100 recursive subroutine track(this, particle, level, tmax)
104 integer(I4B) :: level
105 real(dp),
intent(in) :: tmax
107 logical(LGP) :: advancing
108 integer(I4B) :: nextlevel
113 nextlevel = level + 1
115 call this%load(particle, nextlevel, submethod)
116 call submethod%apply(particle, tmax)
117 call this%try_pass(particle, nextlevel, advancing)
122 subroutine try_pass(this, particle, nextlevel, advancing)
125 integer(I4B) :: nextlevel
126 logical(LGP) :: advancing
129 if (.not. particle%advancing)
then
130 particle%iboundary = 0
135 call this%pass(particle)
136 if (particle%iboundary(nextlevel - 1) .ne. 0) &
142 subroutine load(this, particle, next_level, submethod)
145 integer,
intent(in) :: next_level
146 class(
methodtype),
pointer,
intent(inout) :: submethod
147 call pstop(1,
"load must be overridden")
151 subroutine pass(this, particle)
154 call pstop(1,
"pass must be overridden")
158 subroutine save(this, particle, reason)
163 integer(I4B),
intent(in) :: reason
165 integer(I4B) :: per, stp
175 if (particle%ttrack ==
totimc .and. (per > 1 .or. stp > 1))
then
178 else if (per > 1)
then
184 call this%trackctl%save(particle,
kper=per,
kstp=stp, reason=reason)
This module contains simulation constants.
real(dp), parameter dzero
real constant zero
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
This module defines variable data types.
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Particle tracking strategies.
subroutine load(this, particle, next_level, submethod)
Load the subdomain tracking method (submethod).
subroutine save(this, particle, reason)
Save the particle's state to output files.
recursive subroutine track(this, particle, level, tmax)
Track the particle over domains of the given.
subroutine try_pass(this, particle, nextlevel, advancing)
Try passing the particle to the next subdomain.
subroutine pass(this, particle)
Pass the particle to the next subdomain.
real(dp), pointer, public totimc
simulation time at start of time step
integer(i4b), pointer, public kstp
current time step number
integer(i4b), pointer, public kper
current stress period number
Specify times for some event to occur.
Base grid cell definition.
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Base type for particle tracking methods.
Particle tracked by the PRT model.
Represents a series of instants at which some event should occur.
Manages particle track (i.e. pathline) files.