MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
methodmodule Module Reference

Particle tracking strategies.

Data Types

type  methodtype
 Base type for particle tracking methods. More...
 
interface  apply
 
interface  deallocate
 

Functions/Subroutines

subroutine init (this, fmi, cell, subcell, trackctl, tracktimes, izone, flowja, porosity, retfactor)
 
recursive subroutine track (this, particle, level, tmax)
 Track the particle over domains of the given. More...
 
subroutine try_pass (this, particle, nextlevel, advancing)
 Try passing the particle to the next subdomain. More...
 
subroutine load (this, particle, next_level, submethod)
 Load the subdomain tracking method (submethod). More...
 
subroutine pass (this, particle)
 Pass the particle to the next subdomain. More...
 
subroutine save (this, particle, reason)
 Save the particle's state to output files. More...
 
subroutine check (this, particle, cell_defn)
 Check reporting/terminating conditions before tracking. More...
 

Function/Subroutine Documentation

◆ check()

subroutine methodmodule::check ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
type(celldefntype), intent(inout), pointer  cell_defn 
)

Check a number of conditions determining whether to continue tracking the particle or terminate it, as well as whether to record any output data as per selected reporting conditions.

Definition at line 194 of file Method.f90.

195  ! modules
196  use tdismodule, only: endofsimulation, totim
197  ! dummy
198  class(MethodType), intent(inout) :: this
199  type(ParticleType), pointer, intent(inout) :: particle
200  type(CellDefnType), pointer, intent(inout) :: cell_defn
201  ! local
202  logical(LGP) :: dry_cell, dry_particle, no_exit_face, stop_zone, weak_sink
203 
204  dry_cell = this%fmi%ibdgwfsat0(cell_defn%icell) == 0
205  dry_particle = particle%z > cell_defn%top
206  no_exit_face = cell_defn%inoexitface > 0
207  stop_zone = cell_defn%izone > 0 .and. particle%istopzone == cell_defn%izone
208  weak_sink = cell_defn%iweaksink > 0
209 
210  particle%izone = cell_defn%izone
211  if (stop_zone) then
212  particle%advancing = .false.
213  particle%istatus = 6
214  call this%save(particle, reason=3)
215  return
216  end if
217 
218  if (no_exit_face .and. .not. dry_cell) then
219  particle%advancing = .false.
220  particle%istatus = 5
221  call this%save(particle, reason=3)
222  return
223  end if
224 
225  if (weak_sink) then
226  if (particle%istopweaksink > 0) then
227  particle%advancing = .false.
228  particle%istatus = 3
229  call this%save(particle, reason=3)
230  return
231  else
232  call this%save(particle, reason=4)
233  end if
234  end if
235 
236  if (dry_cell) then
237  if (particle%idrymeth == 0) then
238  no_exit_face = .false.
239  else if (particle%idrymeth == 1) then
240  ! stop
241  particle%advancing = .false.
242  particle%istatus = 7
243  call this%save(particle, reason=3)
244  return
245  else if (particle%idrymeth == 2) then
246  ! stay
247  no_exit_face = .false.
248  particle%advancing = .false.
249  particle%ttrack = totim
250  ! terminate if last period/step
251  if (endofsimulation) then
252  particle%istatus = 5
253  call this%save(particle, reason=3)
254  return
255  end if
256  call this%save(particle, reason=2)
257  end if
258  else if (dry_particle .and. this%name /= "passtobottom") then
259  ! dry particle
260  if (particle%idrymeth == 0) then
261  ! drop to water table
262  particle%z = cell_defn%top
263  call this%save(particle, reason=1)
264  else if (particle%idrymeth == 1) then
265  ! terminate
266  particle%advancing = .false.
267  particle%istatus = 7
268  call this%save(particle, reason=3)
269  return
270  else if (particle%idrymeth == 2) then
271  ! stay
272  no_exit_face = .false.
273  particle%advancing = .false.
274  return
275  end if
276  end if
277 
278  if (no_exit_face) then
279  particle%advancing = .false.
280  particle%istatus = 5
281  call this%save(particle, reason=3)
282  return
283  end if
284 
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32

◆ init()

subroutine methodmodule::init ( class(methodtype), intent(inout)  this,
type(prtfmitype), intent(in), optional, pointer  fmi,
class(celltype), intent(in), optional, pointer  cell,
class(subcelltype), intent(in), optional, pointer  subcell,
type(trackcontroltype), intent(in), optional, pointer  trackctl,
type(timeselecttype), intent(in), optional, pointer  tracktimes,
integer(i4b), dimension(:), intent(in), optional, pointer  izone,
real(dp), dimension(:), intent(in), optional, pointer  flowja,
real(dp), dimension(:), intent(in), optional, pointer  porosity,
real(dp), dimension(:), intent(in), optional, pointer  retfactor 
)
private

Definition at line 75 of file Method.f90.

77  class(MethodType), intent(inout) :: this
78  type(PrtFmiType), intent(in), pointer, optional :: fmi
79  class(CellType), intent(in), pointer, optional :: cell
80  class(SubcellType), intent(in), pointer, optional :: subcell
81  type(TrackControlType), intent(in), pointer, optional :: trackctl
82  type(TimeSelectType), intent(in), pointer, optional :: tracktimes
83  integer(I4B), intent(in), pointer, optional :: izone(:)
84  real(DP), intent(in), pointer, optional :: flowja(:)
85  real(DP), intent(in), pointer, optional :: porosity(:)
86  real(DP), intent(in), pointer, optional :: retfactor(:)
87 
88  if (present(fmi)) this%fmi => fmi
89  if (present(cell)) this%cell => cell
90  if (present(subcell)) this%subcell => subcell
91  if (present(trackctl)) this%trackctl => trackctl
92  if (present(tracktimes)) this%tracktimes => tracktimes
93  if (present(izone)) this%izone => izone
94  if (present(flowja)) this%flowja => flowja
95  if (present(porosity)) this%porosity => porosity
96  if (present(retfactor)) this%retfactor => retfactor

◆ load()

subroutine methodmodule::load ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer, intent(in)  next_level,
class(methodtype), intent(inout), pointer  submethod 
)
private

Definition at line 143 of file Method.f90.

144  class(MethodType), intent(inout) :: this
145  type(ParticleType), pointer, intent(inout) :: particle
146  integer, intent(in) :: next_level
147  class(MethodType), pointer, intent(inout) :: submethod
148  call pstop(1, "load must be overridden")
Here is the call graph for this function:

◆ pass()

subroutine methodmodule::pass ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle 
)
private

Definition at line 152 of file Method.f90.

153  class(MethodType), intent(inout) :: this
154  type(ParticleType), pointer, intent(inout) :: particle
155  call pstop(1, "pass must be overridden")
Here is the call graph for this function:

◆ save()

subroutine methodmodule::save ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b), intent(in)  reason 
)
private

Definition at line 159 of file Method.f90.

160  use tdismodule, only: kper, kstp, totimc
161  ! dummy
162  class(MethodType), intent(inout) :: this
163  type(ParticleType), pointer, intent(inout) :: particle
164  integer(I4B), intent(in) :: reason
165  ! local
166  integer(I4B) :: per, stp
167 
168  per = kper
169  stp = kstp
170 
171  ! If tracking time falls exactly on a boundary between time steps,
172  ! report the previous time step for this datum. This is to follow
173  ! MP7's behavior, and because the particle will have been tracked
174  ! up to this instant under the previous time step's conditions, so
175  ! the time step we're about to start shouldn't get "credit" for it.
176  if (particle%ttrack == totimc .and. (per > 1 .or. stp > 1)) then
177  if (stp > 1) then
178  stp = stp - 1
179  else if (per > 1) then
180  per = per - 1
181  stp = 1
182  end if
183  end if
184 
185  call this%trackctl%save(particle, kper=per, kstp=stp, reason=reason)
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
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ track()

recursive subroutine methodmodule::track ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b)  level,
real(dp), intent(in)  tmax 
)
private

Definition at line 101 of file Method.f90.

102  ! dummy
103  class(MethodType), intent(inout) :: this
104  type(ParticleType), pointer, intent(inout) :: particle
105  integer(I4B) :: level
106  real(DP), intent(in) :: tmax
107  ! local
108  logical(LGP) :: advancing
109  integer(I4B) :: nextlevel
110  class(methodType), pointer :: submethod
111 
112  ! Advance the particle over subdomains
113  advancing = .true.
114  nextlevel = level + 1
115  do while (advancing)
116  call this%load(particle, nextlevel, submethod)
117  call submethod%apply(particle, tmax)
118  call this%try_pass(particle, nextlevel, advancing)
119  end do

◆ try_pass()

subroutine methodmodule::try_pass ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b)  nextlevel,
logical(lgp)  advancing 
)
private

Definition at line 123 of file Method.f90.

124  class(MethodType), intent(inout) :: this
125  type(ParticleType), pointer, intent(inout) :: particle
126  integer(I4B) :: nextlevel
127  logical(LGP) :: advancing
128 
129  ! if the particle is done advancing, reset the domain boundary flag.
130  if (.not. particle%advancing) then
131  particle%iboundary = 0
132  advancing = .false.
133  else
134  ! otherwise pass the particle to the next subdomain.
135  ! if that leaves it on a boundary, stop advancing.
136  call this%pass(particle)
137  if (particle%iboundary(nextlevel - 1) .ne. 0) &
138  advancing = .false.
139  end if