MODFLOW 6  version 6.7.0.dev1
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...
 

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

76  class(MethodType), intent(inout) :: this
77  type(PrtFmiType), intent(in), pointer, optional :: fmi
78  class(CellType), intent(in), pointer, optional :: cell
79  class(SubcellType), intent(in), pointer, optional :: subcell
80  type(TrackControlType), intent(in), pointer, optional :: trackctl
81  type(TimeSelectType), intent(in), pointer, optional :: tracktimes
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(:)
86 
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

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

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

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

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

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

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