BioCMAMC-ST
model_kernel.hpp
1#ifndef __SIMULATION_MC_KERNEL_HPP
2#define __SIMULATION_MC_KERNEL_HPP
3
4#include <Kokkos_Assert.hpp>
5#include <Kokkos_Core.hpp>
6#include <Kokkos_Macros.hpp>
7#include <Kokkos_Printf.hpp>
8#include <Kokkos_Random.hpp>
9#include <biocma_cst_config.hpp>
10#include <cassert>
11#include <common/common.hpp>
12#include <mc/alias.hpp>
13#include <mc/domain.hpp>
14#include <mc/events.hpp>
15#include <mc/particles_container.hpp>
16#include <mc/prng/prng.hpp>
17#include <mc/traits.hpp>
18#include <simulation/probability_leaving.hpp>
19#include <simulation/probe.hpp>
20#include <utility>
21
22#define CHECK_STATUS_OR_RETURN(__idx__) \
23 if (particles.status(__idx__) != MC::Status::Idle) [[unlikely]] \
24 { \
25 return; \
26 }
27
29{
31 {
32 };
34 {
35 };
36 struct TagCycle
37 {
38 };
39
41 {
43 std::size_t dead_total;
44
45 KOKKOS_INLINE_FUNCTION CycleReduceType&
47 {
48 this->waiting_allocation_particle += a.waiting_allocation_particle;
49 this->dead_total += a.dead_total;
50 return *this;
51 }
52 };
53
54 template <class Space> class CycleReducer
55 {
56 public:
57 // Required for Concept
60 using result_view_type = Kokkos::View<value_type, Space>;
61
62 KOKKOS_INLINE_FUNCTION
63 void
64 join(value_type& dest, const value_type& src) const
65 {
66 // dest.dead_total += src.dead_total;
67 // dest.waiting_allocation_particle += src.waiting_allocation_particle;
68 dest += src;
69 }
70
71 [[nodiscard]] KOKKOS_INLINE_FUNCTION value_type&
72 reference() const
73 {
74 return *value.data();
75 }
76
77 KOKKOS_INLINE_FUNCTION
79 view() const
80 {
81 return value;
82 }
83
84 [[nodiscard]] KOKKOS_INLINE_FUNCTION bool
86 {
88 }
89
90 // Optional
91 KOKKOS_INLINE_FUNCTION
92 void
93 init(value_type& val) const
94 {
95 val.dead_total = 0;
97 }
98
99 // KOKKOS_INLINE_FUNCTION
100 // void final(value_type& val) const
101 // {
102 // // NOP
103 // }
104
105 // Part of Build-In reducers for Kokkos
106 KOKKOS_INLINE_FUNCTION
107 explicit CycleReducer(value_type& value_)
108 : value(&value_), references_scalar_v(true)
109 {
110 }
111
112 KOKKOS_INLINE_FUNCTION
113 explicit CycleReducer(const result_view_type& value_)
114 : value(value_), references_scalar_v(false)
115 {
116 }
117
118 private:
121 };
122
123 template <ModelType M> struct CycleFunctor
124 {
125 using TeamPolicy = Kokkos::TeamPolicy<ComputeSpace>;
126 using TeamMember = TeamPolicy::member_type;
128
129 CycleFunctor() = default;
130
131 std::size_t n_p{};
132
133 std::size_t m_p_team;
134
135 KOKKOS_INLINE_FUNCTION
136 CycleFunctor(std::size_t p_per_team,
137 MC::ParticlesContainer<M> _particles,
138 MC::pool_type _random_pool,
139 MC::KernelConcentrationType&& _concentrations,
140 MC::EventContainer _event,
142 : m_p_team(p_per_team), d_t(0.), particles(std::move(_particles)),
143 random_pool(std::move(_random_pool)),
144 concentrations(std::move(_concentrations)), events(std::move(_event)),
145 probes(std::move(_probes))
146 {
147 }
148
149 [[nodiscard]] constexpr bool
151 {
152 return M::n_c > 0;
153 }
154
155 void
156 update(double _d_t, MC::ParticlesContainer<M> _particles)
157 {
158 this->d_t = _d_t;
159 this->particles = std::move(_particles);
160 n_p = this->particles.n_particles();
161 }
162
163 KOKKOS_INLINE_FUNCTION void
165 const TeamMember& team,
166 value_type& reduce_val) const
167 {
168
169 (void)_tag;
170 const std::size_t count = m_p_team;
171 const std::size_t p0 = team.league_rank() * count;
172 const std::size_t n_particle = n_p;
173 const auto _d_t = d_t;
174 const auto& status = particles.status;
175 const auto& ages = particles.ages;
176
177 const auto upper_bound
178 = ((p0 + count) >= n_particle) ? n_particle - p0 : count;
179
180 KOKKOS_ASSERT(upper_bound > 0);
181 KOKKOS_ASSERT(upper_bound < n_particle);
182
183 value_type local;
184 Kokkos::parallel_reduce(
185 Kokkos::TeamThreadRange(team, 0, upper_bound),
186 [&](std::size_t relative_index, value_type& lv)
187 {
188 const std::size_t flatten_index = p0 + relative_index;
189 const bool active = status(flatten_index) == MC::Status::Idle;
190
191 if (active)
192 {
193 ages(flatten_index, 1) += _d_t;
194 exec_per_particle(flatten_index, lv);
195 }
196 },
197 local);
198 // Kokkos::parallel_reduce(
199 // Kokkos::TeamThreadRange(team, count),
200 // [&](std::size_t relative_index, value_type& lv)
201 // {
202 // const std::size_t flatten_index = p0 + relative_index;
203 // const bool active = flatten_index < n_particle
204 // && status(flatten_index) == MC::Status::Idle;
205
206 // if (active)
207 // {
208 // ages(flatten_index, 1) += _d_t;
209
210 // exec_per_particle(flatten_index, lv);
211 // }
212 // },
213 // local);
214 team.team_barrier();
215
216 reduce_val += local;
217 }
218
219 KOKKOS_FORCEINLINE_FUNCTION void
221 const std::size_t idx,
222 value_type& reduce_val) const
223 {
224
225 (void)_tag;
226 (void)reduce_val.dead_total;
227 exec_per_particle(idx, reduce_val);
228 }
229
230 KOKKOS_INLINE_FUNCTION void
231 exec_per_particle(const std::size_t idx, value_type& reduce_val) const
232 {
233 using mem_space = ComputeSpace::memory_space;
234
235 const auto new_status = M::update(random_pool,
236 d_t,
237 idx,
238 particles.model,
239 particles.contribs,
240 particles.position(idx),
242
243 if (new_status == MC::Status::Division)
244 {
245 if constexpr (AutoGenerated::FlagCompileTime::use_probe)
246 {
247 // Register probe here to sample BEFORE division and even if division
248 // procedure fails to spawn new particle, the age is still the same
249 const auto _ = this->probes.template set<mem_space>(
250 particles.ages(idx, 1)); // Skip error
251 }
252
253 if (!particles.handle_division(random_pool, idx)) [[unlikely]]
254 {
255 reduce_val.waiting_allocation_particle += 1;
256 events.wrap_incr<MC::EventType::Overflow>();
257 Kokkos::printf("[KERNEL] Division Overflow\r\n");
258 }
260
261 /*
262 TODO, it seems that after some calculation (example cstr 0d)
263 size(division probes) = number of tally+1 where the last probe is 0
264 where size(leaving time) == tally
265 Is it a bug ?
266 */
267 };
268 }
269
270 M::FloatType d_t;
276 };
277
278} // namespace Simulation::KernelInline
279#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:57
Kokkos::View< value_type, Space > result_view_type
Definition model_kernel.hpp:60
KOKKOS_INLINE_FUNCTION CycleReducer(const result_view_type &value_)
Definition model_kernel.hpp:113
KOKKOS_INLINE_FUNCTION bool references_scalar() const
Definition model_kernel.hpp:85
KOKKOS_INLINE_FUNCTION void init(value_type &val) const
Definition model_kernel.hpp:93
KOKKOS_INLINE_FUNCTION result_view_type view() const
Definition model_kernel.hpp:79
KOKKOS_INLINE_FUNCTION void join(value_type &dest, const value_type &src) const
Definition model_kernel.hpp:64
CycleReduceType value_type
Definition model_kernel.hpp:59
KOKKOS_INLINE_FUNCTION CycleReducer(value_type &value_)
Definition model_kernel.hpp:107
KOKKOS_INLINE_FUNCTION value_type & reference() const
Definition model_kernel.hpp:72
bool references_scalar_v
Definition model_kernel.hpp:120
CycleReducer reducer
Definition model_kernel.hpp:58
result_view_type value
Definition model_kernel.hpp:119
@ Overflow
Definition events.hpp:23
@ NewParticle
Spawn new particle.
Definition events.hpp:19
@ Division
Definition alias.hpp:127
@ Idle
Definition alias.hpp:126
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:100
Kokkos::View< const double **, Kokkos::LayoutLeft, ComputeSpace, Kokkos::MemoryTraits< Kokkos::RandomAccess > > KernelConcentrationType
Definition alias.hpp:165
Definition kernels.hpp:16
Probes< AutoGenerated::probe_buffer_size > ProbeAutogeneratedBuffer
Definition probe.hpp:148
Use to count events that occurs during Monte-Carlo processing cycles.
Definition events.hpp:150
KOKKOS_FORCEINLINE_FUNCTION void operator()(const TagCycle _tag, const std::size_t idx, value_type &reduce_val) const
Definition model_kernel.hpp:220
MC::EventContainer events
Definition model_kernel.hpp:274
M::FloatType d_t
Definition model_kernel.hpp:270
MC::KernelConcentrationType concentrations
Definition model_kernel.hpp:273
MC::ParticlesContainer< M > particles
Definition model_kernel.hpp:271
TeamPolicy::member_type TeamMember
Definition model_kernel.hpp:126
CycleReduceType value_type
Definition model_kernel.hpp:127
KOKKOS_INLINE_FUNCTION CycleFunctor(std::size_t p_per_team, MC::ParticlesContainer< M > _particles, MC::pool_type _random_pool, MC::KernelConcentrationType &&_concentrations, MC::EventContainer _event, ProbeAutogeneratedBuffer _probes)
Definition model_kernel.hpp:136
std::size_t n_p
Definition model_kernel.hpp:131
KOKKOS_INLINE_FUNCTION void operator()(TagCycle _tag, const TeamMember &team, value_type &reduce_val) const
Definition model_kernel.hpp:164
ProbeAutogeneratedBuffer probes
Definition model_kernel.hpp:275
Kokkos::TeamPolicy< ComputeSpace > TeamPolicy
Definition model_kernel.hpp:125
std::size_t m_p_team
Definition model_kernel.hpp:133
KOKKOS_INLINE_FUNCTION void exec_per_particle(const std::size_t idx, value_type &reduce_val) const
Definition model_kernel.hpp:231
void update(double _d_t, MC::ParticlesContainer< M > _particles)
Definition model_kernel.hpp:156
MC::pool_type random_pool
Definition model_kernel.hpp:272
constexpr bool do_contribs() const
Definition model_kernel.hpp:150
Definition model_kernel.hpp:41
std::size_t waiting_allocation_particle
Definition model_kernel.hpp:42
std::size_t dead_total
Definition model_kernel.hpp:43
KOKKOS_INLINE_FUNCTION CycleReduceType & operator+=(const CycleReduceType &a)
Definition model_kernel.hpp:46
Definition model_kernel.hpp:34
Definition model_kernel.hpp:31
Definition model_kernel.hpp:37