MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
ParticleEvents.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
3  use listmodule, only: listtype
6  implicit none
7 
8  private
9 
10  type, public, abstract :: particleeventconsumertype
11  contains
12  procedure(handle_event), deferred :: handle_event
14 
16  type(listtype) :: consumers
17  contains
18  procedure, public :: subscribe
19  procedure, public :: dispatch
20  procedure :: destroy
22 
23  abstract interface
24  subroutine handle_event(this, particle, event)
26  class(particleeventconsumertype), intent(inout) :: this
27  type(particletype), pointer, intent(in) :: particle
28  class(particleeventtype), pointer, intent(in) :: event
29  end subroutine handle_event
30  end interface
31 
32 contains
33  !> @brief Subscribe a consumer to the dispatcher.
34  subroutine subscribe(this, consumer)
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)
40  end subroutine subscribe
41 
42  !> @brief Dispatch an event.
43  subroutine dispatch(this, particle, event)
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
82  end subroutine dispatch
83 
84  !> @brief Destroy the dispatcher.
85  subroutine destroy(this)
86  class(particleeventdispatchertype), intent(inout) :: this
87  call this%consumers%Clear()
88  end subroutine destroy
89 
90 end module particleeventsmodule
This module defines variable data types.
Definition: kind.f90:8
subroutine subscribe(this, consumer)
Subscribe a consumer to the dispatcher.
subroutine dispatch(this, particle, event)
Dispatch an event.
subroutine destroy(this)
Destroy the dispatcher.
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
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
Base type for particle events.
Particle tracked by the PRT model.
Definition: Particle.f90:57