BioCMAMC-ST
post_process.hpp
1#ifndef __CORE_POST_PROCESS_PUBLIC_HPP__
2#define __CORE_POST_PROCESS_PUBLIC_HPP__
3#include "Kokkos_Atomic.hpp"
4#include "Kokkos_Core_fwd.hpp"
5#include "common/common.hpp"
6#include "mc/alias.hpp"
7#include "mc/traits.hpp"
8#include <Kokkos_Core.hpp>
9#include <cstdint>
10#include <mc/particles_container.hpp>
11
12template <typename MemorySpace>
13using ParticlePropertyViewType =
14 Kokkos::View<double**, Kokkos::LayoutRight, MemorySpace>;
15
16using SubViewtype = Kokkos::Subview<ParticlePropertyViewType<ComputeSpace>,
17 decltype(Kokkos::ALL),
18 std::size_t>;
19namespace PostProcessing
20{
22 {
23 ParticlePropertyViewType<HostSpace> particle_values;
24 ParticlePropertyViewType<HostSpace> spatial_values;
25 std::vector<std::string> vnames;
26 };
27
28 template <HasExportPropertiesFull Model, typename ExecutionSpace>
29 void inner(std::size_t n_p,
30 const MC::ParticlePositions& position,
31 const typename Model::SelfParticle& model,
32 ParticlePropertyViewType<ComputeSpace>& particle_values,
33 ParticlePropertyViewType<ComputeSpace>& spatial_values,
34 const MC::ParticleStatus& status)
35 {
36 Kokkos::Experimental::
37 ScatterView<double**, Kokkos::LayoutRight, ComputeSpace>
38 scatter_spatial_values(spatial_values);
39
40 Kokkos::parallel_for(
41 "kernel_get_properties",
42 Kokkos::RangePolicy<ExecutionSpace>(0, n_p),
43 KOKKOS_LAMBDA(const int i_particle) {
44 if (status(i_particle) != MC::Status::Idle)
45 {
46 return;
47 }
48 auto access = scatter_spatial_values.access();
49
50 for (std::size_t i = 0; i < Model::n_var; ++i)
51 {
52 access(i, position(i_particle)) += model(i_particle, i);
53 particle_values(i, i_particle) = model(i_particle, i);
54 }
55 auto mass = Model::mass(i_particle, model);
56 access(Model::n_var, position(i_particle)) += mass;
57 particle_values(Model::n_var, i_particle) = mass;
58 });
59 Kokkos::fence();
60 Kokkos::Experimental::contribute(spatial_values, scatter_spatial_values);
61 }
62
63 template <HasExportPropertiesPartial Model, typename ExecutionSpace>
64 void inner_partial(const std::size_t n_p,
65 const MC::ParticlePositions& position,
66 const typename Model::SelfParticle& model,
67 ParticlePropertyViewType<ComputeSpace>& particle_values,
68 ParticlePropertyViewType<ComputeSpace>& spatial_values,
69 const MC::ParticleStatus& status)
70 {
71 Kokkos::Experimental::
72 ScatterView<double**, Kokkos::LayoutRight, ComputeSpace>
73 scatter_spatial_values(spatial_values);
74 static const auto indices = Model::get_number();
75
76 static Kokkos::View<size_t*, HostSpace> host_index("host_index",
77 indices.size());
78
79 for (std::size_t i = 0; i < indices.size(); ++i)
80 {
81 host_index(i) = indices[i];
82 }
83
84 const Kokkos::View<const size_t*, ComputeSpace> kindices =
85 Kokkos::create_mirror_view_and_copy(ComputeSpace(), host_index);
86 Kokkos::parallel_for(
87 "kernel_get_properties",
88 Kokkos::RangePolicy<ExecutionSpace>(0, n_p),
89 KOKKOS_LAMBDA(const int i_particle) {
90 if (status(i_particle) != MC::Status::Idle)
91 {
92 return;
93 }
94
95 auto access = scatter_spatial_values.access();
96
97 for (size_t k_var = 0; k_var < kindices.extent(0); ++k_var)
98 {
99 const auto index_export = kindices(k_var);
100 const auto current = model(i_particle, index_export);
101
102 access(k_var, position(i_particle)) += current;
103 particle_values(k_var, i_particle) = current;
104 }
105 const auto mass = Model::mass(i_particle, model);
106 access(kindices.extent(0), position(i_particle)) += mass;
107 particle_values(kindices.extent(0), i_particle) = mass;
108 });
109 Kokkos::fence();
110 Kokkos::Experimental::contribute(spatial_values, scatter_spatial_values);
111 }
112
113 template <ModelType M>
114 std::optional<PostProcessing::BonceBuffer>
116 const std::size_t n_compartment)
117 {
118 if constexpr (HasExportProperties<M>)
119 {
120 container.force_remove_dead();
121 BonceBuffer properties;
122 const std::size_t n_p =
123 container.n_particles(); // USE list size not Kokkos View size.
124 // ParticleList allocates more particles than
125 // needed
126 auto ar = M::names();
127 properties.vnames = std::vector<std::string>(ar.begin(), ar.end());
128 properties.vnames.emplace_back("mass");
129
130 std::size_t n_var = M::n_var + 1;
131 if constexpr (HasExportPropertiesPartial<M>)
132 {
133 n_var = ar.size() + 1;
134 }
135 ParticlePropertyViewType<ComputeSpace> spatial_values(
136 "property_spatial", n_var, n_compartment);
137 ParticlePropertyViewType<ComputeSpace> particle_values(
138 "property_values", n_var, n_p);
139
140 if constexpr (HasExportPropertiesPartial<M>)
141 {
143 container.position,
144 container.model,
145 particle_values,
146 spatial_values,
147 container.status);
148 }
149 else
150 {
152 container.position,
153 container.model,
154 particle_values,
155 spatial_values,
156 container.status);
157 }
158
159 properties.particle_values = Kokkos::create_mirror_view_and_copy(
160 Kokkos::HostSpace(), particle_values);
161 properties.spatial_values = Kokkos::create_mirror_view_and_copy(
162 Kokkos::HostSpace(), spatial_values);
163 return properties;
164 }
165 else
166 {
167 return std::nullopt;
168 }
169 }
170
171} // namespace PostProcessing
172#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:143
void force_remove_dead()
Definition particles_container.hpp:304
Model::SelfParticle model
Definition particles_container.hpp:157
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Gets the number of particles in the container.
Definition particles_container.hpp:293
MC::ParticlePositions position
Definition particles_container.hpp:158
MC::ParticleStatus status
Definition particles_container.hpp:159
SFNIAE way to check whether model allow partial value saving.
Definition traits.hpp:170
Model that can export properties.
Definition traits.hpp:177
Kokkos::View< Status *, ComputeSpace, restrict_mt > ParticleStatus
Definition alias.hpp:43
Kokkos::View< uint64_t *, ComputeSpace, restrict_mt > ParticlePositions
Definition alias.hpp:42
@ Idle
Definition alias.hpp:25
Definition impl_post_process.hpp:21
void inner(std::size_t n_p, const MC::ParticlePositions &position, const typename Model::SelfParticle &model, ParticlePropertyViewType< ComputeSpace > &particle_values, ParticlePropertyViewType< ComputeSpace > &spatial_values, const MC::ParticleStatus &status)
Definition post_process.hpp:29
std::optional< PostProcessing::BonceBuffer > get_properties(MC::ParticlesContainer< M > &container, const std::size_t n_compartment)
Definition post_process.hpp:115
void inner_partial(const std::size_t n_p, const MC::ParticlePositions &position, const typename Model::SelfParticle &model, ParticlePropertyViewType< ComputeSpace > &particle_values, ParticlePropertyViewType< ComputeSpace > &spatial_values, const MC::ParticleStatus &status)
Definition post_process.hpp:64
Definition post_process.hpp:22
ParticlePropertyViewType< HostSpace > spatial_values
Definition post_process.hpp:24
std::vector< std::string > vnames
Definition post_process.hpp:25
ParticlePropertyViewType< HostSpace > particle_values
Definition post_process.hpp:23