BioCMAMC-ST
model_kernel.hpp
1#ifndef __SIMULATION_MC_KERNEL_HPP
2#define __SIMULATION_MC_KERNEL_HPP
3
4#include "Kokkos_Macros.hpp"
5#include "mc/alias.hpp"
6#include <Kokkos_Assert.hpp>
7#include <Kokkos_Core.hpp>
8#include <Kokkos_Printf.hpp>
9#include <Kokkos_Random.hpp>
10#include <biocma_cst_config.hpp>
11#include <cassert>
12#include <mc/domain.hpp>
13#include <mc/events.hpp>
14#include <mc/particles_container.hpp>
15#include <mc/prng/prng.hpp>
16#include <mc/traits.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
52 join(value_type& dest, const value_type& src) const
53 {
54 dest.dead_total += src.dead_total;
56 }
57
58 KOKKOS_INLINE_FUNCTION
60 reference() const
61 {
62 return *value.data();
63 }
64
65 KOKKOS_INLINE_FUNCTION
67 view() const
68 {
69 return value;
70 }
71
72 KOKKOS_INLINE_FUNCTION
73 bool
75 {
77 }
78
79 // Optional
80 KOKKOS_INLINE_FUNCTION
81 void
82 init(value_type& val) const
83 {
84 val.dead_total = 0;
86 }
87
88 // KOKKOS_INLINE_FUNCTION
89 // void final(value_type& val) const
90 // {
91 // // NOP
92 // }
93
94 // Part of Build-In reducers for Kokkos
95 KOKKOS_INLINE_FUNCTION
97 {
98 }
99
100 KOKKOS_INLINE_FUNCTION
102 : value(value_), references_scalar_v(false)
103 {
104 }
105
106 private:
109 };
110
111 template <ModelType M> struct CycleFunctor
112 {
113 using TeamPolicy = Kokkos::TeamPolicy<ComputeSpace>;
114 using TeamMember = TeamPolicy::member_type;
116
117 CycleFunctor() = default;
118
119 KOKKOS_INLINE_FUNCTION
121 MC::pool_type _random_pool,
122 MC::KernelConcentrationType&& _concentrations,
123 MC::ContributionView _contribs_scatter,
124 MC::EventContainer _event,
126 : d_t(0.), particles(std::move(_particles)),
127 random_pool(std::move(_random_pool)),
128 concentrations(std::move(_concentrations)),
129 contribs_scatter(std::move(_contribs_scatter)),
130 events(std::move(_event)), probes(std::move(_probes))
131 {
132 }
133
134 void
135 update(double _d_t, MC::ParticlesContainer<M> _particles)
136 {
137 this->d_t = _d_t;
138 this->particles = std::move(_particles);
139 }
140
141 KOKKOS_INLINE_FUNCTION
142 void
143 operator()(const TagContribution _tag, const TeamMember& team) const
144 {
145 (void)_tag;
146 // std ::size_t idx = (team_handle.league_rank() *
147 // team_handle.team_size()) +
148 // team_handle.team_rank();
149 // if (idx >= (particles.n_particles()))
150 // {
151 // return;
152 // }
153 // if (particles.status(idx) != MC::Status::Idle) [[unlikely]]
154 // {
155 // return;
156 // }
157 constexpr int PARTICLES_PER_TEAM = 256;
158 // particles.get_contributions(idx, contribs_scatter);
159
160 const std::size_t p0 = team.league_rank() * PARTICLES_PER_TEAM;
161
162 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, PARTICLES_PER_TEAM),
163 [&](const int i)
164 {
165 const std::size_t p = p0 + i;
166 if (p >= particles.n_particles())
167 {
168 return;
169 }
170
171 if (particles.status(p) != MC::Status::Idle)
172 {
173 return;
174 }
175 particles.get_contributions(p, contribs_scatter);
176 });
177 }
178
179 // KOKKOS_INLINE_FUNCTION
180 // void operator()(const TagContribution _tag, const std::size_t idx) const
181 // {
182 // (void)_tag;
183 // CHECK_STATUS_OR_RETURN(idx);
184
185 // particles.get_contributions(idx, contribs_scatter);
186 // }
187
188 // KOKKOS_INLINE_FUNCTION void operator()(const TagCycle _tag,
189 // const TeamMember& team_handle,
190 // value_type& reduce_val) const
191 // {
192
193 // (void)_tag;
194 // (void)reduce_val.dead_total; // Counter not used currently because
195 // // there is no cell mortality
196 // // GET_INDEX(particles.n_particles());
197
198 // std ::size_t start_idx =
199 // (team_handle.league_rank() * team_handle.team_size());
200 // std ::size_t idx = start_idx + team_handle.team_rank();
201 // // if (particles.status(idx) != MC::Status::Idle) [[unlikely]]
202 // // {
203 // // return;
204 // // }
205 // const std::size_t n_species = concentrations.extent(0);
206 // Kokkos::View<double**, Kokkos::MemoryTraits<Kokkos::Unmanaged>>
207 // shmem_view(
208 // team_handle.team_shmem(), team_handle.team_size(), n_species);
209
210 // Kokkos::parallel_for(
211 // Kokkos::TeamThreadRange(team_handle, 0, team_handle.team_size()),
212 // [&](int i)
213 // {
214 // std::size_t particle_idx = start_idx + i;
215 // const auto position = particles.position(particle_idx);
216 // for (int s = 0; s < n_species; ++s)
217 // {
218 // shmem_view(i, s) = concentrations(s, position);
219 // }
220 // });
221 // // team_handle.team_barrier();
222 // CHECK_STATUS_OR_RETURN(idx);
223 // if (idx > particles.n_particles())
224 // {
225 // return;
226 // }
227 // particles.ages(idx, 1) += d_t;
228
229 // auto local_c =
230 // Kokkos::subview(shmem_view, team_handle.team_rank(), Kokkos::ALL);
231
232 // if (M::update(random_pool, d_t, idx, particles.model, local_c) ==
233 // MC::Status::Division)
234 // {
235 // if (!particles.handle_division(random_pool, idx))
236 // {
237 // reduce_val.waiting_allocation_particle += 1;
238 // events.wrap_incr<MC::EventType::Overflow>();
239 // Kokkos::printf("[KERNEL] Division Overflow\r\n");
240 // }
241 // events.wrap_incr<MC::EventType::NewParticle>();
242 // };
243 // }
244
245 KOKKOS_INLINE_FUNCTION void
247 const Kokkos::TeamPolicy<ComputeSpace>::member_type& team,
248 value_type& reduce_val) const
249 {
250 (void)_tag;
251 const auto team_size = team.team_size();
252
253 const auto work_stride = team_size * team.league_size();
254 const auto work_start = team.league_rank() * team_size + team.team_rank();
255
256 for (auto idx = work_start; idx < particles.n_particles();
257 idx += work_stride)
258 {
259 CHECK_STATUS_OR_RETURN(idx);
260 exec_per_particle(idx, reduce_val);
261 }
262 }
263
264 KOKKOS_FORCEINLINE_FUNCTION void
266 const std::size_t idx,
267 value_type& reduce_val) const
268 {
269
270 (void)_tag;
271 (void)reduce_val.dead_total;
272 exec_per_particle(idx, reduce_val);
273 }
274
275 KOKKOS_INLINE_FUNCTION void
276 exec_per_particle(const std::size_t idx, value_type& reduce_val) const
277 {
278 particles.ages(idx, 1) += d_t;
279 const auto local_c = Kokkos::subview(
280 concentrations, Kokkos::ALL, particles.position(idx));
281 const auto new_status
282 = M::update(random_pool, d_t, idx, particles.model, local_c);
283 if (new_status == MC::Status::Division)
284 {
285 if constexpr (AutoGenerated::FlagCompileTime::use_probe)
286 {
287 // Register probe here to sample BEFORE division and even if division
288 // procedure fails to spawn new particle, the age is still the same
289 const auto _ = this->probes.set(particles.ages(idx, 1)); // Skip error
290 }
291
292 if (!particles.handle_division(random_pool, idx))
293 {
294 reduce_val.waiting_allocation_particle += 1;
295 events.wrap_incr<MC::EventType::Overflow>();
296 Kokkos::printf("[KERNEL] Division Overflow\r\n");
297 }
299 };
300 }
301
302 M::FloatType d_t;
307 // kernelContribution contribs;
311 };
312
313} // namespace Simulation::KernelInline
314#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:56
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:101
KOKKOS_INLINE_FUNCTION bool references_scalar() const
Definition model_kernel.hpp:74
KOKKOS_INLINE_FUNCTION void init(value_type &val) const
Definition model_kernel.hpp:82
KOKKOS_INLINE_FUNCTION result_view_type view() const
Definition model_kernel.hpp:67
KOKKOS_INLINE_FUNCTION void join(value_type &dest, const value_type &src) const
Definition model_kernel.hpp:52
CycleReduceType value_type
Definition model_kernel.hpp:47
KOKKOS_INLINE_FUNCTION CycleReducer(value_type &value_)
Definition model_kernel.hpp:96
KOKKOS_INLINE_FUNCTION value_type & reference() const
Definition model_kernel.hpp:60
bool references_scalar_v
Definition model_kernel.hpp:108
CycleReducer reducer
Definition model_kernel.hpp:46
result_view_type value
Definition model_kernel.hpp:107
@ Overflow
Definition events.hpp:22
@ NewParticle
Spawn new particle.
Definition events.hpp:18
decltype(Kokkos::Experimental::create_scatter_view( kernelContribution())) ContributionView
Definition alias.hpp:88
@ Division
Definition alias.hpp:60
@ Idle
Definition alias.hpp:59
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:39
Kokkos::View< const double **, Kokkos::LayoutLeft, ComputeSpace, Kokkos::MemoryTraits< Kokkos::RandomAccess > > KernelConcentrationType
Definition alias.hpp:91
Definition kernels.hpp:12
Probes< AutoGenerated::probe_buffer_size > ProbeAutogeneratedBuffer
Definition probe.hpp:55
Use to count events that occurs during Monte-Carlo processing cycles.
Definition events.hpp:45
KOKKOS_INLINE_FUNCTION CycleFunctor(MC::ParticlesContainer< M > _particles, MC::pool_type _random_pool, MC::KernelConcentrationType &&_concentrations, MC::ContributionView _contribs_scatter, MC::EventContainer _event, ProbeAutogeneratedBuffer _probes)
Definition model_kernel.hpp:120
KOKKOS_FORCEINLINE_FUNCTION void operator()(const TagCycle _tag, const std::size_t idx, value_type &reduce_val) const
Definition model_kernel.hpp:265
MC::EventContainer events
Definition model_kernel.hpp:309
M::FloatType d_t
Definition model_kernel.hpp:302
MC::ContributionView contribs_scatter
Definition model_kernel.hpp:306
MC::KernelConcentrationType concentrations
Definition model_kernel.hpp:305
MC::ParticlesContainer< M > particles
Definition model_kernel.hpp:303
TeamPolicy::member_type TeamMember
Definition model_kernel.hpp:114
CycleReduceType value_type
Definition model_kernel.hpp:115
KOKKOS_INLINE_FUNCTION void operator()(TagCycle _tag, const Kokkos::TeamPolicy< ComputeSpace >::member_type &team, value_type &reduce_val) const
Definition model_kernel.hpp:246
ProbeAutogeneratedBuffer probes
Definition model_kernel.hpp:310
Kokkos::TeamPolicy< ComputeSpace > TeamPolicy
Definition model_kernel.hpp:113
KOKKOS_INLINE_FUNCTION void exec_per_particle(const std::size_t idx, value_type &reduce_val) const
Definition model_kernel.hpp:276
void update(double _d_t, MC::ParticlesContainer< M > _particles)
Definition model_kernel.hpp:135
MC::KernelConcentrationType limitation_factor
Definition model_kernel.hpp:308
MC::pool_type random_pool
Definition model_kernel.hpp:304
KOKKOS_INLINE_FUNCTION void operator()(const TagContribution _tag, const TeamMember &team) const
Definition model_kernel.hpp:143
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