IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ParticleAttrib.hpp
Go to the documentation of this file.
1//
2// Class ParticleAttrib
3// Templated class for all particle attribute classes.
4//
5// This templated class is used to represent a single particle attribute.
6// An attribute is one data element within a particle object, and is
7// stored as a Kokkos::View. This class stores the type information for the
8// attribute, and provides methods to create and destroy new items, and
9// to perform operations involving this attribute with others.
10//
11// ParticleAttrib is the primary element involved in expressions for
12// particles (just as BareField is the primary element there). This file
13// defines the necessary templated classes and functions to make
14// ParticleAttrib a capable expression-template participant.
15//
16#include "Ippl.h"
17
19
20#include "Utility/IpplTimings.h"
21
22namespace ippl {
23
24 template <typename T, class... Properties>
26 size_type required = *(this->localNum_mp) + n;
27 if (this->size() < required) {
28 int overalloc = Comm->getDefaultOverallocation();
29 this->realloc(required * overalloc);
30 }
31 }
32
33 template <typename T, class... Properties>
35 const hash_type& keepIndex,
36 size_type invalidCount) {
37 // Replace all invalid particles in the valid region with valid
38 // particles in the invalid region
39 using policy_type = Kokkos::RangePolicy<execution_space>;
40 Kokkos::parallel_for(
41 "ParticleAttrib::destroy()", policy_type(0, invalidCount),
42 KOKKOS_CLASS_LAMBDA(const size_t i) {
43 dview_m(deleteIndex(i)) = dview_m(keepIndex(i));
44 });
45 }
46
47 template <typename T, class... Properties>
49 auto size = hash.extent(0);
50 if (buf_m.extent(0) < size) {
51 int overalloc = Comm->getDefaultOverallocation();
52 Kokkos::realloc(buf_m, size * overalloc);
53 }
54
55 using policy_type = Kokkos::RangePolicy<execution_space>;
56 Kokkos::parallel_for(
57 "ParticleAttrib::pack()", policy_type(0, size),
58 KOKKOS_CLASS_LAMBDA(const size_t i) { buf_m(i) = dview_m(hash(i)); });
59 Kokkos::fence();
60 }
61
62 template <typename T, class... Properties>
64 auto size = dview_m.extent(0);
65 size_type required = *(this->localNum_mp) + nrecvs;
66 if (size < required) {
67 int overalloc = Comm->getDefaultOverallocation();
68 this->resize(required * overalloc);
69 }
70
71 size_type count = *(this->localNum_mp);
72 using policy_type = Kokkos::RangePolicy<execution_space>;
73 Kokkos::parallel_for(
74 "ParticleAttrib::unpack()", policy_type(0, nrecvs),
75 KOKKOS_CLASS_LAMBDA(const size_t i) { dview_m(count + i) = buf_m(i); });
76 Kokkos::fence();
77 }
78
79 template <typename T, class... Properties>
80 // KOKKOS_INLINE_FUNCTION
82 using policy_type = Kokkos::RangePolicy<execution_space>;
83 Kokkos::parallel_for(
84 "ParticleAttrib::operator=()", policy_type(0, *(this->localNum_mp)),
85 KOKKOS_CLASS_LAMBDA(const size_t i) { dview_m(i) = x; });
86 return *this;
87 }
88
89 template <typename T, class... Properties>
90 template <typename E, size_t N>
91 // KOKKOS_INLINE_FUNCTION
93 detail::Expression<E, N> const& expr) {
94 using capture_type = detail::CapturedExpression<E, N>;
95 capture_type expr_ = reinterpret_cast<const capture_type&>(expr);
96
97 using policy_type = Kokkos::RangePolicy<execution_space>;
98 Kokkos::parallel_for(
99 "ParticleAttrib::operator=()", policy_type(0, *(this->localNum_mp)),
100 KOKKOS_CLASS_LAMBDA(const size_t i) { dview_m(i) = expr_(i); });
101 return *this;
102 }
103
104 template <typename T, class... Properties>
105 template <typename Field, class PT, typename policy_type>
107 Field& f, const ParticleAttrib<Vector<PT, Field::dim>, Properties...>& pp,
108 policy_type iteration_policy, hash_type hash_array) const {
109 constexpr unsigned Dim = Field::dim;
110 using PositionType = typename Field::Mesh_t::value_type;
112 static IpplTimings::TimerRef scatterTimer = IpplTimings::getTimer("scatter");
113 IpplTimings::startTimer(scatterTimer);
114 using view_type = typename Field::view_type;
115 view_type view = f.getView();
116
117 using mesh_type = typename Field::Mesh_t;
118 const mesh_type& mesh = f.get_mesh();
119
120 using vector_type = typename mesh_type::vector_type;
121 using value_type = typename ParticleAttrib<T, Properties...>::value_type;
122
123 const vector_type& dx = mesh.getMeshSpacing();
124 const vector_type& origin = mesh.getOrigin();
125 const vector_type invdx = 1.0 / dx;
126
127 const FieldLayout<Dim>& layout = f.getLayout();
128 const NDIndex<Dim>& lDom = layout.getLocalNDIndex();
129 const int nghost = f.getNghost();
130
131 //using policy_type = Kokkos::RangePolicy<execution_space>;
132 const bool useHashView = hash_array.extent(0) > 0;
133 if (useHashView && (iteration_policy.end() > hash_array.extent(0))) {
134 Inform m("scatter");
135 m << "Hash array was passed to scatter, but size does not match iteration policy." << endl;
136 ippl::Comm->abort();
137 }
138 Kokkos::parallel_for(
139 "ParticleAttrib::scatter", iteration_policy,
140 KOKKOS_CLASS_LAMBDA(const size_t idx) {
141 // map index to possible hash_map
142 size_t mapped_idx = useHashView ? hash_array(idx) : idx;
143
144 // find nearest grid point
145 vector_type l = (pp(mapped_idx) - origin) * invdx + 0.5;
146 Vector<int, Field::dim> index = l;
147 Vector<PositionType, Field::dim> whi = l - index;
148 Vector<PositionType, Field::dim> wlo = 1.0 - whi;
149
150 Vector<size_t, Field::dim> args = index - lDom.first() + nghost;
151
152 // scatter
153 const value_type& val = dview_m(mapped_idx);
154 detail::scatterToField(std::make_index_sequence<1 << Field::dim>{}, view, wlo, whi,
155 args, val);
156 });
157 IpplTimings::stopTimer(scatterTimer);
158
159 static IpplTimings::TimerRef accumulateHaloTimer = IpplTimings::getTimer("accumulateHalo");
160 IpplTimings::startTimer(accumulateHaloTimer);
161 f.accumulateHalo();
162 IpplTimings::stopTimer(accumulateHaloTimer);
163 }
164
165 template <typename T, class... Properties>
166 template <typename Field, typename P2>
168 Field& f, const ParticleAttrib<Vector<P2, Field::dim>, Properties...>& pp,
169 const bool addToAttribute) {
170 constexpr unsigned Dim = Field::dim;
171 using PositionType = typename Field::Mesh_t::value_type;
172
173 static IpplTimings::TimerRef fillHaloTimer = IpplTimings::getTimer("fillHalo");
175 f.fillHalo();
176 IpplTimings::stopTimer(fillHaloTimer);
177
178 static IpplTimings::TimerRef gatherTimer = IpplTimings::getTimer("gather");
179 IpplTimings::startTimer(gatherTimer);
180 const typename Field::view_type view = f.getView();
181
182 using mesh_type = typename Field::Mesh_t;
183 const mesh_type& mesh = f.get_mesh();
184
185 using vector_type = typename mesh_type::vector_type;
186
187 const vector_type& dx = mesh.getMeshSpacing();
188 const vector_type& origin = mesh.getOrigin();
189 const vector_type invdx = 1.0 / dx;
190
191 const FieldLayout<Dim>& layout = f.getLayout();
192 const NDIndex<Dim>& lDom = layout.getLocalNDIndex();
193 const int nghost = f.getNghost();
194
195 using policy_type = Kokkos::RangePolicy<execution_space>;
196 Kokkos::parallel_for(
197 "ParticleAttrib::gather", policy_type(0, *(this->localNum_mp)),
198 KOKKOS_CLASS_LAMBDA(const size_t idx) {
199 // find nearest grid point
200 vector_type l = (pp(idx) - origin) * invdx + 0.5;
201 Vector<int, Field::dim> index = l;
203 Vector<PositionType, Field::dim> wlo = 1.0 - whi;
204
205 Vector<size_t, Field::dim> args = index - lDom.first() + nghost;
206
207 // gather
208 value_type gathered = detail::gatherFromField(std::make_index_sequence<1 << Field::dim>{},
209 view, wlo, whi, args);
210 if (addToAttribute) {
211 dview_m(idx) += gathered;
212 } else {
213 dview_m(idx) = gathered;
214 }
215 });
216 IpplTimings::stopTimer(gatherTimer);
217 }
218
219 template <typename T, class... Properties>
221 const hash_type& permutation) {
222
223 const auto view = this->getView();
224 const auto size = this->getParticleCount();
225
226 view_type temp("copy", size);
227
228 using policy_type = Kokkos::RangePolicy<execution_space>;
229 Kokkos::parallel_for(
230 "Copy to temp", policy_type(0, size),
231 KOKKOS_LAMBDA(const size_type& i) { temp(permutation(i)) = view(i); });
232
233 Kokkos::fence();
234
235 Kokkos::deep_copy(Kokkos::subview(view, Kokkos::make_pair<size_type, size_type>(0, size)), temp);
236 }
237
238 template<typename T, class... Properties>
240 const hash_type &indices) {
241 auto copySize = indices.size();
242 create(copySize);
243
244 auto view = this->getView();
245 const auto size = this->getParticleCount();
246
247 using policy_type = Kokkos::RangePolicy<execution_space>;
248 Kokkos::parallel_for(
249 "Copy to temp", policy_type(0, copySize),
250 KOKKOS_LAMBDA(const size_type &i) {
251 view(size + i) = view(i);
252 });
253
254 Kokkos::fence();
255 }
256
257 /*
258 * Non-class function
259 *
260 */
261
278 template <typename Attrib1, typename Field, typename Attrib2,
279 typename policy_type = Kokkos::RangePolicy<typename Field::execution_space>>
280 inline void scatter(const Attrib1& attrib, Field& f, const Attrib2& pp) {
281 attrib.scatter(f, pp, policy_type(0, attrib.getParticleCount()));
282 }
283
302 template <typename Attrib1, typename Field, typename Attrib2,
303 typename policy_type = Kokkos::RangePolicy<typename Field::execution_space>>
304 inline void scatter(const Attrib1& attrib, Field& f, const Attrib2& pp,
305 policy_type iteration_policy, typename Attrib1::hash_type hash_array = {}) {
306 attrib.scatter(f, pp, iteration_policy, hash_array);
307 }
308
326 template <typename Attrib1, typename Field, typename Attrib2>
327 inline void gather(Attrib1& attrib, Field& f, const Attrib2& pp,
328 const bool addToAttribute = false) {
329 attrib.gather(f, pp, addToAttribute);
330 }
331
332#define DefineParticleReduction(fun, name, op, MPI_Op) \
333 template <typename T, class... Properties> \
334 T ParticleAttrib<T, Properties...>::name() { \
335 T temp = 0.0; \
336 using policy_type = Kokkos::RangePolicy<execution_space>; \
337 Kokkos::parallel_reduce( \
338 "fun", policy_type(0, *(this->localNum_mp)), \
339 KOKKOS_CLASS_LAMBDA(const size_t i, T& valL) { \
340 T myVal = dview_m(i); \
341 op; \
342 }, \
343 Kokkos::fun<T>(temp)); \
344 T globaltemp = 0.0; \
345 Comm->allreduce(temp, globaltemp, 1, MPI_Op<T>()); \
346 return globaltemp; \
347 }
348
349 DefineParticleReduction(Sum, sum, valL += myVal, std::plus)
350 DefineParticleReduction(Max, max, if (myVal > valL) valL = myVal, std::greater)
351 DefineParticleReduction(Min, min, if (myVal < valL) valL = myVal, std::less)
352 DefineParticleReduction(Prod, prod, valL *= myVal, std::multiplies)
353} // namespace ippl
constexpr unsigned Dim
ippl::detail::size_type size_type
Definition datatypes.h:23
ippl::Field< T, Dim, Mesh_t< Dim >, Centering_t< Dim >, ViewArgs... > Field
Definition datatypes.h:29
#define DefineParticleReduction(fun, name, op, MPI_Op)
Inform & endl(Inform &inf)
Definition Inform.cpp:42
STL namespace.
Definition Archive.h:20
std::greater prod
void gather(Attrib1 &attrib, Field &f, const Attrib2 &pp, const bool addToAttribute=false)
Non-class interface for gathering field data into a particle attribute.
DefineParticleReduction(Sum, sum, valL+=myVal, std::plus) DefineParticleReduction(Max
KOKKOS_INLINE_FUNCTION Vector< T, Dim > min(const Vector< T, Dim > &a, const Vector< T, Dim > &b)
Definition Vector.hpp:217
std::unique_ptr< mpi::Communicator > Comm
Definition Ippl.h:22
void scatter(const Attrib1 &attrib, Field &f, const Attrib2 &pp)
Non-class interface for scattering particle attribute data onto a field.
if(myVal > valL) valL
KOKKOS_INLINE_FUNCTION constexpr View::value_type gatherFromField(const std::index_sequence< GatherPoint... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args)
Definition CIC.hpp:68
KOKKOS_INLINE_FUNCTION constexpr void scatterToField(const std::index_sequence< ScatterPoint... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args, T val=1)
Definition CIC.hpp:47
Layout_t & getLayout() const
Definition BareField.h:134
void accumulateHalo()
int getNghost() const
Definition BareField.h:125
view_type & getView()
Definition BareField.h:169
KOKKOS_INLINE_FUNCTION Mesh_t & get_mesh() const
Definition Field.h:64
const NDIndex_t & getLocalNDIndex() const
KOKKOS_INLINE_FUNCTION Vector< int, Dim > first() const
Definition NDIndex.hpp:170
KOKKOS_INLINE_FUNCTION vector_type getOrigin() const
Definition Mesh.hpp:10
virtual KOKKOS_INLINE_FUNCTION const vector_type & getMeshSpacing() const =0
void destroy(const hash_type &deleteIndex, const hash_type &keepIndex, size_type invalidCount) override
void create(size_type) override
void unpack(size_type) override
void internalCopy(const hash_type &indices) override
detail::size_type size_type
void scatter(Field &f, const ParticleAttrib< Vector< P2, Field::dim >, Properties... > &pp, policy_type iteration_policy, hash_type hash_array={}) const
typename detail::ViewType< ippl::Vector< double, 2 >, 1, Properties... >::view_type view_type
void gather(Field &f, const ParticleAttrib< Vector< P2, Field::dim >, Properties... > &pp, const bool addToAttribute=false)
view_type & getView()
typename Base::hash_type hash_type
void realloc(size_type n)
void pack(const hash_type &) override
void applyPermutation(const hash_type &permutation) override
ParticleAttrib< T, Properties... > & operator=(T x)
size_type size() const override
Timing::TimerRef TimerRef
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)