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

Particle tracking strategies.

Data Types

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

Enumerations

enum  
 Tracking method level enumeration. More...
 

Functions/Subroutines

subroutine init (this, fmi, cell, subcell, events, tracktimes, izone, flowja, porosity, retfactor)
 
recursive subroutine track (this, particle, level, tmax)
 Track the particle over subdomains 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...
 
integer(i4b) function get_level (this)
 Get the tracking method's level. More...
 
subroutine pass (this, particle)
 Pass the particle to the next subdomain. More...
 
subroutine release (this, particle)
 Particle is released. More...
 
subroutine terminate (this, particle, status)
 Particle terminates. More...
 
subroutine timestep (this, particle)
 Time step ends. More...
 
subroutine weaksink (this, particle)
 Particle leaves a weak sink. More...
 
subroutine usertime (this, particle)
 User-defined tracking time occurs. More...
 

Variables

@, public level_model = 1
 
@, public level_feature = 2
 
@, public level_subfeature = 3
 

Enumeration Type Documentation

◆ anonymous enum

anonymous enum

Tracking levels: 1: model, 2: grid feature, 3: grid subfeature. A tracking level identifies the domain through which a tracking method is responsible for moving a particle. Methods operate on a particular level and delegate to submethods for levels higher than (i.e. below the scope of) their own.

Definition at line 35 of file Method.f90.

Function/Subroutine Documentation

◆ get_level()

integer(i4b) function methodmodule::get_level ( class(methodtype), intent(in)  this)
private

Definition at line 191 of file Method.f90.

192  class(MethodType), intent(in) :: this
193  integer(I4B) :: level
194  level = -1 ! suppress compiler warning
195  call pstop(1, "get_level must be overridden")
Here is the call graph for this function:

◆ 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(particleeventdispatchertype), intent(in), optional, pointer  events,
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 114 of file Method.f90.

116  class(MethodType), intent(inout) :: this
117  type(PrtFmiType), intent(in), pointer, optional :: fmi
118  class(CellType), intent(in), pointer, optional :: cell
119  class(SubcellType), intent(in), pointer, optional :: subcell
120  type(ParticleEventDispatcherType), intent(in), pointer, optional :: events
121  type(TimeSelectType), intent(in), pointer, optional :: tracktimes
122  integer(I4B), intent(in), pointer, optional :: izone(:)
123  real(DP), intent(in), pointer, optional :: flowja(:)
124  real(DP), intent(in), pointer, optional :: porosity(:)
125  real(DP), intent(in), pointer, optional :: retfactor(:)
126 
127  if (present(fmi)) this%fmi => fmi
128  if (present(cell)) this%cell => cell
129  if (present(subcell)) this%subcell => subcell
130  if (present(events)) this%events => events
131  if (present(tracktimes)) this%tracktimes => tracktimes
132  if (present(izone)) this%izone => izone
133  if (present(flowja)) this%flowja => flowja
134  if (present(porosity)) this%porosity => porosity
135  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 182 of file Method.f90.

183  class(MethodType), intent(inout) :: this
184  type(ParticleType), pointer, intent(inout) :: particle
185  integer, intent(in) :: next_level
186  class(MethodType), pointer, intent(inout) :: submethod
187  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 199 of file Method.f90.

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

◆ release()

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

Definition at line 206 of file Method.f90.

207  class(MethodType), intent(inout) :: this
208  type(ParticleType), pointer, intent(inout) :: particle
209  class(ParticleEventType), pointer :: event
210  allocate (releaseeventtype :: event)
211  call this%events%dispatch(particle, event)
212  deallocate (event)

◆ terminate()

subroutine methodmodule::terminate ( class(methodtype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
integer(i4b), intent(in), optional  status 
)
private

Definition at line 216 of file Method.f90.

217  class(MethodType), intent(inout) :: this
218  type(ParticleType), pointer, intent(inout) :: particle
219  integer(I4B), intent(in), optional :: status
220  class(ParticleEventType), pointer :: event
221  particle%advancing = .false.
222  if (present(status)) particle%istatus = status
223  allocate (terminationeventtype :: event)
224  call this%events%dispatch(particle, event)
225  deallocate (event)

◆ timestep()

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

Definition at line 229 of file Method.f90.

230  class(MethodType), intent(inout) :: this
231  type(ParticleType), pointer, intent(inout) :: particle
232  class(ParticleEventType), pointer :: event
233  allocate (timestepeventtype :: event)
234  call this%events%dispatch(particle, event)
235  deallocate (event)

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

141  ! dummy
142  class(MethodType), intent(inout) :: this
143  type(ParticleType), pointer, intent(inout) :: particle
144  integer(I4B) :: level
145  real(DP), intent(in) :: tmax
146  ! local
147  logical(LGP) :: advancing
148  integer(I4B) :: nextlevel
149  class(methodType), pointer :: submethod
150 
151  advancing = .true.
152  nextlevel = level + 1
153  do while (advancing)
154  call this%load(particle, nextlevel, submethod)
155  call submethod%apply(particle, tmax)
156  call this%try_pass(particle, nextlevel, advancing)
157  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 161 of file Method.f90.

162  class(MethodType), intent(inout) :: this
163  type(ParticleType), pointer, intent(inout) :: particle
164  integer(I4B) :: nextlevel
165  logical(LGP) :: advancing
166 
167  if (particle%advancing) then
168  ! if still advancing, pass to the next subdomain.
169  ! if that puts us on a boundary, then we're done.
170  call this%pass(particle)
171  if (particle%iboundary(nextlevel - 1) .ne. 0) &
172  advancing = .false.
173  else
174  ! otherwise we're already done so
175  ! reset the domain boundary value.
176  advancing = .false.
177  particle%iboundary = 0
178  end if

◆ usertime()

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

Definition at line 249 of file Method.f90.

250  class(MethodType), intent(inout) :: this
251  type(ParticleType), pointer, intent(inout) :: particle
252  class(ParticleEventType), pointer :: event
253  allocate (usertimeeventtype :: event)
254  call this%events%dispatch(particle, event)
255  deallocate (event)

◆ weaksink()

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

Definition at line 239 of file Method.f90.

240  class(MethodType), intent(inout) :: this
241  type(ParticleType), pointer, intent(inout) :: particle
242  class(ParticleEventType), pointer :: event
243  allocate (weaksinkeventtype :: event)
244  call this%events%dispatch(particle, event)
245  deallocate (event)

Variable Documentation

◆ level_feature

@, public methodmodule::level_feature = 2

Definition at line 37 of file Method.f90.

37  enumerator :: LEVEL_FEATURE = 2

◆ level_model

@, public methodmodule::level_model = 1

Definition at line 36 of file Method.f90.

36  enumerator :: LEVEL_MODEL = 1

◆ level_subfeature

@, public methodmodule::level_subfeature = 3

Definition at line 38 of file Method.f90.

38  enumerator :: LEVEL_SUBFEATURE = 3