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>
9#include <mc/particles_container.hpp>
10#include <mc/traits.hpp>
13template <
typename MemorySpace>
14using ParticlePropertyViewType
15 = Kokkos::View<double**, Kokkos::LayoutRight, MemorySpace>;
17using SubViewtype = Kokkos::Subview<ParticlePropertyViewType<ComputeSpace>,
18 decltype(Kokkos::ALL),
26 std::optional<ParticlePropertyViewType<HostSpace>>
ages;
27 std::optional<std::vector<std::string>>
vnames;
32 template <
typename Model,
typename ExecutionSpace>
33 struct GetPropertiesFunctor
35 using ComputeSpace =
typename ExecutionSpace::memory_space;
40 typename Model::SelfParticle model,
42 ParticlePropertyViewType<ComputeSpace> particle_values,
43 ParticlePropertyViewType<ComputeSpace> spatial_values,
44 ParticlePropertyViewType<ComputeSpace> ages_value,
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))
53 KOKKOS_INLINE_FUNCTION
55 operator()(
const int i_particle)
const
63 auto access = scatter_spatial_values.access();
64 for (
size_t k_var = 0; k_var < Model::n_var; ++k_var)
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;
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);
78 KOKKOS_INLINE_FUNCTION
80 operator()(
const int i_particle)
const
88 auto access = scatter_spatial_values.access();
89 for (
size_t k_var = 0; k_var < kindices.extent(0); ++k_var)
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;
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);
107 scatter_spatial_values
108 = Kokkos::Experimental::create_scatter_view(spatial_values);
114 static const auto indices = Model::get_number();
116 Kokkos::View<size_t*, HostSpace> host_index(
"host_index",
118 for (
size_t i = 0; i < indices.size(); ++i)
120 host_index(i) = indices[i];
123 = Kokkos::create_mirror_view_and_copy(ComputeSpace(), host_index);
126 Kokkos::parallel_for(
"get_properties",
127 Kokkos::RangePolicy<ExecutionSpace>(0, n_p),
131 Kokkos::Experimental::contribute(spatial_values,
132 scatter_spatial_values);
137 typename Model::SelfParticle model;
139 ParticlePropertyViewType<ComputeSpace> particle_values;
140 ParticlePropertyViewType<ComputeSpace> spatial_values;
141 ParticlePropertyViewType<ComputeSpace> ages_value;
143 Kokkos::View<const size_t*, ComputeSpace> kindices;
144 Kokkos::Experimental::
145 ScatterView<double**, Kokkos::LayoutRight, ComputeSpace>
146 scatter_spatial_values;
150 template <ModelType M>
151 ParticlePropertyViewType<ComputeSpace>
158 auto p_age = container.
ages;
159 ParticlePropertyViewType<ComputeSpace> ages_values(
"ages_values", 2, n_p);
161 Kokkos::parallel_for(
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);
171 template <ModelType M>
172 std::optional<PostProcessing::BonceBuffer>
174 const std::size_t n_compartment,
181 properties.
ages = std::nullopt;
182 const std::size_t n_p
186 auto ar = M::names();
187 properties.
vnames = std::vector<std::string>(ar.begin(), ar.end());
188 properties.
vnames->emplace_back(
"mass");
190 std::size_t n_var = M::n_var + 1;
193 n_var = ar.size() + 1;
194 if (ar.size() != M::get_number().size())
196 throw std::invalid_argument(
"Partial export Model need to have same "
197 "number of name and indices");
201 ParticlePropertyViewType<ComputeSpace> spatial_values(
202 "property_spatial", n_var, n_compartment);
203 ParticlePropertyViewType<ComputeSpace> particle_values(
204 "property_values", n_var, n_p);
208 ParticlePropertyViewType<ComputeSpace> ages_values(
"ages_values", 2, n_p);
210 GetPropertiesFunctor<M, Kokkos::DefaultExecutionSpace>(n_p,
221 Kokkos::HostSpace(), particle_values);
223 Kokkos::HostSpace(), spatial_values);
227 ParticlePropertyViewType<ComputeSpace> ages_values(
228 "ages_values", 2, n_p);
230 properties.
ages = Kokkos::create_mirror_view_and_copy(
231 Kokkos::HostSpace(), ages_values);
242 properties.
ages = Kokkos::create_mirror_view_and_copy(
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