BioCMAMC-ST
post_process.hpp
1#ifndef __CORE_POST_PROCESS_PUBLIC_HPP__
2#define __CORE_POST_PROCESS_PUBLIC_HPP__
3#include "Kokkos_Macros.hpp"
4#include <Kokkos_Core.hpp>
5#include <Kokkos_Core_fwd.hpp>
6#include <Kokkos_ScatterView.hpp>
7#include <common/common.hpp>
8#include <mc/alias.hpp>
9#include <mc/particles_container.hpp>
10#include <mc/traits.hpp>
11#include <stdexcept>
12
13template <typename MemorySpace>
14using ParticlePropertyViewType
15 = Kokkos::View<double**, Kokkos::LayoutRight, MemorySpace>;
16
17using SubViewtype = Kokkos::Subview<ParticlePropertyViewType<ComputeSpace>,
18 decltype(Kokkos::ALL),
19 std::size_t>;
20namespace PostProcessing
22 struct BonceBuffer
23 {
24 std::optional<ParticlePropertyViewType<HostSpace>> particle_values;
25 std::optional<ParticlePropertyViewType<HostSpace>> spatial_values;
26 std::optional<ParticlePropertyViewType<HostSpace>> ages;
27 std::optional<std::vector<std::string>> vnames;
28 };
29
30 namespace
31 {
32 template <typename Model, typename ExecutionSpace>
33 struct GetPropertiesFunctor
34 {
35 using ComputeSpace = typename ExecutionSpace::memory_space;
36
37 GetPropertiesFunctor(
38 std::size_t n_p,
39 MC::ParticlePositions position,
40 typename Model::SelfParticle model,
42 ParticlePropertyViewType<ComputeSpace> particle_values,
43 ParticlePropertyViewType<ComputeSpace> spatial_values,
44 ParticlePropertyViewType<ComputeSpace> ages_value,
45 MC::ParticleStatus status)
46 : n_p(n_p), position(std::move(position)), model(std::move(model)),
47 ages(std::move(ages)), particle_values(std::move(particle_values)),
48 spatial_values(std::move(spatial_values)),
49 ages_value(std::move(ages_value)), status(std::move(status))
50 {
51 }
52
53 KOKKOS_INLINE_FUNCTION
54 void
55 operator()(const int i_particle) const
57 {
58 if (status(i_particle) != MC::Status::Idle)
59 {
60 return;
61 }
62
63 auto access = scatter_spatial_values.access();
64 for (size_t k_var = 0; k_var < Model::n_var; ++k_var)
65 {
66 const auto current = model(i_particle, k_var);
67 access(k_var, position(i_particle)) += current;
68 particle_values(k_var, i_particle) = current;
69 }
70
71 const auto mass = Model::mass(i_particle, model);
72 access(Model::n_var, position(i_particle)) += mass;
73 particle_values(Model::n_var, i_particle) = mass;
74 ages_value(0, i_particle) = ages(i_particle, 0);
75 ages_value(1, i_particle) = ages(i_particle, 1);
76 }
77
78 KOKKOS_INLINE_FUNCTION
79 void
80 operator()(const int i_particle) const
82 {
83 if (status(i_particle) != MC::Status::Idle)
84 {
85 return;
86 }
87
88 auto access = scatter_spatial_values.access();
89 for (size_t k_var = 0; k_var < kindices.extent(0); ++k_var)
90 {
91 const auto current = model(i_particle, kindices(k_var));
92 access(k_var, position(i_particle)) += current;
93 particle_values(k_var, i_particle) = current;
94 }
95
96 const auto mass = Model::mass(i_particle, model);
97 access(kindices.extent(0), position(i_particle)) += mass;
98 particle_values(kindices.extent(0), i_particle) = mass;
99 ages_value(0, i_particle) = ages(i_particle, 0);
100 ages_value(1, i_particle) = ages(i_particle, 1);
101 }
102
103 void
104 run()
105 {
106
107 scatter_spatial_values
108 = Kokkos::Experimental::create_scatter_view(spatial_values);
109
110 // Use kindices for partial to map with correct vector position
111 // Warning UB if len(kindex)!=len(vnames)
113 {
114 static const auto indices = Model::get_number();
115
116 Kokkos::View<size_t*, HostSpace> host_index("host_index",
117 indices.size());
118 for (size_t i = 0; i < indices.size(); ++i)
119 {
120 host_index(i) = indices[i];
121 }
122 kindices
123 = Kokkos::create_mirror_view_and_copy(ComputeSpace(), host_index);
124 }
125
126 Kokkos::parallel_for("get_properties",
127 Kokkos::RangePolicy<ExecutionSpace>(0, n_p),
128 *this);
129
130 Kokkos::fence();
131 Kokkos::Experimental::contribute(spatial_values,
132 scatter_spatial_values);
133 }
134
135 std::size_t n_p;
136 MC::ParticlePositions position;
137 typename Model::SelfParticle model;
138 MC::ParticleAges ages;
139 ParticlePropertyViewType<ComputeSpace> particle_values;
140 ParticlePropertyViewType<ComputeSpace> spatial_values;
141 ParticlePropertyViewType<ComputeSpace> ages_value;
142 MC::ParticleStatus status;
143 Kokkos::View<const size_t*, ComputeSpace> kindices;
144 Kokkos::Experimental::
145 ScatterView<double**, Kokkos::LayoutRight, ComputeSpace>
146 scatter_spatial_values;
147 };
148 } // namespace
149
150 template <ModelType M>
151 ParticlePropertyViewType<ComputeSpace>
153 {
154 // USE list size not Kokkos View size.
155 // bcause container allocates more
156 // particles than needed
157 const std::size_t n_p = container.n_particles();
158 auto p_age = container.ages;
159 ParticlePropertyViewType<ComputeSpace> ages_values("ages_values", 2, n_p);
160
161 Kokkos::parallel_for(
162 "get_age",
163 Kokkos::RangePolicy<ComputeSpace>(0, n_p),
164 KOKKOS_LAMBDA(const std::size_t i_particle) {
165 ages_values(0, i_particle) = p_age(i_particle, 0);
166 ages_values(1, i_particle) = p_age(i_particle, 1);
167 });
168 return ages_values;
169 }
170
171 template <ModelType M>
172 std::optional<PostProcessing::BonceBuffer>
174 const std::size_t n_compartment,
175 const bool with_age)
176 {
177 if constexpr (HasExportProperties<M>)
178 {
179 container.force_remove_dead();
180 BonceBuffer properties;
181 properties.ages = std::nullopt;
182 const std::size_t n_p
183 = container.n_particles(); // USE list size not Kokkos View size.
184 // bcause container allocates more
185 // particles than needed
186 auto ar = M::names();
187 properties.vnames = std::vector<std::string>(ar.begin(), ar.end());
188 properties.vnames->emplace_back("mass");
189
190 std::size_t n_var = M::n_var + 1;
191 if constexpr (HasExportPropertiesPartial<M>)
192 {
193 n_var = ar.size() + 1;
194 if (ar.size() != M::get_number().size())
195 {
196 throw std::invalid_argument("Partial export Model need to have same "
197 "number of name and indices");
198 }
199 }
200
201 ParticlePropertyViewType<ComputeSpace> spatial_values(
202 "property_spatial", n_var, n_compartment);
203 ParticlePropertyViewType<ComputeSpace> particle_values(
204 "property_values", n_var, n_p);
205
206 // TODO Currently view is allocated and filled even if with_age flag is
207 // turned of, with_age controlls export but not postprocessing
208 ParticlePropertyViewType<ComputeSpace> ages_values("ages_values", 2, n_p);
209
210 GetPropertiesFunctor<M, Kokkos::DefaultExecutionSpace>(n_p,
211 container.position,
212 container.model,
213 container.ages,
214 particle_values,
215 spatial_values,
216 ages_values,
217 container.status)
218 .run();
219
220 properties.particle_values = Kokkos::create_mirror_view_and_copy(
221 Kokkos::HostSpace(), particle_values);
222 properties.spatial_values = Kokkos::create_mirror_view_and_copy(
223 Kokkos::HostSpace(), spatial_values);
224
225 if (with_age)
226 {
227 ParticlePropertyViewType<ComputeSpace> ages_values(
228 "ages_values", 2, n_p);
229
230 properties.ages = Kokkos::create_mirror_view_and_copy(
231 Kokkos::HostSpace(), ages_values);
232 }
233
234 return properties;
235 }
236 else
237 {
238
239 if (with_age)
240 {
241 BonceBuffer properties;
242 properties.ages = Kokkos::create_mirror_view_and_copy(
243 Kokkos::HostSpace(), get_particle_age_only(container));
244 return properties;
245 }
246 return std::nullopt;
247 }
248 }
249
250} // namespace PostProcessing
251#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:57
void force_remove_dead()
Clean all the non-idle particle even if number smaller than threshold.
Definition particles_container.hpp:464
Model::SelfParticle model
Definition particles_container.hpp:83
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Gets the number of particles in the container.
Definition particles_container.hpp:450
MC::ParticlePositions position
Definition particles_container.hpp:84
ParticleAges ages
Definition particles_container.hpp:87
MC::ParticleStatus status
Definition particles_container.hpp:85
SFNIAE way to check whether model allow all value saving.
Definition traits.hpp:152
SFNIAE way to check whether model allow partial value saving.
Definition traits.hpp:158
Model that can export properties.
Definition traits.hpp:166
Kokkos::View< Status *, ComputeSpace > ParticleStatus
Definition alias.hpp:74
Kokkos::View< uint64_t *, ComputeSpace > ParticlePositions
Definition alias.hpp:73
@ Idle
Definition alias.hpp:59
ParticleAgesBase< ComputeSpace > ParticleAges
Definition alias.hpp:76
Definition impl_post_process.hpp:21
std::optional< PostProcessing::BonceBuffer > get_properties(MC::ParticlesContainer< M > &container, const std::size_t n_compartment, const bool with_age)
Definition post_process.hpp:172
ParticlePropertyViewType< ComputeSpace > get_particle_age_only(MC::ParticlesContainer< M > &container)
Definition post_process.hpp:151
Definition post_process.hpp:22
std::optional< ParticlePropertyViewType< HostSpace > > ages
Definition post_process.hpp:25
std::optional< ParticlePropertyViewType< HostSpace > > particle_values
Definition post_process.hpp:23
std::optional< std::vector< std::string > > vnames
Definition post_process.hpp:26
std::optional< ParticlePropertyViewType< HostSpace > > spatial_values
Definition post_process.hpp:24