IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ParticleBase.hpp
Go to the documentation of this file.
1//
2// Class ParticleBase
3// Base class for all user-defined particle classes.
4//
5// ParticleBase is a container and manager for a set of particles.
6// The user must define a class derived from ParticleBase which describes
7// what specific data attributes the particle has (e.g., mass or charge).
8// Each attribute is an instance of a ParticleAttribute<T> class; ParticleBase
9// keeps a list of pointers to these attributes, and performs particle creation
10// and destruction.
11//
12// ParticleBase is templated on the ParticleLayout mechanism for the particles.
13// This template parameter should be a class derived from ParticleLayout.
14// ParticleLayout-derived classes maintain the info on which particles are
15// located on which processor, and performs the specific communication
16// required between processors for the particles. The ParticleLayout is
17// templated on the type and dimension of the atom position attribute, and
18// ParticleBase uses the same types for these items as the given
19// ParticleLayout.
20//
21// ParticleBase and all derived classes have the following common
22// characteristics:
23// - The spatial positions of the N particles are stored in the
24// particle_position_type variable R
25// - The global index of the N particles are stored in the
26// particle_index_type variable ID
27// - A pointer to an allocated layout class. When you construct a
28// ParticleBase, you must provide a layout instance, and ParticleBase
29// will delete this instance when it (the ParticleBase) is deleted.
30//
31// To use this class, the user defines a derived class with the same
32// structure as in this example:
33//
34// class UserParticles :
35// public ParticleBase< ParticleSpatialLayout<double,3> > {
36// public:
37// // attributes for this class
38// ParticleAttribute<double> rad; // radius
39// particle_position_type vel; // velocity, same storage type as R
40//
41// // constructor: add attributes to base class
42// UserParticles(ParticleSpatialLayout<double,2>* L) : ParticleBase(L) {
43// addAttribute(rad);
44// addAttribute(vel);
45// }
46// };
47//
48// This example defines a user class with 3D position and two extra
49// attributes: a radius rad (double), and a velocity vel (a 3D Vector).
50//
51
52namespace ippl {
53
54 template <class PLayout, typename... IP>
56 : layout_m(nullptr)
57 , localNum_m(0)
58 , totalNum_m(0)
59 , nextID_m(Comm->rank())
60 , numNodes_m(Comm->size()) {
61 if constexpr (EnableIDs) {
63 }
65 ID.set_name("ID");
66 R.set_name("position");
67 }
68
69 template <class PLayout, typename... IP>
71 : ParticleBase() {
72 initialize(layout);
73 }
74
75 template <class PLayout, typename... IP>
76 template <typename MemorySpace>
81
82 template <class PLayout, typename... IP>
84 // PAssert(layout_m == nullptr);
85
86 // save the layout, and perform setup tasks
87 layout_m = &layout;
88 }
89
90 template <class PLayout, typename... IP>
92 PAssert(layout_m != nullptr);
93
94 if (nLocal > 0) {
95 forAllAttributes([&]<typename Attribute>(Attribute& attribute) {
96 attribute->create(nLocal);
97 });
98
99 if constexpr (EnableIDs) {
100 // set the unique ID value for these new particles
101 using policy_type =
102 Kokkos::RangePolicy<size_type, typename particle_index_type::execution_space>;
103 auto pIDs = ID.getView();
104 auto nextID = this->nextID_m;
105 auto numNodes = this->numNodes_m;
106 Kokkos::parallel_for(
107 "ParticleBase<...>::create(size_t)", policy_type(localNum_m, nLocal),
108 KOKKOS_LAMBDA(const std::int64_t i) { pIDs(i) = nextID + numNodes * i; });
109 // nextID_m += numNodes_m * (nLocal - localNum_m);
110 nextID_m += numNodes_m * nLocal;
111 }
112
113 // remember that we're creating these new particles
114 localNum_m += nLocal;
115 }
116
117 Comm->allreduce(localNum_m, totalNum_m, 1, std::plus<size_type>());
118 }
119
120 template <class PLayout, typename... IP>
122 PAssert(layout_m != nullptr);
123
124 size_type n = (id > -1) ? 1 : 0;
125
126 // temporary change
127 index_type tmpNextID = nextID_m;
128 nextID_m = id;
129 numNodes_m = 0;
130
131 create(n);
132
133 nextID_m = tmpNextID;
134 numNodes_m = Comm->size();
135 }
136
137 template <class PLayout, typename... IP>
139 PAssert(layout_m != nullptr);
140
141 // Compute the number of particles local to each processor
142 size_type nLocal = nTotal / numNodes_m;
143
144 const size_t rank = Comm->rank();
145
146 size_type rest = nTotal - nLocal * rank;
147 if (rank < rest) {
148 ++nLocal;
149 }
151 create(nLocal);
152 }
153
154 template <class PLayout, typename... IP>
155 template <typename... Properties>
156 void ParticleBase<PLayout, IP...>::destroy(const Kokkos::View<bool*, Properties...>& invalid,
157 const size_type destroyNum) {
158 this->internalDestroy(invalid, destroyNum);
159
160 Comm->allreduce(localNum_m, totalNum_m, 1, std::plus<size_type>());
161 }
162
163 template <class PLayout, typename... IP>
164 template <typename... Properties>
166 const Kokkos::View<bool*, Properties...>& invalid, const size_type destroyNum) {
167 PAssert(destroyNum <= localNum_m);
168
169 // If there aren't any particles to delete, do nothing
170 if (destroyNum == 0) {
171 return;
172 }
173
174 // If we're deleting all the particles, there's no point in doing
175 // anything because the valid region will be empty; we only need to
176 // update the particle count
177 if (destroyNum == localNum_m) {
178 localNum_m = 0;
179 return;
180 }
181
182 using view_type = Kokkos::View<bool*, Properties...>;
183 using memory_space = typename view_type::memory_space;
184 using execution_space = typename view_type::execution_space;
185 using policy_type = Kokkos::RangePolicy<execution_space>;
186 auto& locDeleteIndex = deleteIndex_m.get<memory_space>();
187 auto& locKeepIndex = keepIndex_m.get<memory_space>();
188
189 // Resize buffers, if necessary
190 detail::runForAllSpaces([&]<typename MemorySpace>() {
191 if (attributes_m.template get<memory_space>().size() > 0) {
192 int overalloc = Comm->getDefaultOverallocation();
193 auto& del = deleteIndex_m.get<memory_space>();
194 auto& keep = keepIndex_m.get<memory_space>();
195 if (del.size() < destroyNum) {
196 Kokkos::realloc(del, destroyNum * overalloc);
197 Kokkos::realloc(keep, destroyNum * overalloc);
198 }
199 }
200 });
201
202 // Reset index buffer
203 Kokkos::deep_copy(locDeleteIndex, -1);
204
205 // Find the indices of the invalid particles in the valid region
206 Kokkos::parallel_scan(
207 "Scan in ParticleBase::destroy()", policy_type(0, localNum_m - destroyNum),
208 KOKKOS_LAMBDA(const size_t i, int& idx, const bool final) {
209 if (final && invalid(i)) {
210 locDeleteIndex(idx) = i;
211 }
212 if (invalid(i)) {
213 idx += 1;
214 }
215 });
216 Kokkos::fence();
217
218 // Determine the total number of invalid particles in the valid region
219 size_type maxDeleteIndex = 0;
220 Kokkos::parallel_reduce(
221 "Reduce in ParticleBase::destroy()", policy_type(0, destroyNum),
222 KOKKOS_LAMBDA(const size_t i, size_t& maxIdx) {
223 if (locDeleteIndex(i) >= 0 && i > maxIdx) {
224 maxIdx = i;
225 }
226 },
227 Kokkos::Max<size_type>(maxDeleteIndex));
228
229 // Find the indices of the valid particles in the invalid region
230 Kokkos::parallel_scan(
231 "Second scan in ParticleBase::destroy()",
232 Kokkos::RangePolicy<size_type, execution_space>(localNum_m - destroyNum, localNum_m),
233 KOKKOS_LAMBDA(const size_t i, int& idx, const bool final) {
234 if (final && !invalid(i)) {
235 locKeepIndex(idx) = i;
236 }
237 if (!invalid(i)) {
238 idx += 1;
239 }
240 });
241
242 Kokkos::fence();
243
244 localNum_m -= destroyNum;
245
246 // We need to delete particles in all memory spaces. If there are any attributes not stored
247 // in the memory space we've already been using, we need to copy the index views to the
248 // other spaces.
249 auto filter = [&]<typename MemorySpace>() {
250 return attributes_m.template get<MemorySpace>().size() > 0;
251 };
252 deleteIndex_m.copyToOtherSpaces<memory_space>(filter);
253 keepIndex_m.copyToOtherSpaces<memory_space>(filter);
254
255 // Partition the attributes into valid and invalid regions
256 // NOTE: The vector elements are pointers, but we want to extract
257 // the memory space from the class type, so we explicitly
258 // make the lambda argument a pointer to the template parameter
259 forAllAttributes([&]<typename Attribute>(Attribute*& attribute) {
260 using att_memory_space = typename Attribute::memory_space;
261 auto& del = deleteIndex_m.get<att_memory_space>();
262 auto& keep = keepIndex_m.get<att_memory_space>();
263 attribute->destroy(del, keep, maxDeleteIndex + 1);
264 });
265 }
266
267 template <class PLayout, typename... IP>
268 template <typename HashType>
270 std::vector<MPI_Request>& requests,
271 const HashType& hash) {
272 size_type nSends = hash.size();
273 requests.resize(requests.size() + 1);
274
275 auto hashes = hash_container_type(hash, [&]<typename MemorySpace>() {
276 return attributes_m.template get<MemorySpace>().size() > 0;
277 });
278 pack(hashes);
279 detail::runForAllSpaces([&]<typename MemorySpace>() {
280 size_type bufSize = packedSize<MemorySpace>(nSends);
281 if (bufSize == 0) {
282 return;
283 }
284
285 auto buf = Comm->getBuffer<MemorySpace>(bufSize);
286
287 Comm->isend(rank, tag++, *this, *buf, requests.back(), nSends);
288 buf->resetWritePos();
289 });
291
292 template <class PLayout, typename... IP>
294 detail::runForAllSpaces([&]<typename MemorySpace>() {
295 size_type bufSize = packedSize<MemorySpace>(nRecvs);
296 if (bufSize == 0) {
297 return;
298 }
299
300 auto buf = Comm->getBuffer<MemorySpace>(bufSize);
301
302 Comm->recv(rank, tag++, *this, *buf, bufSize, nRecvs);
303 buf->resetReadPos();
304 });
305 unpack(nRecvs);
306 }
307
308 template <class PLayout, typename... IP>
309 template <typename Archive>
311 using memory_space = typename Archive::buffer_type::memory_space;
312 forAllAttributes<memory_space>([&]<typename Attribute>(Attribute& att) {
313 att->serialize(ar, nsends);
314 });
315 }
316
317 template <class PLayout, typename... IP>
318 template <typename Archive>
320 using memory_space = typename Archive::buffer_type::memory_space;
321 forAllAttributes<memory_space>([&]<typename Attribute>(Attribute& att) {
322 att->deserialize(ar, nrecvs);
323 });
324 }
325
326 template <class PLayout, typename... IP>
327 template <typename MemorySpace>
329 size_type total = 0;
330 forAllAttributes<MemorySpace>([&]<typename Attribute>(const Attribute& att) {
331 total += att->packedSize(count);
332 });
333 return total;
334 }
335
336 template <class PLayout, typename... IP>
338 detail::runForAllSpaces([&]<typename MemorySpace>() {
339 auto& att = attributes_m.template get<MemorySpace>();
340 for (unsigned j = 0; j < att.size(); j++) {
341 att[j]->pack(hash.template get<MemorySpace>());
342 }
343 });
344 }
346 template <class PLayout, typename... IP>
348 detail::runForAllSpaces([&]<typename MemorySpace>() {
349 auto& att = attributes_m.template get<MemorySpace>();
350 for (unsigned j = 0; j < att.size(); j++) {
351 att[j]->unpack(nrecvs);
352 }
353 });
354 localNum_m += nrecvs;
355 }
356} // namespace ippl
typename ippl::detail::ViewType< ippl::Vector< double, Dim >, 1 >::view_type view_type
ippl::detail::size_type size_type
Definition datatypes.h:23
void unpack(const ippl::NDIndex< 3 > intersect, const Kokkos::View< Tf *** > &view, ippl::detail::FieldBufferData< Tb > &fd, int nghost, const ippl::NDIndex< 3 > ldom, bool x=false, bool y=false, bool z=false)
void pack(const ippl::NDIndex< 3 > intersect, Kokkos::View< Tf *** > &view, ippl::detail::FieldBufferData< Tb > &fd, int nghost, const ippl::NDIndex< 3 > ldom, ippl::mpi::Communicator::size_type &nsends)
#define PAssert(c)
Definition PAssert.h:116
Definition Archive.h:20
std::unique_ptr< mpi::Communicator > Comm
Definition Ippl.h:22
KOKKOS_INLINE_FUNCTION auto & get(Tuple< Ts... > &t)
Accessor function to get an element mutable reference at a specific index from a Tuple.
Definition Tuple.h:314
std::size_t size_type
Definition IpplTypes.h:13
void runForAllSpaces(Functor &&f)
Definition TypeUtils.h:428
void create(size_type nLocal)
void addAttribute(detail::ParticleAttribBase< MemorySpace > &pa)
typename detail::ContainerForAllSpaces< detail::hash_type >::type hash_container_type
index_type nextID_m
next unique particle ID
void pack(const hash_container_type &hash)
size_type localNum_m
processor local number of particles
void serialize(Archive &ar, size_type nsends)
attribute_container_type attributes_m
all attributes
void internalDestroy(const Kokkos::View< bool *, Properties... > &invalid, const size_type destroyNum)
detail::size_type size_type
void globalCreate(size_type nTotal)
void forAllAttributes(Functor &&f) const
Layout_t * layout_m
particle layout
void initialize(Layout_t &layout)
void unpack(size_type nrecvs)
void destroy(const Kokkos::View< bool *, Properties... > &invalid, const size_type destroyNum)
size_type packedSize(const size_type count) const
static constexpr bool EnableIDs
void recvFromRank(int rank, int tag, size_type nRecvs)
void createWithID(index_type id)
typename PLayout::index_type index_type
index_type numNodes_m
number of MPI ranks
hash_container_type keepIndex_m
void sendToRank(int rank, int tag, std::vector< MPI_Request > &requests, const HashType &hash)
size_type totalNum_m
total number of particles (across all processes)
hash_container_type deleteIndex_m
buffers for particle partitioning
particle_position_type R
view of particle positions
particle_index_type ID
view of particle IDs
void deserialize(Archive &ar, size_type nrecvs)