BioCMAMC-ST
kernels.hpp
1#ifndef __SIMULATION_KERNELS_HPP__
2#define __SIMULATION_KERNELS_HPP__
3
4#include <Kokkos_Core_fwd.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>
10
12{
13 template <typename Space>
16 template <typename Space>
17 using move_reducer_view_type = Kokkos::View<std::size_t, Space>;
19 template <typename Space, ModelType Model> struct CycleFunctors
20 {
21 using FModel = Model;
26 ComputeSpace move_space;
28 ComputeSpace model_space;
35 CycleFunctors() = default;
36
37 void
38 update(const double d_t,
41 {
42 const bool enable_move = new_move.liquid_volume.size() > 1;
43 const bool enable_leave = new_move.leaving_flow.size() != 0;
44
45 auto n_sample = enable_move ? 2 : 0;
46 n_sample += enable_leave ? 1 : 0;
47 container.change_nsample(n_sample);
48
49 // FIXME: cycle_kernel: need to update container because we change:
50 // n_used_element counter which is size_t outside of kernel. As functor
51 // owns a copy of container this counter is not update between iterations.
52 // 1. Use n_used_element as a reference counter type
53 // 2. Manually update counter in update function (dirty way)
54
55 cycle_kernel.update(d_t, container);
56
57 // TODO: Why need to update all views (where did we lost the refcount ? )
59 d_t,
60 container.n_particles(),
61 std::move(new_move),
62 container.position,
63 container.status,
64 container.ages,
65 container.random,
68 }
69
70 auto
72 {
73 const auto host_red
74 = Kokkos::create_mirror_view_and_copy(HostSpace(), cycle_reducer)();
75
76 const auto host_out_counter
77 = Kokkos::create_mirror_view_and_copy(HostSpace(), move_reducer)();
78
79 return std::tuple(host_red, host_out_counter);
80 }
83 MC::pool_type _random_pool,
84 MC::KernelConcentrationType&& _concentrations,
85 MC::ContributionView _contribs_scatter,
86 MC::EventContainer _event,
89 ProbeAutogeneratedBuffer _probes_div)
90 : cycle_reducer("cycle_reducer"), move_reducer("move_reducer"),
91 cycle_kernel(container,
92 _random_pool,
93 std::move(_concentrations),
94 _contribs_scatter,
95 _event,
96 _probes_div),
97 move_kernel(container.position,
98 container.status,
99 std::move(m),
100 _random_pool,
101 _event,
102 _probes,
103 container.ages,
104 container.random)
105 {
106 }
107
108 void
109 launch_move(const std::size_t n_particle) const
110 {
111
112 if (move_kernel.enable_move)
113 {
114 // const auto _policy_move =
115 // Common::get_policy<move_kernel_type, KernelInline::TagMove>(
116 // move_kernel, n_particle, false);
117 // const auto _policy_move =
118 // Common::get_policy<KernelInline::TagMove>(n_particle, false);
119
120 const auto _policy_move = Kokkos::RangePolicy<KernelInline::TagMove>(
121 move_space, 0, n_particle);
122
123 Kokkos ::parallel_for("cycle_move", _policy_move, move_kernel);
124 }
125 if (move_kernel.enable_leave)
126 {
127
128 // const auto _policy_leave =
129 // Kokkos::TeamPolicy<ComputeSpace, KernelInline::TagLeave>(
130 // 256, Kokkos::AUTO);
131
132 // Kokkos ::parallel_reduce(
133 // "cycle_move_leave", _policy_leave, move_kernel, move_reducer);
134
135 const auto _policy_leave = Kokkos::RangePolicy<KernelInline::TagLeave>(
136 move_space, 0, n_particle);
137 // const auto _policy_move = Kokkos::TeamPolicy<TagLeave>(56, 256);
138 Kokkos ::parallel_reduce(
139 "cycle_move_leave", _policy_leave, move_kernel, move_reducer);
140 }
141 }
142
143 void
144 launch_model(const std::size_t n_particle) const
145 {
146
147 auto _policy = Kokkos::RangePolicy<TagCycle>(model_space, 0, n_particle);
148 Kokkos::parallel_reduce(
149 "cycle_model",
150 _policy,
153 Kokkos::fence(); // TODO needed ?
154
155 constexpr int PARTICLES_PER_TEAM = 256;
156 int league_size = Kokkos::ceil(n_particle / PARTICLES_PER_TEAM);
157
158 auto _policy2 = Kokkos::TeamPolicy<TagContribution>(
159 model_space, league_size, Kokkos::AUTO(), Kokkos::AUTO());
160 Kokkos::parallel_for("cycle_model2", _policy2, cycle_kernel);
161 }
162 };
163
164} // namespace Simulation::KernelInline
165
166#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:56
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Gets the number of particles in the container.
Definition particles_container.hpp:445
ParticleSamples random
Definition particles_container.hpp:87
void change_nsample(std::size_t new_n_sample)
Definition particles_container.hpp:540
MC::ParticlePositions position
Definition particles_container.hpp:83
ParticleAges ages
Definition particles_container.hpp:86
MC::ParticleStatus status
Definition particles_container.hpp:84
Definition model_kernel.hpp:43
Kokkos::View< value_type, Space > result_view_type
Definition model_kernel.hpp:48
decltype(Kokkos::Experimental::create_scatter_view( kernelContribution())) ContributionView
Definition alias.hpp:88
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:39
Kokkos::View< const double **, Kokkos::LayoutLeft, ComputeSpace, Kokkos::MemoryTraits< Kokkos::RandomAccess > > KernelConcentrationType
Definition alias.hpp:91
Definition kernels.hpp:12
constexpr bool enable_leave
Definition move_kernel.hpp:25
KernelInline::CycleReducer< Space >::result_view_type cycle_reducer_view_type
Definition kernels.hpp:14
constexpr bool enable_move
Definition move_kernel.hpp:28
Kokkos::View< std::size_t, Space > move_reducer_view_type
Definition kernels.hpp:16
Probes< AutoGenerated::probe_buffer_size > ProbeAutogeneratedBuffer
Definition probe.hpp:55
Structure to store information about domain needed during MC cycle data is likely to change between e...
Definition domain.hpp:27
Use to count events that occurs during Monte-Carlo processing cycles.
Definition events.hpp:45
Definition model_kernel.hpp:112
CycleFunctor< Model > cycle_kernel_type
Definition kernels.hpp:22
cycle_kernel_type cycle_kernel
Definition kernels.hpp:31
ComputeSpace model_space
Definition kernels.hpp:27
MoveFunctor move_kernel_type
Definition kernels.hpp:23
move_kernel_type move_kernel
Definition kernels.hpp:32
auto get_host_reduction()
Definition kernels.hpp:70
void launch_move(const std::size_t n_particle) const
Definition kernels.hpp:108
ComputeSpace move_space
Definition kernels.hpp:25
cycle_reducer_view_type< Space > cycle_reducer
Definition kernels.hpp:29
move_reducer_view_type< Space > move_reducer
Definition kernels.hpp:30
void launch_model(const std::size_t n_particle) const
Definition kernels.hpp:143
Model FModel
Definition kernels.hpp:20
void update(const double d_t, MC::ParticlesContainer< Model > container, MC::DomainState< ComputeSpace > &&new_move)
Definition kernels.hpp:37
Definition move_kernel.hpp:134