BioCMAMC-ST
contribution_kernel.hpp
1#ifndef __CONTRIBUTION_KERNEL_HPP__
2#define __CONTRIBUTION_KERNEL_HPP__
3#include "Kokkos_Assert.hpp"
4#include <Kokkos_Core.hpp>
5#include <mc/particles_container.hpp>
6#include <mc/traits.hpp>
7
8template <ModelType M> struct ContributionFunctor
9{
10 struct Tag0D
11 {
12 };
13 struct Tag3D
14 {
15 };
16
17 using TeamPolicy = Kokkos::TeamPolicy<ComputeSpace>;
18 using TeamMember = TeamPolicy::member_type;
19
21 = Kokkos::View<float_t*,
22 TeamPolicy::execution_space::scratch_memory_space>;
23 using float_t = float;
25 ContributionFunctor(std::size_t particle_per_team,
26 MC::ContributionView contribution_scatter,
28 : m_particle_per_team(particle_per_team),
29 m_contribution_scatter(std::move(contribution_scatter)),
30 m_particles(std::move(particles))
31 {
32 }
34 size_t np{};
35
36 void
38
39 {
40 this->m_particles = std::move(_particles);
41 np = m_particles.n_particles();
42 }
44 std::size_t m_particle_per_team;
48
49 KOKKOS_INLINE_FUNCTION
50 void
51 operator()(Tag0D _tag, const TeamMember& team) const
52 {
53 (void)_tag;
54 const std::size_t p0 = team.league_rank() * m_particle_per_team;
55 const std::size_t n_particle = m_particles.n_particles();
56 const auto& status = m_particles.status;
57 const auto& contribs = m_particles.contribs;
58 static const auto n_c = M::n_c;
59
60 ScratchView scratch(team.team_scratch(0), n_c);
61
62 Kokkos::parallel_for(Kokkos::TeamVectorRange(team, n_c),
63 [&](const std::size_t j)
64 { scratch(j) = float_t{ 0 }; });
65
66 const auto upper_bound = ((p0 + m_particle_per_team) >= n_particle)
67 ? n_particle - p0
69
70 KOKKOS_ASSERT(upper_bound > 0 && upper_bound <= m_particle_per_team
71 && (p0 + upper_bound) <= n_particle);
72
73 team.team_barrier();
74
75 Kokkos::parallel_for(
76 Kokkos::TeamThreadRange(team, 0, upper_bound),
77 [&](const std::size_t relative_index)
78 {
79 const std::size_t flatten_index = p0 + relative_index;
80 if (status(flatten_index) != MC::Status::Idle)
81 {
82 return;
83 }
84 const auto weight = m_particles.get_weight(flatten_index);
85 for (std::size_t j = 0; j < n_c; ++j)
86 {
87 Kokkos::atomic_add(&scratch(j),
88 weight * contribs(flatten_index, j));
89 }
90 });
91
92 team.team_barrier();
93 const auto& cs = m_contribution_scatter;
94 Kokkos::single(Kokkos::PerTeam(team),
95 [=]()
96 {
97 auto access = cs.access();
98 for (std::size_t j = 0; j < M::n_c; ++j)
99 {
100 access(j, 0) += scratch(j);
101 }
102 });
103 }
104
105 // KOKKOS_INLINE_FUNCTION
106 // void
107 // operator()(const Tag3D _tag, const TeamMember& team) const
108 // {
109 // (void)_tag;
110
111 // const std::size_t p0 = team.league_rank() * m_particle_per_team;
112 // const auto& c = m_particles.contribs;
113 // const auto& status = m_particles.status;
114 // const auto n_particle = m_particles.n_particles();
115 // const auto& positions = m_particles.position;
116
117 // const auto upper_bound = ((p0 + m_particle_per_team) >= n_particle)
118 // ? n_particle - p0 - 1
119 // : m_particle_per_team;
120 // KOKKOS_ASSERT(upper_bound >= 0 && upper_bound < n_particle);
121
122 // constexpr std::size_t work_per_thread = 32;
123 // KOKKOS_ASSERT(m_particle_per_team%work_per_thread==0);
124
125 // Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, upper_bound),
126 // [&](const std::size_t i)
127 // {
128 // const std::size_t p = p0 + i;
129 // if (status(p) != MC::Status::Idle)
130 // {
131 // return;
132 // }
133 // auto access = m_contribution_scatter.access();
134 // const double weight = m_particles.get_weight(p);
135 // const auto pos = positions(p);
136
137 // Kokkos::parallel_for(
138 // Kokkos::ThreadVectorRange(team, 0, M::n_c),
139 // [&](const int j)
140 // { access(j, pos) += weight * c(p, j); });
141 // });
142
143 // }
144
145 KOKKOS_INLINE_FUNCTION
146 void
147 operator()(const Tag3D _tag, const TeamMember& team) const
148 {
149 (void)_tag;
150
151 const std::size_t p0 = team.league_rank() * m_particle_per_team;
152 const auto& c = m_particles.contribs;
153 const auto& status = m_particles.status;
154 const auto n_particle = m_particles.n_particles();
155 const auto& positions = m_particles.position;
156 constexpr std::size_t work_per_thread = 32;
157
158 const auto upper_bound = ((p0 + m_particle_per_team) >= n_particle)
159 ? n_particle - p0
161 KOKKOS_ASSERT(upper_bound >= 0 && upper_bound < n_particle);
162 const auto range = (upper_bound + work_per_thread - 1) / work_per_thread;
163
164 KOKKOS_ASSERT(m_particle_per_team % work_per_thread == 0);
165
166 Kokkos::parallel_for(
167 Kokkos::TeamThreadRange(team, 0, range),
168 [&](const std::size_t i)
169 {
170 auto access = m_contribution_scatter.access();
171 for (std::size_t k = 0; k < work_per_thread; ++k)
172 {
173 const std::size_t p = p0 + i * work_per_thread + k;
174 if (status(p) != MC::Status::Idle || i >= n_particle)
175 {
176 return;
177 }
178 const double weight = m_particles.get_weight(p);
179 const auto pos = positions(p);
180 Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, 0, M::n_c),
181 [&](const int j)
182 { access(j, pos) += weight * c(p, j); });
183 }
184 });
185 }
186};
187
188#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:57
decltype(Kokkos::Experimental::create_scatter_view( kernelContribution())) ContributionView
Definition alias.hpp:162
@ Idle
Definition alias.hpp:126
Definition contribution_kernel.hpp:11
Definition contribution_kernel.hpp:14
TeamPolicy::member_type TeamMember
Definition contribution_kernel.hpp:18
KOKKOS_INLINE_FUNCTION void operator()(Tag0D _tag, const TeamMember &team) const
Definition contribution_kernel.hpp:50
ContributionFunctor(std::size_t particle_per_team, MC::ContributionView contribution_scatter, MC::ParticlesContainer< M > particles)
Definition contribution_kernel.hpp:24
Kokkos::View< float_t *, TeamPolicy::execution_space::scratch_memory_space > ScratchView
Definition contribution_kernel.hpp:20
Kokkos::TeamPolicy< ComputeSpace > TeamPolicy
Definition contribution_kernel.hpp:17
void update(MC::ParticlesContainer< M > _particles)
Definition contribution_kernel.hpp:36
float float_t
Definition contribution_kernel.hpp:22
size_t np
Definition contribution_kernel.hpp:33
MC::ParticlesContainer< M > m_particles
Definition contribution_kernel.hpp:46
MC::ContributionView m_contribution_scatter
Definition contribution_kernel.hpp:45
std::size_t m_particle_per_team
Definition contribution_kernel.hpp:43