MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
particleeventsmodule Module Reference

Data Types

type  particleeventsubscriptiontype
 Subscription to particle events: a procedure to handle the event with an unlimited pointer for storing arbitrary context the handling procedure may reference. More...
 
type  particleeventdispatchertype
 Dispatcher for particle events. Consumers subscribe handlers to the dispatcher. Events may be dispatched, with the first handler to handle the event stopping propagation, or broadcast, with all handlers receiving the event. More...
 
interface  handle_event
 Event handler interface. Handlers may signal to the dispatching caller whether they have handled the event, but the signal is ignored for broadcasts. More...
 

Functions/Subroutines

subroutine subscribe (this, handler, context)
 Add a subscription to the dispatcher. More...
 
subroutine prep_event (this, particle, event)
 Prepare an event for dispatching, loading it with the current state of the particle. For internal use only. More...
 
subroutine dispatch (this, particle, event)
 Dispatch an event for handling. The first subscriber to handle the event stops propagation. More...
 
subroutine broadcast (this, particle, event)
 Broadcast an event to all subscribers so all receive the event and a chance to handle it. More...
 
subroutine destroy (this)
 Destroy the dispatcher. More...
 

Function/Subroutine Documentation

◆ broadcast()

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

Definition at line 123 of file ParticleEvents.f90.

124  class(ParticleEventDispatcherType), intent(inout) :: this
125  type(ParticleType), pointer, intent(inout) :: particle
126  class(ParticleEventType), pointer, intent(in) :: event
127  ! local
128  logical(LGP) :: handled
129  integer(I4B) :: i
130  class(*), pointer :: p
131 
132  call this%prep_event(particle, event)
133 
134  do i = 1, this%subscriptions%Count()
135  p => this%subscriptions%GetItem(i)
136  select type (subscription => p)
137  type is (particleeventsubscriptiontype)
138  handled = subscription%handler(subscription%context, particle, event)
139  end select
140  end do

◆ destroy()

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

Definition at line 144 of file ParticleEvents.f90.

145  class(ParticleEventDispatcherType), intent(inout) :: this
146  call this%subscriptions%Clear(destroy=.true.)

◆ dispatch()

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

Definition at line 97 of file ParticleEvents.f90.

98  class(ParticleEventDispatcherType), intent(inout) :: this
99  type(ParticleType), pointer, intent(inout) :: particle
100  class(ParticleEventType), pointer, intent(in) :: event
101  ! local
102  logical(LGP) :: handled
103  integer(I4B) :: i
104  class(*), pointer :: p
105 
106  call this%prep_event(particle, event)
107 
108  do i = 1, this%subscriptions%Count()
109  p => this%subscriptions%GetItem(i)
110  select type (subscription => p)
111  type is (particleeventsubscriptiontype)
112  handled = subscription%handler( &
113  subscription%context, &
114  particle, &
115  event)
116  if (handled) exit
117  end select
118  end do

◆ prep_event()

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

Definition at line 64 of file ParticleEvents.f90.

65  use tdismodule, only: kper, kstp
66  ! dummy
67  class(ParticleEventDispatcherType), intent(inout) :: this
68  type(ParticleType), pointer, intent(inout) :: particle
69  class(ParticleEventType), pointer, intent(in) :: event
70  ! local
71  real(DP) :: x, y, z
72 
73  x = particle%x
74  y = particle%y
75  z = particle%z
76  ! switch to model coordinates if needed
77  call particle%get_model_coords(x, y, z)
78 
79  event%kper = kper
80  event%kstp = kstp
81  event%imdl = particle%imdl
82  event%iprp = particle%iprp
83  event%irpt = particle%irpt
84  event%ilay = particle%ilay
85  event%icu = particle%icu
86  event%izone = particle%izone
87  event%trelease = particle%trelease
88  event%ttrack = particle%ttrack
89  event%x = x
90  event%y = y
91  event%z = z
92  event%istatus = particle%istatus
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,
procedure(handle_event handler,
class(*), pointer  context 
)
private

Definition at line 47 of file ParticleEvents.f90.

48  class(ParticleEventDispatcherType), intent(inout) :: this
49  procedure(handle_event) :: handler
50  class(*), pointer :: context
51  ! local
52  type(ParticleEventSubscriptionType), pointer :: subscription
53  class(*), pointer :: p
54 
55  allocate (subscription)
56  subscription%handler => handler
57  subscription%context => context
58  p => subscription
59  call this%subscriptions%Add(p)