MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
particleeventsmodule Module Reference

Data Types

type  particleeventconsumertype
 
type  particleeventdispatchertype
 
interface  handle_event
 

Functions/Subroutines

subroutine subscribe (this, consumer)
 Subscribe a consumer to the dispatcher. More...
 
subroutine dispatch (this, particle, event)
 Dispatch an event. More...
 
subroutine destroy (this)
 Destroy the dispatcher. More...
 

Function/Subroutine Documentation

◆ destroy()

subroutine particleeventsmodule::destroy ( class(particleeventdispatchertype), intent(inout)  this)

Definition at line 85 of file ParticleEvents.f90.

86  class(ParticleEventDispatcherType), intent(inout) :: this
87  call this%consumers%Clear()

◆ dispatch()

subroutine particleeventsmodule::dispatch ( class(particleeventdispatchertype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
class(particleeventtype), intent(inout), pointer  event 
)
private

Definition at line 43 of file ParticleEvents.f90.

44  use tdismodule, only: kper, kstp
45  ! dummy
46  class(ParticleEventDispatcherType), intent(inout) :: this
47  type(ParticleType), pointer, intent(inout) :: particle
48  class(ParticleEventType), pointer, intent(inout) :: event
49  ! local
50  integer(I4B) :: i
51  real(DP) :: x, y, z
52  class(*), pointer :: p
53 
54  ! Convert to model coordinates if we need to
55  x = particle%x
56  y = particle%y
57  z = particle%z
58  call particle%get_model_coords(x, y, z)
59 
60  event%kper = kper
61  event%kstp = kstp
62  event%imdl = particle%imdl
63  event%iprp = particle%iprp
64  event%irpt = particle%irpt
65  event%ilay = particle%ilay
66  event%icu = particle%icu
67  event%izone = particle%izone
68  event%trelease = particle%trelease
69  event%ttrack = particle%ttrack
70  event%x = x
71  event%y = y
72  event%z = z
73  event%istatus = particle%istatus
74 
75  do i = 1, this%consumers%Count()
76  p => this%consumers%GetItem(i)
77  select type (consumer => p)
78  class is (particleeventconsumertype)
79  call consumer%handle_event(particle, event)
80  end select
81  end do
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

◆ subscribe()

subroutine particleeventsmodule::subscribe ( class(particleeventdispatchertype), intent(inout)  this,
class(particleeventconsumertype), intent(inout), target  consumer 
)
private

Definition at line 34 of file ParticleEvents.f90.

35  class(ParticleEventDispatcherType), intent(inout) :: this
36  class(ParticleEventConsumerType), target, intent(inout) :: consumer
37  class(*), pointer :: p
38  p => consumer
39  call this%consumers%Add(p)