MODFLOW 6  version 6.6.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 particle through subdomains. 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 subdomain tracking method (submethod) More...
 
subroutine pass (this, particle)
 Pass a particle to the next subdomain, internal use only. More...
 
subroutine save (this, particle, reason)
 Save a particle's current state. More...
 
subroutine update (this, particle, cell_defn)
 Update particle state and check termination conditions. More...
 

Function/Subroutine Documentation

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

75  class(MethodType), intent(inout) :: this
76  type(PrtFmiType), intent(in), pointer, optional :: fmi
77  class(CellType), intent(in), pointer, optional :: cell
78  class(SubcellType), intent(in), pointer, optional :: subcell
79  type(TrackControlType), intent(in), pointer, optional :: trackctl
80  type(TimeSelectType), intent(in), pointer, optional :: tracktimes
81  integer(I4B), intent(in), pointer, optional :: izone(:)
82  real(DP), intent(in), pointer, optional :: flowja(:)
83  real(DP), intent(in), pointer, optional :: porosity(:)
84  real(DP), intent(in), pointer, optional :: retfactor(:)
85 
86  if (present(fmi)) this%fmi => fmi
87  if (present(cell)) this%cell => cell
88  if (present(subcell)) this%subcell => subcell
89  if (present(trackctl)) this%trackctl => trackctl
90  if (present(tracktimes)) this%tracktimes => tracktimes
91  if (present(izone)) this%izone => izone
92  if (present(flowja)) this%flowja => flowja
93  if (present(porosity)) this%porosity => porosity
94  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 140 of file Method.f90.

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

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

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

99  ! dummy
100  class(MethodType), intent(inout) :: this
101  type(ParticleType), pointer, intent(inout) :: particle
102  integer(I4B) :: level
103  real(DP), intent(in) :: tmax
104  ! local
105  logical(LGP) :: advancing
106  integer(I4B) :: nextlevel
107  class(methodType), pointer :: submethod
108 
109  advancing = .true.
110  nextlevel = level + 1
111  do while (advancing)
112  call this%load(particle, nextlevel, submethod)
113  call submethod%apply(particle, tmax)
114  call this%try_pass(particle, nextlevel, advancing)
115  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 119 of file Method.f90.

120  class(MethodType), intent(inout) :: this
121  type(ParticleType), pointer, intent(inout) :: particle
122  integer(I4B) :: nextlevel
123  logical(LGP) :: advancing
124 
125  ! tracking submethod marked tracking complete?
126  ! reset domain boundary flag and don't advance
127  if (.not. particle%advancing) then
128  particle%iboundary = 0
129  advancing = .false.
130  else
131  ! otherwise pass particle to next subdomain
132  ! and if it's on a boundary, stop advancing
133  call this%pass(particle)
134  if (particle%iboundary(nextlevel - 1) .ne. 0) &
135  advancing = .false.
136  end if

◆ update()

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

Update the particle's properties (e.g. advancing flag, zone number, status). If any termination conditions apply, the particle's status will be set to the appropriate termination value. If any reporting conditions apply, save particle state with the proper reason code.

Definition at line 194 of file Method.f90.

195  ! dummy
196  class(MethodType), intent(inout) :: this
197  type(ParticleType), pointer, intent(inout) :: particle
198  type(CellDefnType), pointer, intent(inout) :: cell_defn
199 
200  particle%izone = cell_defn%izone
201  if (cell_defn%izone .ne. 0) then
202  if (particle%istopzone .eq. cell_defn%izone) then
203  particle%advancing = .false.
204  particle%istatus = 6
205  call this%save(particle, reason=3) ! particle terminated
206  return
207  end if
208  end if
209  if (cell_defn%inoexitface .ne. 0) then
210  particle%advancing = .false.
211  particle%istatus = 5
212  call this%save(particle, reason=3) ! particle terminated
213  return
214  end if
215  if (cell_defn%iweaksink .ne. 0) then
216  if (particle%istopweaksink .ne. 0) then
217  particle%advancing = .false.
218  particle%istatus = 3
219  call this%save(particle, reason=3) ! particle terminated
220  return
221  else
222  call this%save(particle, reason=4) ! particle exited weak sink
223  return
224  end if
225  end if