BioCMAMC-ST
model_kernel.hpp
1#ifndef __SIMULATION_MC_KERNEL_HPP
2#define __SIMULATION_MC_KERNEL_HPP
3
4#include "mc/alias.hpp"
5#include <Kokkos_Assert.hpp>
6#include <Kokkos_Core.hpp>
7#include <Kokkos_Printf.hpp>
8#include <Kokkos_Random.hpp>
9#include <biocma_cst_config.hpp>
10#include <cassert>
11#include <mc/domain.hpp>
12#include <mc/events.hpp>
13#include <mc/particles_container.hpp>
14#include <mc/prng/prng.hpp>
15#include <mc/traits.hpp>
16#include <simulation/alias.hpp>
17#include <simulation/probability_leaving.hpp>
18#include <simulation/probe.hpp>
19#include <utility>
20
21#define CHECK_STATUS_OR_RETURN(__idx__) \
22 if (particles.status(__idx__) != MC::Status::Idle) [[unlikely]] \
23 { \
24 return; \
25 }
26
28{
30 {
31 };
32 struct TagCycle
33 {
34 };
35
37 {
39 std::size_t dead_total;
40 };
41
42 template <class Space> class CycleReducer
43 {
44 public:
45 // Required for Concept
48 using result_view_type = Kokkos::View<value_type, Space>;
49
50 KOKKOS_INLINE_FUNCTION
51 void join(value_type& dest, const value_type& src) const
52 {
53 dest.dead_total += src.dead_total;
55 }
56
57 KOKKOS_INLINE_FUNCTION
59 {
60 return *value.data();
61 }
62
63 KOKKOS_INLINE_FUNCTION
65 {
66 return value;
67 }
68
69 KOKKOS_INLINE_FUNCTION
70 bool references_scalar() const
71 {
73 }
74
75 // Optional
76 KOKKOS_INLINE_FUNCTION
77 void init(value_type& val) const
78 {
79 val.dead_total = 0;
81 }
82
83 // KOKKOS_INLINE_FUNCTION
84 // void final(value_type& val) const
85 // {
86 // // NOP
87 // }
88
89 // Part of Build-In reducers for Kokkos
90 KOKKOS_INLINE_FUNCTION
92 {
93 }
94
95 KOKKOS_INLINE_FUNCTION
97 : value(value_), references_scalar_v(false)
98 {
99 }
100
101 private:
104 };
105 template <ModelType M> struct CycleFunctor
106 {
107 using TeamPolicy = Kokkos::TeamPolicy<ComputeSpace>;
108 using TeamMember = TeamPolicy::member_type;
110
111 CycleFunctor() = default;
112
113 KOKKOS_INLINE_FUNCTION
115 MC::pool_type _random_pool,
116 MC::KernelConcentrationType&& _concentrations,
117 MC::ContributionView _contribs_scatter,
118 MC::EventContainer _event)
119 : d_t(0.), particles(std::move(_particles)), random_pool(_random_pool),
120 concentrations(std::move(_concentrations)),
121 contribs_scatter(std::move(_contribs_scatter)),
122 events(std::move(_event))
123 {
124 }
125
126 void update(double _d_t, MC::ParticlesContainer<M> _particles)
127 {
128 this->d_t = _d_t;
129 this->particles = std::move(_particles);
130 }
131
132 // KOKKOS_INLINE_FUNCTION
133 // void operator()(const TagFirstPass _tag,
134 // const TeamMember& team_handle) const
135 // {
136 // (void)_tag;
137 // GET_INDEX(particles.n_particles());
138 // if (particles.status(idx) != MC::Status::Idle) [[unlikely]]
139 // {
140 // return;
141 // }
142 // particles.get_contributions(idx, contribs_scatter);
143 // }
144 //
145 KOKKOS_INLINE_FUNCTION
146 void operator()(const TagContribution _tag, const std::size_t idx) const
147 {
148 (void)_tag;
149 CHECK_STATUS_OR_RETURN(idx);
150
151 particles.get_contributions(idx, contribs_scatter);
152 }
153
154 // KOKKOS_INLINE_FUNCTION void operator()(const TagCycle _tag,
155 // const TeamMember& team_handle,
156 // value_type& reduce_val) const
157 // {
158
159 // (void)_tag;
160 // (void)reduce_val.dead_total; // Counter not used currently because
161 // // there is no cell mortality
162 // // GET_INDEX(particles.n_particles());
163
164 // std ::size_t start_idx =
165 // (team_handle.league_rank() * team_handle.team_size());
166 // std ::size_t idx = start_idx + team_handle.team_rank();
167 // // if (particles.status(idx) != MC::Status::Idle) [[unlikely]]
168 // // {
169 // // return;
170 // // }
171 // const std::size_t n_species = concentrations.extent(0);
172 // Kokkos::View<double**, Kokkos::MemoryTraits<Kokkos::Unmanaged>>
173 // shmem_view(
174 // team_handle.team_shmem(), team_handle.team_size(), n_species);
175
176 // Kokkos::parallel_for(
177 // Kokkos::TeamThreadRange(team_handle, 0, team_handle.team_size()),
178 // [&](int i)
179 // {
180 // std::size_t particle_idx = start_idx + i;
181 // const auto position = particles.position(particle_idx);
182 // for (int s = 0; s < n_species; ++s)
183 // {
184 // shmem_view(i, s) = concentrations(s, position);
185 // }
186 // });
187 // // team_handle.team_barrier();
188 // CHECK_STATUS_OR_RETURN(idx);
189 // if (idx > particles.n_particles())
190 // {
191 // return;
192 // }
193 // particles.ages(idx, 1) += d_t;
194
195 // auto local_c =
196 // Kokkos::subview(shmem_view, team_handle.team_rank(), Kokkos::ALL);
197
198 // if (M::update(random_pool, d_t, idx, particles.model, local_c) ==
199 // MC::Status::Division)
200 // {
201 // if (!particles.handle_division(random_pool, idx))
202 // {
203 // reduce_val.waiting_allocation_particle += 1;
204 // events.wrap_incr<MC::EventType::Overflow>();
205 // Kokkos::printf("[KERNEL] Division Overflow\r\n");
206 // }
207 // events.wrap_incr<MC::EventType::NewParticle>();
208 // };
209 // }
210
211 KOKKOS_INLINE_FUNCTION void
213 const Kokkos::TeamPolicy<ComputeSpace>::member_type& team,
214 value_type& reduce_val) const
215 {
216 (void)_tag;
217 const auto team_size = team.team_size();
218
219 const auto work_stride = team_size * team.league_size();
220 const auto work_start = team.league_rank() * team_size + team.team_rank();
221
222 for (auto idx = work_start; idx < particles.n_particles();
223 idx += work_stride)
224 {
225 CHECK_STATUS_OR_RETURN(idx);
226 exec_per_particle(idx, reduce_val);
227 }
228 }
229
230 KOKKOS_INLINE_FUNCTION void operator()(const TagCycle _tag,
231 const std::size_t idx,
232 value_type& reduce_val) const
233 {
234
235 (void)_tag;
236 (void)reduce_val.dead_total;
237 exec_per_particle(idx, reduce_val);
238 }
239
240 KOKKOS_INLINE_FUNCTION void exec_per_particle(const std::size_t idx,
241 value_type& reduce_val) const
242 {
243 particles.ages(idx, 1) += d_t;
244 auto local_c =
245 Kokkos::subview(concentrations, Kokkos::ALL, particles.position(idx));
246
247 if (M::update(random_pool, d_t, idx, particles.model, local_c) ==
249 {
250 if (!particles.handle_division(random_pool, idx))
251 {
252 reduce_val.waiting_allocation_particle += 1;
253 events.wrap_incr<MC::EventType::Overflow>();
254 Kokkos::printf("[KERNEL] Division Overflow\r\n");
255 }
257 };
258
259 particles.get_contributions(idx, contribs_scatter);
260 }
261
262 M::FloatType d_t;
267 // kernelContribution contribs;
270 };
271
272} // namespace Simulation::KernelInline
273#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:53
Kokkos::View< value_type, Space > result_view_type
Definition model_kernel.hpp:48
KOKKOS_INLINE_FUNCTION CycleReducer(const result_view_type &value_)
Definition model_kernel.hpp:96
KOKKOS_INLINE_FUNCTION bool references_scalar() const
Definition model_kernel.hpp:70
KOKKOS_INLINE_FUNCTION void init(value_type &val) const
Definition model_kernel.hpp:77
KOKKOS_INLINE_FUNCTION result_view_type view() const
Definition model_kernel.hpp:64
KOKKOS_INLINE_FUNCTION void join(value_type &dest, const value_type &src) const
Definition model_kernel.hpp:51
CycleReduceType value_type
Definition model_kernel.hpp:47
KOKKOS_INLINE_FUNCTION CycleReducer(value_type &value_)
Definition model_kernel.hpp:91
KOKKOS_INLINE_FUNCTION value_type & reference() const
Definition model_kernel.hpp:58
bool references_scalar_v
Definition model_kernel.hpp:103
CycleReducer reducer
Definition model_kernel.hpp:46
result_view_type value
Definition model_kernel.hpp:102
@ Overflow
Definition events.hpp:22
@ NewParticle
Spawn new particle.
Definition events.hpp:18
@ Division
Definition alias.hpp:59
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:38
decltype(Kokkos::Experimental::create_scatter_view(kernelContribution())) ContributionView
Definition alias.hpp:87
Kokkos::View< const double **, Kokkos::LayoutLeft, ComputeSpace, Kokkos::MemoryTraits< Kokkos::RandomAccess > > KernelConcentrationType
Definition alias.hpp:90
Definition kernels.hpp:11
Use to count events that occurs during Monte-Carlo processing cycles.
Definition events.hpp:44
KOKKOS_INLINE_FUNCTION void operator()(const TagContribution _tag, const std::size_t idx) const
Definition model_kernel.hpp:146
MC::EventContainer events
Definition model_kernel.hpp:269
M::FloatType d_t
Definition model_kernel.hpp:262
MC::ContributionView contribs_scatter
Definition model_kernel.hpp:266
MC::KernelConcentrationType concentrations
Definition model_kernel.hpp:265
MC::ParticlesContainer< M > particles
Definition model_kernel.hpp:263
TeamPolicy::member_type TeamMember
Definition model_kernel.hpp:108
CycleReduceType value_type
Definition model_kernel.hpp:109
KOKKOS_INLINE_FUNCTION void operator()(const TagCycle _tag, const std::size_t idx, value_type &reduce_val) const
Definition model_kernel.hpp:230
KOKKOS_INLINE_FUNCTION CycleFunctor(MC::ParticlesContainer< M > _particles, MC::pool_type _random_pool, MC::KernelConcentrationType &&_concentrations, MC::ContributionView _contribs_scatter, MC::EventContainer _event)
Definition model_kernel.hpp:114
KOKKOS_INLINE_FUNCTION void operator()(TagCycle _tag, const Kokkos::TeamPolicy< ComputeSpace >::member_type &team, value_type &reduce_val) const
Definition model_kernel.hpp:212
Kokkos::TeamPolicy< ComputeSpace > TeamPolicy
Definition model_kernel.hpp:107
KOKKOS_INLINE_FUNCTION void exec_per_particle(const std::size_t idx, value_type &reduce_val) const
Definition model_kernel.hpp:240
void update(double _d_t, MC::ParticlesContainer< M > _particles)
Definition model_kernel.hpp:126
MC::KernelConcentrationType limitation_factor
Definition model_kernel.hpp:268
MC::pool_type random_pool
Definition model_kernel.hpp:264
Definition model_kernel.hpp:37
std::size_t waiting_allocation_particle
Definition model_kernel.hpp:38
std::size_t dead_total
Definition model_kernel.hpp:39
Definition model_kernel.hpp:30
Definition model_kernel.hpp:33