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 195 of file Method.f90.

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

145  class(MethodType), intent(inout) :: this
146  type(ParticleType), pointer, intent(inout) :: particle
147  integer, intent(in) :: next_level
148  class(MethodType), pointer, intent(inout) :: submethod
149  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 153 of file Method.f90.

154  class(MethodType), intent(inout) :: this
155  type(ParticleType), pointer, intent(inout) :: particle
156  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 160 of file Method.f90.

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