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