BioCMAMC-ST
kernels.hpp
1#ifndef __SIMULATION_KERNELS_HPP__
2#define __SIMULATION_KERNELS_HPP__
3
4#include "Kokkos_Printf.hpp"
5#include "common/common.hpp"
6#include <common/kokkos_getpolicy.hpp>
7#include <mc/unit.hpp>
8#include <simulation/kernels/model_kernel.hpp>
9#include <simulation/kernels/move_kernel.hpp>
11{
12 template <typename Space>
15 template <typename Space>
16 using move_reducer_view_type = Kokkos::View<std::size_t, Space>;
17
18 template <typename Space, ModelType Model> struct CycleFunctors
19 {
20 using FModel = Model;
21
24
25 ComputeSpace move_space{};
26
31
32 CycleFunctors() = default;
33
34 void update(double d_t,
37 {
38 bool enable_move = new_move.liquid_volume.size() > 1;
39 bool enable_leave = new_move.leaving_flow.size() != 0;
40
41 auto n_sample = enable_move ? 2 : 0;
42 n_sample += enable_leave ? 1 : 0;
43 // std::cout<<n_sample<<std::endl;
44 container.change_nsample(n_sample);
45
46 // FIXME: cycle_kernel: need to update container because we change:
47 // n_used_element counter which is size_t outside of kernel. As functor
48 // owns a copy of container this counter is not update between iterations.
49 // 1. Use n_used_element as a reference counter type
50 // 2. Manually update counter in update function (dirty way)
51
52 cycle_kernel.update(d_t, container);
53
54 // TODO: Why need to update all views (where did we lost the refcount ? )
56 d_t,
57 container.n_particles(),
58 std::move(new_move),
59 container.position,
60 container.status,
61 container.ages,
62 container.random,
65 }
66
68 {
69 const auto host_red =
70 Kokkos::create_mirror_view_and_copy(HostSpace(), cycle_reducer)();
71
72 const auto host_out_counter =
73 Kokkos::create_mirror_view_and_copy(HostSpace(), move_reducer)();
74
75 return std::tuple(host_red, host_out_counter);
76 }
77
79 MC::pool_type _random_pool,
80 MC::KernelConcentrationType&& _concentrations,
81 MC::ContributionView _contribs_scatter,
82 MC::EventContainer _event,
85 : cycle_reducer("cycle_reducer"), move_reducer("move_reducer"),
86 cycle_kernel(container,
87 _random_pool,
88 std::move(_concentrations),
89 _contribs_scatter,
90 _event),
91 move_kernel(container.position,
92 container.status,
93 std::move(m),
94 _random_pool,
95 _event,
96 _probes,
97 container.ages,
98 container.random)
99 {
100 }
101
102 void launch_move(const std::size_t n_particle) const
103 {
104
105 if (move_kernel.enable_move)
106 {
107 // const auto _policy_move =
108 // Common::get_policy<move_kernel_type, KernelInline::TagMove>(
109 // move_kernel, n_particle, false);
110 const auto _policy_move =
112
113 Kokkos ::parallel_for("cycle_move", _policy_move, move_kernel);
114 }
115 if (move_kernel.enable_leave)
116 {
117
118 // const auto _policy_leave =
119 // Kokkos::TeamPolicy<ComputeSpace, KernelInline::TagLeave>(
120 // 256, Kokkos::AUTO);
121
122 // Kokkos ::parallel_reduce(
123 // "cycle_move_leave", _policy_leave, move_kernel, move_reducer);
124
125 const auto _policy_leave =
126 Kokkos::RangePolicy<KernelInline::TagLeave>(0, n_particle);
127 // const auto _policy_move = Kokkos::TeamPolicy<TagLeave>(56, 256);
128 Kokkos ::parallel_reduce(
129 "cycle_move_leave", _policy_leave, move_kernel, move_reducer);
130 }
131 }
132
133 void launch_model(const std::size_t n_particle) const
134 {
135 constexpr bool is_reduce = true;
136 // const auto _policy =
137 // Common::get_policy<cycle_kernel_type, KernelInline::TagSecondPass>(
138 // cycle_kernel, n_particle, is_reduce);
139
140 // auto _policy =
141 // Common::get_policy<KernelInline::TagCycle>(n_particle, is_reduce);
142
143 // _policy.set_scratch_size(
144 // 0,
145 // Kokkos::PerTeam(cycle_kernel.concentrations.extent(0) *
146 // _policy.team_size() * sizeof(double)));
147
148 // const auto scatter_policy =
149 // Kokkos::RangePolicy<KernelInline::TagContribution>(0, n_particle);
150
151 // auto _policy = Kokkos::RangePolicy<KernelInline::TagCycle>(0,
152 // n_particle);
153
154 // auto _policy = Kokkos::TeamPolicy<TagCycle>(53, 256);
155 // Kokkos::parallel_reduce(
156 // "cycle_model",
157 // _policy,
158 // cycle_kernel,
159 // KernelInline::CycleReducer<ComputeSpace>(cycle_reducer));
160 //
161
162 auto _policy = Kokkos::RangePolicy<TagCycle>(0, n_particle);
163 Kokkos::parallel_reduce(
164 "cycle_model",
165 _policy,
168
169 // Kokkos::parallel_reduce(
170 // "cycle_model",
171 // _policy,
172 // cycle_kernel,
173 // KernelInline::CycleReducer<ComputeSpace>(cycle_reducer));
174 // Kokkos::fence(); // TODO needed ?
175 // Kokkos::parallel_for("cycle_scatter", scatter_policy, cycle_kernel);
176 }
177 };
178
179} // namespace Simulation::KernelInline
180
181#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:53
void change_nsample(const std::size_t new_n_sample)
Definition particles_container.hpp:172
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Gets the number of particles in the container.
Definition particles_container.hpp:464
ParticleSamples random
Definition particles_container.hpp:84
MC::ParticlePositions position
Definition particles_container.hpp:80
ParticleAges ages
Definition particles_container.hpp:83
MC::ParticleStatus status
Definition particles_container.hpp:81
Definition model_kernel.hpp:43
Kokkos::View< value_type, Space > result_view_type
Definition model_kernel.hpp:48
Kokkos::TeamPolicy< ComputeSpace, Tag > get_policy(const FunctorType &f, std::size_t range, bool reduce=false)
Definition kokkos_getpolicy.hpp:37
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:38
decltype(Kokkos::Experimental::create_scatter_view(kernelContribution())) ContributionView
Definition alias.hpp:87
Kokkos::View< const double **, Kokkos::LayoutLeft, ComputeSpace, Kokkos::MemoryTraits< Kokkos::RandomAccess > > KernelConcentrationType
Definition alias.hpp:90
Definition kernels.hpp:11
constexpr bool enable_leave
Definition move_kernel.hpp:26
constexpr bool enable_move
Definition move_kernel.hpp:29
KernelInline::CycleReducer< Space >::result_view_type cycle_reducer_view_type
Definition kernels.hpp:13
Kokkos::View< std::size_t, Space > move_reducer_view_type
Definition kernels.hpp:16
Probes< AutoGenerated::probe_buffer_size > ProbeAutogeneratedBuffer
Definition probe.hpp:61
Use to count events that occurs during Monte-Carlo processing cycles.
Definition events.hpp:44
Definition model_kernel.hpp:106
CycleFunctor< Model > cycle_kernel_type
Definition kernels.hpp:22
void update(double d_t, MC::ParticlesContainer< Model > container, MoveInfo< ComputeSpace > new_move)
Definition kernels.hpp:34
cycle_kernel_type cycle_kernel
Definition kernels.hpp:29
MoveFunctor move_kernel_type
Definition kernels.hpp:23
CycleFunctors(MC::ParticlesContainer< Model > container, MC::pool_type _random_pool, MC::KernelConcentrationType &&_concentrations, MC::ContributionView _contribs_scatter, MC::EventContainer _event, MoveInfo< ComputeSpace > m, ProbeAutogeneratedBuffer _probes)
Definition kernels.hpp:78
move_kernel_type move_kernel
Definition kernels.hpp:30
auto get_host_reduction()
Definition kernels.hpp:67
void launch_move(const std::size_t n_particle) const
Definition kernels.hpp:102
ComputeSpace move_space
Definition kernels.hpp:25
cycle_reducer_view_type< Space > cycle_reducer
Definition kernels.hpp:27
move_reducer_view_type< Space > move_reducer
Definition kernels.hpp:28
void launch_model(const std::size_t n_particle) const
Definition kernels.hpp:133
Model FModel
Definition kernels.hpp:20
Definition move_kernel.hpp:133
Definition move_info.hpp:19
Kokkos::View< LeavingFlow *, Kokkos::SharedHostPinnedSpace, Kokkos::MemoryTraits< Kokkos::MemoryTraitsFlags::Aligned|Kokkos::MemoryTraitsFlags::Restrict > > leaving_flow
Definition move_info.hpp:29
Kokkos::View< double *, ExecSpace > liquid_volume
Definition move_info.hpp:30