MODFLOW 6  version 6.7.0.dev2
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, totimc
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, per, stp
51  real(DP) :: x, y, z
52  class(*), pointer :: p
53 
54  ! If tracking time falls exactly on a boundary between time steps,
55  ! report the previous time step for this datum. This is to follow
56  ! MP7's behavior, and because the particle will have been tracked
57  ! up to this instant under the previous time step's conditions, so
58  ! the time step we're about to start shouldn't get "credit" for it.
59  per = kper
60  stp = kstp
61  if (particle%ttrack == totimc .and. (per > 1 .or. stp > 1)) then
62  if (stp > 1) then
63  stp = stp - 1
64  else if (per > 1) then
65  per = per - 1
66  stp = 1
67  end if
68  end if
69 
70  ! Convert to model coordinates if we need to
71  x = particle%x
72  y = particle%y
73  z = particle%z
74  call particle%get_model_coords(x, y, z)
75 
76  event%kper = per
77  event%kstp = stp
78  event%imdl = particle%imdl
79  event%iprp = particle%iprp
80  event%irpt = particle%irpt
81  event%ilay = particle%ilay
82  event%icu = particle%icu
83  event%izone = particle%izone
84  event%trelease = particle%trelease
85  event%ttrack = particle%ttrack
86  event%x = x
87  event%y = y
88  event%z = z
89  event%istatus = particle%istatus
90 
91  do i = 1, this%consumers%Count()
92  p => this%consumers%GetItem(i)
93  select type (consumer => p)
94  class is (particleeventconsumertype)
95  call consumer%handle_event(particle, event)
96  end select
97  end do
98  end subroutine dispatch
99 
100  !> @brief Destroy the dispatcher.
101  subroutine destroy(this)
102  class(particleeventdispatchertype), intent(inout) :: this
103  call this%consumers%Clear()
104  end subroutine destroy
105 
106 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.
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
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
Base type for particle events.
Particle tracked by the PRT model.
Definition: Particle.f90:56