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/domain.hpp>
8#include <mc/unit.hpp>
9#include <simulation/kernels/contribution_kernel.hpp>
10#include <simulation/kernels/model_kernel.hpp>
11#include <simulation/kernels/move_kernel.hpp>
12
13#include <common/execinfo.hpp>
14
16{
17
18 template <typename Space>
21 template <typename Space>
22 using move_reducer_view_type = Kokkos::View<std::size_t, Space>;
24 template <typename Space, ModelType Model> struct CycleFunctors
25 {
26 using FModel = Model;
31 ComputeSpace move_space;
33 ComputeSpace model_space;
42 CycleFunctors() = default;
47
48 void
49 update(const double d_t,
52 {
53 f_multi_compartment = new_move.liquid_volume.size() > 1;
55 const bool enable_leave = new_move.leaving_flow.size() != 0;
56
57 // FIXME: cycle_kernel: need to update container because we change:
58 // n_used_element counter which is size_t outside of kernel. As functor
59 // owns a copy of container this counter is not update between iterations.
60 // 1. Use n_used_element as a reference counter type
61 // 2. Manually update counter in update function (dirty way)
62
63 cycle_kernel.update(d_t, container);
64
65 contribution_kernel.update(container);
66
67 // TODO: Why need to update all views (where did we lost the refcount ? )
68 move_kernel.update(d_t,
69 container.n_particles(),
70 std::move(new_move),
71 container.position,
72 container.status,
73 container.ages,
76 }
77
78 auto
80 {
81 const auto host_red
82 = Kokkos::create_mirror_view_and_copy(HostSpace(), cycle_reducer)();
83
84 const auto host_out_counter
85 = Kokkos::create_mirror_view_and_copy(HostSpace(), move_reducer)();
86
87 return std::tuple(host_red, host_out_counter);
88 }
92 MC::pool_type _random_pool,
93 MC::KernelConcentrationType _concentrations,
94 MC::ContributionView _contribs_scatter,
95 MC::EventContainer _event,
98 ProbeAutogeneratedBuffer _probes_div)
99 : cycle_reducer("cycle_reducer"), move_reducer("move_reducer"),
100 cycle_kernel(options.m_p_p_team_model,
101 container,
102 _random_pool,
103 std::move(_concentrations),
104 _event,
105 _probes_div),
106 move_kernel(options.m_p_p_team_move,
107 options.m_p_p_team_leave,
108 container.position,
109 container.status,
110 std::move(m),
111 _random_pool,
112 _event,
113 _probes,
114 container.ages),
116 options.m_p_p_team_contribs, _contribs_scatter, container),
117 m_options(options)
118
119 {
120 }
121
122 void
123 launch_move(const std::size_t n_particle) const
124 {
125
126 if (move_kernel.enable_move)
127 {
128 const auto npt = m_options.m_p_p_team_move;
129
130 const std::size_t league_size = Common::c_league_size(n_particle, npt);
131
132 auto cycle_policy
133 = Kokkos::TeamPolicy<TagMove>(model_space,
134 static_cast<int>(league_size),
135 Kokkos::AUTO(),
136 Kokkos::AUTO());
137
138 cycle_policy.set_scratch_size(
139 0, Kokkos::PerTeam(sizeof(float_t) * npt * 2));
140
141 Kokkos ::parallel_for("cycle_move", cycle_policy, move_kernel);
142 }
143
144 if (move_kernel.enable_leave)
145 {
146
147 const auto _policy_leave = Kokkos::RangePolicy<KernelInline::TagLeave>(
148 move_space, 0, n_particle);
149 Kokkos ::parallel_reduce(
150 "cycle_move_leave", _policy_leave, move_kernel, move_reducer);
151 }
152 }
153
154 void
155 launch_model(const std::size_t n_particle) const
156 {
157 std::size_t league_size
158 = Common::c_league_size(n_particle, m_options.m_p_p_team_model);
159
160 const auto cycle_policy
161 = Kokkos::TeamPolicy<TagCycle>(model_space,
162 static_cast<int>(league_size),
163 Kokkos::AUTO(),
164 Kokkos::AUTO());
165
166 Kokkos::parallel_reduce(
167 "cycle_model",
168 cycle_policy,
171 Kokkos::fence(); // TODO needed ?
172
173 // Assumptions
174 // Newborn cells do not contribte in current time step
175 // Mother cell doesn´t exist but
176 // Contribution array is not changed during division and
177
178 if (cycle_kernel.do_contribs())
179 {
180 league_size
181 = Common::c_league_size(n_particle, m_options.m_p_p_team_contribs);
182 static_assert(ConstWeightModelType<Model>,
183 "ModelType:Constapply_weight()");
184
186 {
187
188 const auto policy_contribs
189 = Kokkos::TeamPolicy<typename ContributionFunctor<Model>::Tag3D>(
190 model_space, league_size, Kokkos::AUTO(), Kokkos::AUTO());
191 Kokkos::parallel_for(
192 "cycle_model_contribs", policy_contribs, contribution_kernel);
193 }
194 else
195 {
196 auto policy_contribs
197 = Kokkos::TeamPolicy<typename ContributionFunctor<Model>::Tag0D>(
198 model_space, league_size, Kokkos::AUTO(), Kokkos::AUTO());
199 policy_contribs.set_scratch_size(
200 0, Kokkos::PerTeam(sizeof(float_t) * Model::n_c));
201 Kokkos::parallel_for(
202 "cycle_model_contribs_0d", policy_contribs, contribution_kernel);
203 }
204 }
205 }
206 };
207
208} // namespace Simulation::KernelInline
209
210#endif
Main owning object for Monte-Carlo particles.
Definition particles_container.hpp:57
KOKKOS_INLINE_FUNCTION std::size_t n_particles() const
Gets the number of particles in the container.
Definition particles_container.hpp:449
MC::ParticlePositions position
Definition particles_container.hpp:85
ParticleAges ages
Definition particles_container.hpp:88
MC::ParticleStatus status
Definition particles_container.hpp:86
Definition model_kernel.hpp:55
Kokkos::View< value_type, Space > result_view_type
Definition model_kernel.hpp:60
Concept to check if a model type has uniform_weight
Definition traits.hpp:189
std::size_t c_league_size(std::size_t n_tot, std::size_t n_per_team) noexcept
Definition common.cpp:17
decltype(Kokkos::Experimental::create_scatter_view( kernelContribution())) ContributionView
Definition alias.hpp:162
gen_pool_type< Kokkos::DefaultExecutionSpace > pool_type
Definition alias.hpp:100
Kokkos::View< const double **, Kokkos::LayoutLeft, ComputeSpace, Kokkos::MemoryTraits< Kokkos::RandomAccess > > KernelConcentrationType
Definition alias.hpp:165
Definition kernels.hpp:16
constexpr bool enable_leave
Definition move_kernel.hpp:22
KernelInline::CycleReducer< Space >::result_view_type cycle_reducer_view_type
Definition kernels.hpp:19
constexpr bool enable_move
Definition move_kernel.hpp:25
Kokkos::View< std::size_t, Space > move_reducer_view_type
Definition kernels.hpp:21
Probes< AutoGenerated::probe_buffer_size > ProbeAutogeneratedBuffer
Definition probe.hpp:148
Definition contribution_kernel.hpp:9
Definition execinfo.hpp:12
Structure to store information about domain needed during MC cycle data is likely to change between e...
Definition domain.hpp:28
Use to count events that occurs during Monte-Carlo processing cycles.
Definition events.hpp:150
Definition model_kernel.hpp:124
CycleFunctor< Model > cycle_kernel_type
Definition kernels.hpp:27
cycle_kernel_type cycle_kernel
Definition kernels.hpp:36
ComputeSpace model_space
Definition kernels.hpp:32
MoveFunctor move_kernel_type
Definition kernels.hpp:28
bool f_multi_compartment
Definition kernels.hpp:45
move_kernel_type move_kernel
Definition kernels.hpp:37
auto get_host_reduction()
Definition kernels.hpp:78
void launch_move(const std::size_t n_particle) const
Definition kernels.hpp:122
KernelDispatchOptions m_options
Definition kernels.hpp:43
ComputeSpace move_space
Definition kernels.hpp:30
cycle_reducer_view_type< Space > cycle_reducer
Definition kernels.hpp:34
move_reducer_view_type< Space > move_reducer
Definition kernels.hpp:35
void launch_model(const std::size_t n_particle) const
Definition kernels.hpp:154
Model FModel
Definition kernels.hpp:25
ContributionFunctor< Model > contribution_kernel
Definition kernels.hpp:39
void update(const double d_t, MC::ParticlesContainer< Model > container, MC::DomainState< ComputeSpace > &&new_move)
Definition kernels.hpp:48
Definition move_kernel.hpp:115