Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleContainerI.H
Go to the documentation of this file.
2
3#include <algorithm>
4#include <iterator>
5#include <set>
6#include <stdexcept>
7#include <string>
8#include <type_traits>
9#include <utility>
10#include <vector>
11
12namespace amrex {
13
14template <typename ParticleType, int NArrayReal, int NArrayInt,
15 template<class> class Allocator, class CellAssignor>
16void
18{
19 num_real_comm_comps = 0;
20 int comm_comps_start = 0;
21 if constexpr (!ParticleType::is_soa_particle) {
22 comm_comps_start += AMREX_SPACEDIM + NStructReal;
23 }
24 for (int i = comm_comps_start; i < comm_comps_start + NumRealComps(); ++i) {
25 if (h_redistribute_real_comp[i]) {++num_real_comm_comps;}
26 }
27
28 num_int_comm_comps = 0;
29 comm_comps_start = 2 + NStructInt;
30 for (int i = comm_comps_start; i < comm_comps_start + NumIntComps(); ++i) {
31 if (h_redistribute_int_comp[i]) {++num_int_comm_comps;}
32 }
33
34 if constexpr (ParticleType::is_soa_particle) {
35 particle_size = sizeof(uint64_t); // idcpu
36 } else {
37 particle_size = sizeof(ParticleType);
38 }
39 superparticle_size = particle_size +
40 num_real_comm_comps*sizeof(ParticleReal) + num_int_comm_comps*sizeof(int);
41}
42
43template <typename ParticleType, int NArrayReal, int NArrayInt,
44 template<class> class Allocator, class CellAssignor>
45void
47{
48 levelDirectoriesCreated = false;
49 usePrePost = false;
50 doUnlink = true;
51
52 // by default communicate all components
53 if constexpr (ParticleType::is_soa_particle)
54 {
55 h_redistribute_real_comp.resize(NArrayReal, true);
56 } else {
57 h_redistribute_real_comp.resize(AMREX_SPACEDIM + NStructReal + NArrayReal, true);
58 }
59 h_redistribute_int_comp.resize(2 + NStructInt + NArrayInt, true);
60
61 SetParticleSize();
62
63 // add default names for SoA Real and Int compile-time arguments
64 m_soa_rdata_names.clear();
65 for (int i=0; i<NArrayReal; ++i)
66 {
67 m_soa_rdata_names.push_back(getDefaultCompNameReal<ParticleType>(i));
68 }
69 m_soa_idata_names.clear();
70 for (int i=0; i<NArrayInt; ++i)
71 {
72 m_soa_idata_names.push_back(getDefaultCompNameInt<ParticleType>(i));
73 }
74
75 static bool initialized = false;
76 if ( ! initialized)
77 {
78 static_assert(sizeof(ParticleType)%sizeof(RealType) == 0,
79 "sizeof ParticleType is not a multiple of sizeof RealType");
80
81 ParmParse pp("particles");
82 pp.queryAdd("do_tiling", do_tiling);
83 Vector<int> tilesize(AMREX_SPACEDIM);
84 if (pp.queryarr("tile_size", tilesize, 0, AMREX_SPACEDIM)) {
85 for (int i=0; i<AMREX_SPACEDIM; ++i) { tile_size[i] = tilesize[i]; }
86 }
87
88 static_assert(std::is_standard_layout_v<ParticleType>,
89 "Particle type must be standard layout");
90 // && std::is_trivial<ParticleType>::value,
91 // "Particle type must be standard layout and trivial.");
92
93 pp.query("use_prepost", usePrePost);
94 pp.query("do_unlink", doUnlink);
95 pp.queryAdd("do_mem_efficient_sort", memEfficientSort);
96 pp.queryAdd("use_comms_arena", use_comms_arena);
97
98 initialized = true;
99 }
100}
101
102template<
103 typename ParticleType,
104 int NArrayReal,
105 int NArrayInt,
106 template<class> class Allocator,
107 class CellAssignor
108>
109void
111 std::vector<std::string> const & rdata_name, std::vector<std::string> const & idata_name
112)
113{
114 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(std::ssize(rdata_name) == NArrayReal, "rdata_name must be equal to NArrayReal");
115 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(std::ssize(idata_name) == NArrayInt, "idata_name must be equal to NArrayInt");
116
117 // ensure names for components are unique
118 std::set<std::string> const unique_r_names(rdata_name.begin(), rdata_name.end());
119 std::set<std::string> const unique_i_names(idata_name.begin(), idata_name.end());
120 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(rdata_name.size() == unique_r_names.size(), "SetSoACompileTimeNames: Provided names in rdata_name are not unique!");
121 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(idata_name.size() == unique_i_names.size(), "SetSoACompileTimeNames: Provided names in idata_name are not unique!");
122
123 for (int i=0; i<NArrayReal; ++i)
124 {
125 m_soa_rdata_names.at(i) = rdata_name.at(i);
126 }
127 for (int i=0; i<NArrayInt; ++i)
128 {
129 m_soa_idata_names.at(i) = idata_name.at(i);
130 }
131}
132
133template<
134 typename ParticleType,
135 int NArrayReal,
136 int NArrayInt,
137 template<class> class Allocator,
138 class CellAssignor
139>
140bool
142{
143 return std::find(m_soa_rdata_names.begin(), m_soa_rdata_names.end(), name) != std::end(m_soa_rdata_names);
144}
145
146template <typename ParticleType, int NArrayReal, int NArrayInt,
147 template<class> class Allocator, class CellAssignor>
148bool
150{
151 return std::find(m_soa_idata_names.begin(), m_soa_idata_names.end(), name) != std::end(m_soa_idata_names);
152}
153
154template<
155 typename ParticleType,
156 int NArrayReal,
157 int NArrayInt,
158 template<class> class Allocator,
159 class CellAssignor
160>
161int
163{
164 const auto it = std::find(m_soa_rdata_names.begin(), m_soa_rdata_names.end(), name);
165
166 if (it == m_soa_rdata_names.end())
167 {
168 throw std::runtime_error("GetRealCompIndex: Component " + name + " does not exist!");
169 }
170 else
171 {
172 return std::distance(m_soa_rdata_names.begin(), it);
173 }
174}
175
176template<
177 typename ParticleType,
178 int NArrayReal,
179 int NArrayInt,
180 template<class> class Allocator,
181 class CellAssignor
182>
183int
185{
186 const auto it = std::find(m_soa_idata_names.begin(), m_soa_idata_names.end(), name);
187
188 if (it == m_soa_idata_names.end())
189 {
190 throw std::runtime_error("GetIntCompIndex: Component " + name + " does not exist!");
191 }
192 else
193 {
194 return std::distance(m_soa_idata_names.begin(), it);
195 }
196}
197
198template <typename ParticleType, int NArrayReal, int NArrayInt,
199 template<class> class Allocator, class CellAssignor>
200template <typename P, typename Assignor>
203{
204 const Geometry& geom = Geom(lev);
205 const auto& domain = geom.Domain();
206 const auto& plo = geom.ProbLoArray();
207 const auto& dxi = geom.InvCellSizeArray();
208
209 return Assignor{}(p, plo, dxi, domain);
210}
211
212template <typename ParticleType, int NArrayReal, int NArrayInt,
213 template<class> class Allocator, class CellAssignor>
214template <typename P>
215bool
217::Where (const P& p,
218 ParticleLocData& pld,
219 int lev_min,
220 int lev_max,
221 int nGrow,
222 int local_grid) const
223{
224
225 AMREX_ASSERT(m_gdb != nullptr);
226
227 if (lev_max == -1) {
228 lev_max = finestLevel();
229 }
230
231 AMREX_ASSERT(lev_max <= finestLevel());
232
233 AMREX_ASSERT(nGrow == 0 || (nGrow >= 0 && lev_min == lev_max));
234
235 std::vector< std::pair<int, Box> > isects;
236
237 for (int lev = lev_max; lev >= lev_min; lev--) {
238 const IntVect& iv = Index(p, lev);
239 if (lev == pld.m_lev) {
240 // The fact that we are here means this particle does not belong to any finer grids.
241 if (pld.m_grid >= 0) {
242 if (pld.m_grown_gridbox.contains(iv)) {
243 pld.m_cell = iv;
244 if (!pld.m_tilebox.contains(iv)) {
245 pld.m_tile = getTileIndex(iv, pld.m_gridbox, do_tiling, tile_size, pld.m_tilebox);
246 }
247 return true;
248 }
249 }
250 }
251
252 int grid;
253 const BoxArray& ba = ParticleBoxArray(lev);
255
256 if (local_grid < 0) {
257 bool findfirst = (nGrow == 0) ? true : false;
258 ba.intersections(Box(iv, iv), isects, findfirst, nGrow);
259 grid = isects.empty() ? -1 : isects[0].first;
260 if (nGrow > 0 && std::ssize(isects) > 1) {
261 for (auto & isect : isects) {
262 Box bx = ba.getCellCenteredBox(isect.first);
263 for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
264 Box gbx = bx;
266 gr[dir] = nGrow;
267 gbx.grow(gr);
268 if (gbx.contains(iv)) {
269 grid = isect.first;
270 }
271 }
272 }
273 }
274 } else {
275 grid = (*redistribute_mask_ptr)[local_grid](iv, 0);
276 }
277
278 if (grid >= 0) {
279 const Box& bx = ba.getCellCenteredBox(grid);
280 pld.m_lev = lev;
281 pld.m_grid = grid;
282 pld.m_tile = getTileIndex(iv, bx, do_tiling, tile_size, pld.m_tilebox);
283 pld.m_cell = iv;
284 pld.m_gridbox = bx;
285 pld.m_grown_gridbox = amrex::grow(bx, nGrow);
286 return true;
287 }
288 }
289
290 return false;
291}
292
293template <typename ParticleType, int NArrayReal, int NArrayInt,
294 template<class> class Allocator, class CellAssignor>
295template <typename P>
296bool
299 ParticleLocData& pld,
300 int lev_min,
301 int lev_max,
302 int local_grid) const
303{
304
305 AMREX_ASSERT(m_gdb != nullptr);
306
307 if (!Geom(0).isAnyPeriodic()) { return false; }
308
309 if (lev_max == -1) {
310 lev_max = finestLevel();
311 }
312
313 AMREX_ASSERT(lev_max <= finestLevel());
314
315 // Create a copy "dummy" particle to check for periodic outs.
316 Particle<0, 0> p_prime;
317 AMREX_D_TERM(p_prime.pos(0) = p.pos(0);,
318 p_prime.pos(1) = p.pos(1);,
319 p_prime.pos(2) = p.pos(2));
320 if (PeriodicShift(p_prime)) {
321 std::vector< std::pair<int,Box> > isects;
322 for (int lev = lev_max; lev >= lev_min; lev--) {
323
324 int grid;
325 IntVect iv;
326 const BoxArray& ba = ParticleBoxArray(lev);
328
329 if (local_grid < 0) {
330 iv = Index<amrex::Particle<0, 0>, DefaultAssignor>(p_prime, lev);
331 ba.intersections(Box(iv, iv), isects, true, 0);
332 grid = isects.empty() ? -1 : isects[0].first;
333 } else {
334 iv = Index<amrex::Particle<0, 0>, DefaultAssignor>(p_prime, lev);
335 if (ba[local_grid].contains(iv))
336 {
337 grid = local_grid;
338 }
339 else
340 {
341 ba.intersections(Box(iv, iv), isects, true, 0);
342 grid = isects.empty() ? -1 : isects[0].first;
343 if(grid == -1)
344 {
345 grid = (*redistribute_mask_ptr)[local_grid](Index(p, lev), 0);
346 }
347 }
348 }
349
350 if (grid >= 0) {
351 AMREX_D_TERM(p.pos(0) = p_prime.pos(0);,
352 p.pos(1) = p_prime.pos(1);,
353 p.pos(2) = p_prime.pos(2););
354
355 const Box& bx = ba.getCellCenteredBox(grid);
356
357 pld.m_lev = lev;
358 pld.m_grid = grid;
359 pld.m_tile = getTileIndex(iv, bx, do_tiling, tile_size, pld.m_tilebox);
360 pld.m_cell = iv;
361 pld.m_gridbox = bx;
362 pld.m_grown_gridbox = bx;
363 return true;
364 }
365 }
366 }
367
368 return false;
369}
370
371
372template <typename ParticleType, int NArrayReal, int NArrayInt,
373 template<class> class Allocator, class CellAssignor>
374template <typename P>
375bool
377::PeriodicShift (P& p) const
378{
379 const auto& geom = Geom(0);
380 const auto plo = geom.ProbLoArray();
381 const auto phi = geom.ProbHiArray();
382 const auto rlo = geom.ProbLoArrayInParticleReal();
383 const auto rhi = geom.ProbHiArrayInParticleReal();
384 const auto is_per = geom.isPeriodicArray();
386 return enforcePeriodic(p, plo, phi, rlo, rhi, is_per);
387}
388
389template <typename ParticleType, int NArrayReal, int NArrayInt,
390 template<class> class Allocator, class CellAssignor>
394 bool /*update*/,
395 bool verbose,
396 ParticleLocData pld) const
397{
398 AMREX_ASSERT(m_gdb != nullptr);
399
400 bool ok = Where(p, pld);
401
402 if (!ok && Geom(0).isAnyPeriodic())
403 {
404 // Attempt to shift the particle back into the domain if it
405 // crossed a periodic boundary.
406 PeriodicShift(p);
407 ok = Where(p, pld);
408 }
409
410 if (!ok) {
411 // invalidate the particle.
412 if (verbose) {
413 amrex::AllPrint()<< "Invalidating out-of-domain particle: " << p << '\n';
414 }
415
416 AMREX_ASSERT(p.id().is_valid());
417
418 p.id().make_invalid();
419 }
420
421 return pld;
422}
423
424template <typename ParticleType, int NArrayReal, int NArrayInt,
425 template<class> class Allocator, class CellAssignor>
426void
432
433template <typename ParticleType, int NArrayReal, int NArrayInt,
434 template<class> class Allocator, class CellAssignor>
435void
437{
439 int nlevs = std::max(0, finestLevel()+1);
440 m_particles.resize(nlevs);
441}
442
443template <typename ParticleType, int NArrayReal, int NArrayInt,
444 template<class> class Allocator, class CellAssignor>
445template <typename P>
446void
448 int lev_min, int lev_max, int nGrow, int local_grid) const
449{
450 bool success;
451 if (Geom(0).outsideRoundoffDomain(AMREX_D_DECL(p.pos(0), p.pos(1), p.pos(2))))
452 {
453 // Note that EnforcePeriodicWhere may shift the particle if it is successful.
454 success = EnforcePeriodicWhere(p, pld, lev_min, lev_max, local_grid);
455 if (!success && lev_min == 0)
456 {
457 // The particle has left the domain; invalidate it.
458 p.id().make_invalid();
459 success = true;
460 }
461 }
462 else
463 {
464 success = Where(p, pld, lev_min, lev_max, 0, local_grid);
465 }
466
467 if (!success)
468 {
469 success = (nGrow > 0) && Where(p, pld, lev_min, lev_min, nGrow);
470 pld.m_grown_gridbox = pld.m_gridbox; // reset grown box for subsequent calls.
471 }
473 if (!success)
475 amrex::Abort("ParticleContainer::locateParticle(): invalid particle.");
476 }
477}
478
479template <typename ParticleType, int NArrayReal, int NArrayInt,
480 template<class> class Allocator, class CellAssignor>
481Long
483{
484 Long nparticles = 0;
485 for (int lev = 0; lev <= finestLevel(); lev++) {
486 nparticles += NumberOfParticlesAtLevel(lev,only_valid,true);
487 }
488 if (!only_local) {
490 }
491 return nparticles;
492}
493
494template <typename ParticleType, int NArrayReal, int NArrayInt,
495 template<class> class Allocator, class CellAssignor>
498{
499 AMREX_ASSERT(lev >= 0 && lev < std::ssize(m_particles));
500
501 LayoutData<Long> np_per_grid_local(ParticleBoxArray(lev),
502 ParticleDistributionMap(lev));
503
504 for (ParConstIterType pti(*this, lev); pti.isValid(); ++pti)
505 {
506 int gid = pti.index();
507 if (only_valid)
508 {
509 const auto& ptile = ParticlesAt(lev, pti);
510 const int np = ptile.numParticles();
511 auto const ptd = ptile.getConstParticleTileData();
512
514 ReduceData<int> reduce_data(reduce_op);
515 using ReduceTuple = typename decltype(reduce_data)::Type;
516
517 reduce_op.eval(np, reduce_data,
518 [=] AMREX_GPU_DEVICE (int i) -> ReduceTuple
519 {
520 return (ptd.id(i).is_valid()) ? 1 : 0;
521 });
522
523 int np_valid = amrex::get<0>(reduce_data.value(reduce_op));
524 np_per_grid_local[gid] += np_valid;
525 } else
526 {
527 np_per_grid_local[gid] += pti.numParticles();
528 }
530
531 Vector<Long> nparticles(np_per_grid_local.size(), 0);
532 if (only_local)
533 {
534 for (ParConstIterType pti(*this, lev); pti.isValid(); ++pti)
535 {
536 nparticles[pti.index()] = np_per_grid_local[pti.index()];
537 }
538 }
539 else
540 {
541 ParallelDescriptor::GatherLayoutDataToVector(np_per_grid_local, nparticles,
543 ParallelDescriptor::Bcast(nparticles.data(), nparticles.size(),
546
547 return nparticles;
548}
549
550template <typename ParticleType, int NArrayReal, int NArrayInt,
551 template<class> class Allocator, class CellAssignor>
553{
554 Long nparticles = 0;
555
556 if (level < 0 || level >= std::ssize(m_particles)) { return nparticles; }
557
558 if (only_valid) {
559 ReduceOps<ReduceOpSum> reduce_op;
560 ReduceData<unsigned long long> reduce_data(reduce_op);
561 using ReduceTuple = typename decltype(reduce_data)::Type;
562
563 for (const auto& kv : GetParticles(level)) {
564 const auto& ptile = kv.second;
565 auto const ptd = ptile.getConstParticleTileData();
566
567 reduce_op.eval(ptile.numParticles(), reduce_data,
568 [=] AMREX_GPU_DEVICE (int i) -> ReduceTuple
569 {
570 return (ptd.id(i).is_valid()) ? 1 : 0;
571 });
573
574 nparticles = static_cast<Long>(amrex::get<0>(reduce_data.value(reduce_op)));
575 }
576 else {
577 for (const auto& kv : GetParticles(level)) {
578 const auto& ptile = kv.second;
579 nparticles += ptile.numParticles();
580 }
581 }
582
583 if (!only_local) {
585 }
586
587 return nparticles;
589
590template <typename ParticleType, int NArrayReal, int NArrayInt,
591 template<class> class Allocator, class CellAssignor>
592template <std::integral I>
593requires (sizeof(I) >= sizeof(Long))
594void
596{
597 AMREX_ASSERT(lev >= 0 && lev < std::ssize(m_particles));
598
599 AMREX_ASSERT(BoxArray::SameRefs(mem.boxArray(), ParticleBoxArray(lev)) &&
600 mem.DistributionMap() == ParticleDistributionMap(lev));
601
602 [[maybe_unused]] Gpu::NoSyncRegion no_sync{};
603 for (ParConstIterType pti(*this, lev); pti.isValid(); ++pti)
604 {
605 int gid = pti.index();
606 mem[gid] += static_cast<I>(pti.capacity());
608}
609
610//
611// This includes both valid and invalid particles since invalid particles still take up space.
612//
613
614template <typename ParticleType, int NArrayReal, int NArrayInt,
615 template<class> class Allocator, class CellAssignor>
616std::array<Long, 3>
620 Long cnt = 0;
621
622 for (unsigned lev = 0; lev < m_particles.size(); lev++) {
623 const auto& pmap = m_particles[lev];
624 for (const auto& kv : pmap) {
625 const auto& ptile = kv.second;
626 cnt += ptile.numParticles();
627 }
629
630 Long mn = cnt, mx = mn;
631
632 const int IOProc = ParallelContext::IOProcessorNumberSub();
633 const Long sz = sizeof(ParticleType)+NumRealComps()*sizeof(ParticleReal)+NumIntComps()*sizeof(int);
634
635#ifdef AMREX_LAZY
636 Lazy::QueueReduction( [=] () mutable {
637#endif
641
642 amrex::Print() << "ParticleContainer spread across MPI nodes - bytes (num particles): [Min: "
643 << mn*sz
644 << " (" << mn << ")"
645 << ", Max: "
646 << mx*sz
647 << " (" << mx << ")"
648 << ", Total: "
649 << cnt*sz
650 << " (" << cnt << ")]\n";
651#ifdef AMREX_LAZY
652 });
653#endif
654
655 return {mn*sz, mx*sz, cnt*sz};
657
658template <typename ParticleType, int NArrayReal, int NArrayInt,
659 template<class> class Allocator, class CellAssignor>
660std::array<Long, 3>
663{
664 Long cnt = 0;
665
666 for (unsigned lev = 0; lev < m_particles.size(); lev++) {
667 const auto& pmap = m_particles[lev];
668 for (const auto& kv : pmap) {
669 const auto& ptile = kv.second;
670 cnt += ptile.capacity();
671 }
672 }
673
674 Long mn = cnt, mx = mn;
676 const int IOProc = ParallelContext::IOProcessorNumberSub();
677
678#ifdef AMREX_LAZY
679 Lazy::QueueReduction( [=] () mutable {
680#endif
684
685 amrex::Print() << "ParticleContainer spread across MPI nodes - bytes: [Min: "
686 << mn
687 << ", Max: "
688 << mx
689 << ", Total: "
690 << cnt
691 << "]\n";
692#ifdef AMREX_LAZY
693 });
694#endif
695
696 return {mn, mx, cnt};
697}
698
699template <typename ParticleType, int NArrayReal, int NArrayInt,
700 template<class> class Allocator, class CellAssignor>
701void
703{
704 for (unsigned lev = 0; lev < m_particles.size(); lev++) {
705 auto& pmap = m_particles[lev];
706 for (auto& kv : pmap) {
707 auto& ptile = kv.second;
708 ptile.shrink_to_fit();
709 }
710 }
711}
712
719template <typename ParticleType, int NArrayReal, int NArrayInt,
720 template<class> class Allocator, class CellAssignor>
721void
724 BL_PROFILE("ParticleContainer::Increment");
725
726 AMREX_ASSERT(OK());
727 if (m_particles.empty()) { return; }
728 AMREX_ASSERT(lev >= 0 && lev < std::ssize(m_particles));
729 AMREX_ASSERT(numParticlesOutOfRange(*this, 0) == 0);
730
731 const auto& geom = Geom(lev);
732 const auto plo = geom.ProbLoArray();
733 const auto dxi = geom.InvCellSizeArray();
734 const auto domain = geom.Domain();
735 amrex::ParticleToMesh(*this, mf, lev,
736 [=] AMREX_GPU_DEVICE (const typename ParticleTileType::ConstParticleTileDataType& ptd, int ip,
737 amrex::Array4<amrex::Real> const& count)
738 {
739 const auto p = ptd[ip];
740 CellAssignor assignor;
741 IntVect iv = assignor(p, plo, dxi, domain);
742 amrex::Gpu::Atomic::AddNoRet(&count(iv), 1.0_rt);
743 }, false);
744}
745
746template <typename ParticleType, int NArrayReal, int NArrayInt,
747 template<class> class Allocator, class CellAssignor>
748Long
750{
751 BL_PROFILE("ParticleContainer::IncrementWithTotal(lev)");
752 Increment(mf, lev);
753 return TotalNumberOfParticles(true, local);
754}
755
756template <typename ParticleType, int NArrayReal, int NArrayInt,
757 template<class> class Allocator, class CellAssignor>
758void
760{
761 BL_PROFILE("ParticleContainer::RemoveParticlesAtLevel()");
762 if (level >= std::ssize(this->m_particles)) { return; }
763
764 if (!this->m_particles[level].empty())
765 {
766 ParticleLevel().swap(this->m_particles[level]);
767 }
768}
769
770template <typename ParticleType, int NArrayReal, int NArrayInt,
771 template<class> class Allocator, class CellAssignor>
772void
774{
775 BL_PROFILE("ParticleContainer::RemoveParticlesNotAtFinestLevel()");
776 AMREX_ASSERT(this->finestLevel()+1 == std::ssize(this->m_particles));
777
778 Long cnt = 0;
779
780 for (unsigned lev = 0; lev < m_particles.size() - 1; ++lev) {
781 auto& pmap = m_particles[lev];
782 if (!pmap.empty()) {
783 for (auto& kv : pmap) {
784 const auto& pbx = kv.second;
785 cnt += pbx.numParticles();
786 }
787 ParticleLevel().swap(pmap);
788 }
789 }
790
791 //
792 // Print how many particles removed on each processor if any were removed.
793 //
794 if (this->m_verbose > 1 && cnt > 0) {
795 amrex::AllPrint() << "Processor " << ParallelContext::MyProcSub() << " removed " << cnt
796 << " particles not in finest level\n";
797 }
798}
799
801{
802
806
808 const GpuArray<Real, AMREX_SPACEDIM> & dxi, const Box& domain)
809 : m_assign_buffer_grid(assign_buffer_grid), m_plo(plo), m_dxi(dxi), m_domain(domain)
810 {}
811
812 template <typename SrcData>
814 int operator() (const SrcData& src, int src_i) const noexcept
815 {
816 auto iv = getParticleCell(src, src_i, m_plo, m_dxi, m_domain);
817 return (m_assign_buffer_grid(iv).first != -1);
818 }
819};
820
822{
823 template <typename DstData, typename SrcData>
825 void operator() (DstData& dst, const SrcData& src,
826 int src_i, int dst_i) const noexcept
827 {
828 copyParticle(dst, src, src_i, dst_i);
829
831 dst.cpu(dst_i) = 0;
832 }
833};
834
835
836template <typename ParticleType, int NArrayReal, int NArrayInt,
837 template<class> class Allocator, class CellAssignor>
838void
839ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
840::CreateVirtualParticles (int level, AoS& virts) const
841{
842 ParticleTileType ptile;
843 CreateVirtualParticles(level, ptile);
844 ptile.GetArrayOfStructs().swap(virts);
845}
846
847template <typename ParticleType, int NArrayReal, int NArrayInt,
848 template<class> class Allocator, class CellAssignor>
849void
852{
853 BL_PROFILE("ParticleContainer::CreateVirtualParticles()");
854 AMREX_ASSERT(level > 0);
855 AMREX_ASSERT(virts.empty());
856
857 if (level >= std::ssize(m_particles)) {
858 return;
859 }
860
861 std::string const& aggregation_type = AggregationType();
862 int aggregation_buffer = AggregationBuffer();
863
864 if (aggregation_type == "None")
865 {
866 auto virts_offset = virts.numParticles();
867 for(ParConstIterType pti(*this, level); pti.isValid(); ++pti)
868 {
869 const auto& src_tile = ParticlesAt(level, pti);
870
871 auto np = src_tile.numParticles();
872 virts.resize(virts_offset+np);
873 transformParticles(virts, src_tile, 0, virts_offset, np, TransformerVirt());
874 virts_offset += np;
875 }
876 }
877 if (aggregation_type == "Cell")
878 {
879 //Components would be based on
880 int nComp = AMREX_SPACEDIM + NStructReal + NArrayReal;
881 // NArrayReal, NStructInt, NArrayInt behavior as before
882 int nGhost = 0;
883 MultiFab mf(ParticleBoxArray(level), ParticleDistributionMap(level), nComp, nGhost);
884
885 nComp = 1 + NStructInt + NArrayInt;
886 iMultiFab imf(ParticleBoxArray(level), ParticleDistributionMap(level), nComp, nGhost);
887
888 const auto& geom = Geom(level);
889 const auto plo = geom.ProbLoArray();
890 const auto dxi = geom.InvCellSizeArray();
891 const auto domain = geom.Domain();
892
893 BoxList bl_buffer;
894 bl_buffer.complementIn(Geom(level).Domain(), ParticleBoxArray(level));
895 BoxArray buffer(std::move(bl_buffer));
896 buffer.grow(aggregation_buffer);
897
899 locator.build(buffer, geom, do_tiling, tile_size);
900 AssignGrid<DenseBinIteratorFactory<Box>> assign_buffer_grid = locator.getGridAssignor();
901
902 amrex::ParticleToMesh(*this, mf, level,
903 [=] AMREX_GPU_DEVICE (const ParticleType& p,
904 amrex::Array4<amrex::Real> const& partData,
907 {
908 auto iv = getParticleCell(p, plo_loc, dxi_loc, domain);
909 if(assign_buffer_grid(iv).first == -1)
910 {
911 //Ordering may make this not unique
912 for (int i = 0; i < NArrayReal; ++i)
913 {
914 amrex::Gpu::Atomic::AddNoRet(&partData(iv,AMREX_SPACEDIM+NStructReal+i), partData(iv,AMREX_SPACEDIM)!=0.0 ? static_cast<Real>(0) : static_cast<Real>(p.rdata(NStructReal+i)));
915 }
916 //Add the rdata(0)-weighted sum of position
917 for (int i = 0; i < AMREX_SPACEDIM; ++i)
918 {
919 amrex::Gpu::Atomic::AddNoRet(&partData(iv,i), static_cast<Real>((p.rdata(0)*p.pos(i))));
920 }
921 //Add the rdata(0)-weighted sum of other rdata fields
922 for (int i = 1; i < NStructReal; ++i)
923 {
924 amrex::Gpu::Atomic::AddNoRet(&partData(iv,AMREX_SPACEDIM+i), static_cast<Real>((p.rdata(0)*p.rdata(i))));
925 }
926 //Add the rdata(0) sum
927 for (int i = 0; i < 1; ++i)
928 {
929 amrex::Gpu::Atomic::AddNoRet(&partData(iv,AMREX_SPACEDIM+i), static_cast<Real>(p.rdata(0)));
930 }
931 }
932
933 }); //skipping extra false argument, doing mf.setVal(0) at beginning
934
935 amrex::ParticleToMesh(*this, imf, level,
936 [=] AMREX_GPU_DEVICE (const ParticleType& p,
937 amrex::Array4<int> const& partData,
940 {
941
942 auto iv = getParticleCell(p, plo_loc, dxi_loc, domain);
943 if(assign_buffer_grid(iv).first == -1)
944 {
945 //if this cell has no particle id info, do a straight copy to store idata
946 if(partData(iv,0)==0)
947 {
948 //Add 1 to indicate at least 1 particle at cell iv
949 amrex::Gpu::Atomic::AddNoRet(&partData(iv,0), 1);
950 for (int i = 0; i < NStructInt; ++i)
951 {
952 amrex::Gpu::Atomic::AddNoRet(&partData(iv,1+i), p.idata(i));
953 }
954 for (int i = 0; i < NArrayInt; ++i)
955 {
956 amrex::Gpu::Atomic::AddNoRet(&partData(iv,1+NStructInt+i), p.idata(NStructInt+i));
957 }
958 }
959 }
960 });
961
962 //There may be a better way to ensure virts is the right length
963 virts.resize(imf.sum(0));
964
965 int last_offset = 0;
966 for (MFIter mfi(mf); mfi.isValid(); ++mfi)
967 {
968 const auto bx = mfi.tilebox();
969 const auto partData = mf.array(mfi);
970 const auto imf_arr = imf.array(mfi);
971
972 Gpu::DeviceVector<int> offsets(bx.numPts());
973 auto *offsets_ptr = offsets.dataPtr();
974 int next_offset = Scan::ExclusiveSum((int) bx.numPts(),(imf_arr.ptr(bx.smallEnd(),0)),(offsets.dataPtr()),Scan::retSum);
975 auto dst = virts.getParticleTileData();
976 ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
977 {
978 if(imf_arr(i,j,k,0)!=0)
979 {
980 const auto idx = last_offset + offsets_ptr[
981 imf_arr.get_offset(IntVectND<3>{i,j,k})
982 ];
983
984 dst.cpu(idx) = 0;
986
987 auto& p = dst[idx];
988 //Set rdata(0) first so we can normalize the weighted fields
989 //Note that this does not work for soa PC
990 p.rdata(0) = static_cast<ParticleReal>(partData(i,j,k,AMREX_SPACEDIM));;
991 //Set pos with the normalized weighted field
992 for (int n = 0; n < AMREX_SPACEDIM; ++n)
993 {
994 p.pos(n) = static_cast<ParticleReal>(partData(i,j,k,n) / p.rdata(0));
995 }
996 //Set rdata(n>0) with the normalized weighted field for NStructReal
997 for (int n = 1; n < NStructReal; ++n)
998 {
999 p.rdata(n) = static_cast<ParticleReal>(partData(i,j,k,AMREX_SPACEDIM+n) / p.rdata(0));
1000 }
1001 //Set rdata(n>0) with the normalized weighted field for NArrayReal
1002 for (int n = 0; n < NArrayReal; ++n)
1003 {
1004 dst.rdata(n)[idx] = static_cast<ParticleReal>(partData(i,j,k,AMREX_SPACEDIM+NStructReal+n));
1005 }
1006 //Set idata with the "first" particles idata field for NStructInt
1007 for (int n = 0; n < NStructInt; ++n)
1008 {
1009 p.idata(n) = imf_arr(i,j,k,1+n);
1010 }
1011 //Set idata with the "first" particles idata field for NArrayInt
1012 for (int n = 0; n < NArrayInt; ++n)
1013 {
1014 dst.idata(n)[idx] = imf_arr(i,j,k,1+NStructInt+n);
1015 }
1016 }
1017
1018 });
1019 last_offset+=next_offset;
1021 }
1022
1023 // last_offset should equal virts.numParticles()
1024 auto virts_offset = last_offset;
1025 for(ParConstIterType pti(*this, level); pti.isValid(); ++pti)
1026 {
1027 const auto& src_tile = ParticlesAt(level, pti);
1028
1029 auto np = src_tile.numParticles();
1030 virts.resize(virts_offset+np);
1031 virts_offset += filterAndTransformParticles(virts, src_tile, FilterVirt(assign_buffer_grid,plo,dxi,domain), TransformerVirt(),0,virts_offset);
1033 }
1034 virts.resize(virts_offset);
1036 }
1037}
1038
1040{
1041
1044
1050 : m_lev_min(level), m_lev_max(level+1), m_nGrow(nGrow), m_gid(gid), m_assign_grid(assign_grid)
1051 {}
1052
1053 template <typename SrcData>
1055 int operator() (const SrcData& src, int src_i) const noexcept
1056 {
1057 const auto tup_min = (m_assign_grid)(src[src_i], m_lev_min, m_lev_max, m_nGrow, DefaultAssignor{});
1058 const auto tup_max = (m_assign_grid)(src[src_i], m_lev_max, m_lev_max, m_nGrow, DefaultAssignor{});
1059 const auto p_boxes = amrex::get<0>(tup_min);
1060 const auto p_boxes_max = amrex::get<0>(tup_max);
1061 const auto p_levs_max = amrex::get<2>(tup_max);
1062 return p_boxes_max >=0 && p_boxes == m_gid && p_levs_max == m_lev_max;
1063 }
1064};
1065
1067{
1068
1069 template <typename DstData, typename SrcData>
1071 void operator() (DstData& dst, const SrcData& src,
1072 int src_i, int dst_i) const noexcept
1073 {
1074 copyParticle(dst, src, src_i, dst_i);
1075
1076 dst.id(dst_i) = LongParticleIds::GhostParticleID;
1077 dst.cpu(dst_i) = 0;
1078 }
1079};
1080
1081template <typename ParticleType, int NArrayReal, int NArrayInt,
1082 template<class> class Allocator, class CellAssignor>
1083void
1084ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1085::CreateGhostParticles (int level, int nGrow, AoS& ghosts) const
1086{
1087 ParticleTileType ptile;
1088 CreateGhostParticles(level, nGrow, ptile);
1089 ptile.GetArrayOfStructs().swap(ghosts);
1090}
1091
1092template <typename ParticleType, int NArrayReal, int NArrayInt,
1093 template<class> class Allocator, class CellAssignor>
1094void
1096::CreateGhostParticles (int level, int nGrow, ParticleTileType& ghosts) const
1097{
1098 BL_PROFILE("ParticleContainer::CreateGhostParticles()");
1099 AMREX_ASSERT(ghosts.empty());
1100 AMREX_ASSERT(level < finestLevel());
1101
1102 if (level >= std::ssize(m_particles)) {
1103 return;
1104 }
1105
1106 if (! m_particle_locator.isValid(GetParGDB(), do_tiling, tile_size)) {
1107 m_particle_locator.build(GetParGDB(), do_tiling, tile_size);
1108 }
1109
1110 m_particle_locator.setGeometry(GetParGDB());
1111 AmrAssignGrid<DenseBinIteratorFactory<Box>> assign_grid = m_particle_locator.getGridAssignor();
1112 auto ghost_offset = ghosts.numParticles();
1113 for(ParConstIterType pti(*this, level); pti.isValid(); ++pti)
1114 {
1115 const auto& src_tile = ParticlesAt(level, pti);
1116 int gid = pti.index();
1117
1118 auto np = src_tile.numParticles();
1119 ghosts.resize(ghost_offset+np);
1120 ghost_offset += filterAndTransformParticles(ghosts, src_tile, AssignGridFilter(assign_grid,gid,level,nGrow), TransformerGhost(),0,ghost_offset);
1121 }
1122 ghosts.resize(ghost_offset);
1124}
1125
1126template <typename ParticleType, int NArrayReal, int NArrayInt,
1127 template<class> class Allocator, class CellAssignor>
1128void
1131{
1132 BL_PROFILE("ParticleContainer::clearParticles()");
1133
1134 for (int lev = 0; lev < std::ssize(m_particles); ++lev)
1135 {
1136 for (auto& kv : m_particles[lev]) { kv.second.resize(0); }
1137 particle_detail::clearEmptyEntries(m_particles[lev]);
1138 }
1139}
1140
1141template <typename ParticleType, int NArrayReal, int NArrayInt,
1142 template<class> class Allocator, class CellAssignor>
1143template <class PCType>
1145void
1147copyParticles (const PCType& other, bool local)
1148{
1149 using PData = typename ParticleTileType::ConstParticleTileDataType;
1150 copyParticles(other, [] AMREX_GPU_HOST_DEVICE (const PData& /*data*/, int /*i*/) { return 1; }, local);
1151}
1152
1153template <typename ParticleType, int NArrayReal, int NArrayInt,
1154 template<class> class Allocator, class CellAssignor>
1155template <class PCType>
1157void
1159addParticles (const PCType& other, bool local)
1160{
1161 using PData = typename ParticleTileType::ConstParticleTileDataType;
1162 addParticles(other, [] AMREX_GPU_HOST_DEVICE (const PData& /*data*/, int /*i*/) { return 1; }, local);
1163}
1164
1165template <typename ParticleType, int NArrayReal, int NArrayInt,
1166 template<class> class Allocator, class CellAssignor>
1167template <class F, class PCType>
1168requires ((IsParticleContainer<PCType>::value) && (!std::is_integral_v<F>))
1169void
1171copyParticles (const PCType& other, F&& f, bool local)
1172{
1173 BL_PROFILE("ParticleContainer::copyParticles");
1174 clearParticles();
1175 addParticles(other, std::forward<F>(f), local);
1176}
1177
1178template <typename ParticleType, int NArrayReal, int NArrayInt,
1179 template<class> class Allocator, class CellAssignor>
1180template <class F, class PCType>
1181requires ((IsParticleContainer<PCType>::value) && (!std::is_integral_v<F>))
1182void
1184addParticles (const PCType& other, F const& f, bool local)
1185{
1186 BL_PROFILE("ParticleContainer::addParticles");
1187
1188 // touch all tiles in serial
1189 for (int lev = 0; lev < other.numLevels(); ++lev)
1190 {
1191 [[maybe_unused]] Gpu::NoSyncRegion no_sync{};
1192 const auto& plevel_other = other.GetParticles(lev);
1193 for(MFIter mfi = other.MakeMFIter(lev); mfi.isValid(); ++mfi)
1194 {
1195 auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
1196 if(!plevel_other.contains(index)) { continue; }
1197
1198 DefineAndReturnParticleTile(lev, mfi.index(), mfi.LocalTileIndex());
1199 }
1200 }
1201
1202#ifdef AMREX_USE_OMP
1203#pragma omp parallel if (Gpu::notInLaunchRegion())
1204#endif
1205 for (int lev = 0; lev < other.numLevels(); ++lev)
1206 {
1207 const auto& plevel_other = other.GetParticles(lev);
1208 for(MFIter mfi = other.MakeMFIter(lev); mfi.isValid(); ++mfi)
1209 {
1210 auto index = std::make_pair(mfi.index(), mfi.LocalTileIndex());
1211 if(!plevel_other.contains(index)) { continue; }
1212
1213 // this has already had define() called above
1214 auto& ptile = ParticlesAt(lev, mfi.index(), mfi.LocalTileIndex());
1215 const auto& ptile_other = plevel_other.at(index);
1216 auto np = ptile_other.numParticles();
1217 if (np == 0) { continue; }
1218
1219 auto dst_index = ptile.numParticles();
1220 ptile.resize(dst_index + np);
1221
1222 auto count = amrex::filterParticles(ptile, ptile_other, f,
1223 static_cast<decltype(np)>(0), dst_index, np);
1224
1225 ptile.resize(dst_index + count);
1226 }
1227 }
1228
1229 if (! local) { Redistribute(); }
1230}
1231
1232//
1233// This redistributes valid particles and discards invalid ones.
1235template <typename ParticleType, int NArrayReal, int NArrayInt,
1236 template<class> class Allocator, class CellAssignor>
1237void
1239::Redistribute (int lev_min, int lev_max, int nGrow, int local, bool remove_negative)
1240{
1241 BL_PROFILE_SYNC_START_TIMED("SyncBeforeComms: Redist");
1242
1244 !is_rtsoa_pc || NumRuntimeRealComps() >= AMREX_SPACEDIM,
1245 "ParticleContainer with RTSoA requires at least AMREX_SPACEDIM "
1246 "runtime real components for positions"
1247 );
1248
1249 Redistribute_impl(lev_min, lev_max, nGrow, local, remove_negative);
1250
1252}
1253
1254template <typename ParticleType, int NArrayReal, int NArrayInt,
1255 template<class> class Allocator, class CellAssignor>
1256template <class index_type>
1257void
1259::ReorderParticles (int lev, const MFIter& mfi, const index_type* permutations)
1260{
1261 auto& ptile = ParticlesAt(lev, mfi);
1262 const size_t np = ptile.numParticles();
1263 const size_t np_total = np + ptile.numNeighborParticles();
1264
1265 if (memEfficientSort) {
1266 amrex::ReorderParticles(ptile, permutations);
1267 } else {
1268 ParticleTileType ptile_tmp;
1269 ptile_tmp.define(m_num_runtime_real, m_num_runtime_int,
1270 &m_soa_rdata_names, &m_soa_idata_names, arena());
1271 ptile_tmp.resize(np_total);
1272 // copy re-ordered particles
1273 gatherParticles(ptile_tmp, ptile, np, permutations);
1274 // copy neighbor particles
1275 amrex::copyParticles(ptile_tmp, ptile, np, np, np_total-np);
1276 ptile.swap(ptile_tmp);
1277 }
1278}
1279
1280template <typename ParticleType, int NArrayReal, int NArrayInt,
1281 template<class> class Allocator, class CellAssignor>
1282void
1284{
1285 SortParticlesByBin(IntVect(AMREX_D_DECL(1, 1, 1)));
1286}
1287
1288template <typename ParticleType, int NArrayReal, int NArrayInt,
1289 template<class> class Allocator, class CellAssignor>
1290void
1293{
1294 BL_PROFILE("ParticleContainer::SortParticlesByBin()");
1295
1296 if (bin_size == IntVect::TheZeroVector()) { return; }
1297
1298 for (int lev = 0; lev < numLevels(); ++lev)
1299 {
1300 const Geometry& geom = Geom(lev);
1301 const auto dxi = geom.InvCellSizeArray();
1302 const auto plo = geom.ProbLoArray();
1303 const auto domain = geom.Domain();
1304
1305 for(MFIter mfi = MakeMFIter(lev); mfi.isValid(); ++mfi)
1306 {
1307 auto& ptile = ParticlesAt(lev, mfi);
1308 const size_t np = ptile.numParticles();
1309
1310 const Box& box = mfi.validbox();
1311
1312 int ntiles = numTilesInBox(box, true, bin_size);
1313
1314 m_bins.build(np, ptile.getParticleTileData(), ntiles,
1315 GetParticleBin{.plo = plo, .dxi = dxi, .domain = domain,
1316 .bin_size = bin_size, .box = box});
1317 ReorderParticles(lev, mfi, m_bins.permutationPtr());
1318 }
1319 }
1320}
1321
1322template <typename ParticleType, int NArrayReal, int NArrayInt,
1323 template<class> class Allocator, class CellAssignor>
1324void
1327{
1328 BL_PROFILE("ParticleContainer::SortParticlesForDeposition()");
1329
1330 for (int lev = 0; lev < numLevels(); ++lev)
1331 {
1332 const Geometry& geom = Geom(lev);
1333
1334 for(MFIter mfi = MakeMFIter(lev); mfi.isValid(); ++mfi)
1335 {
1336 const auto& ptile = ParticlesAt(lev, mfi);
1337 const size_t np = ptile.numParticles();
1338
1339 const Box& box = mfi.validbox();
1340
1341 using index_type = typename decltype(m_bins)::index_type;
1343 PermutationForDeposition<index_type>(perm, np, ptile, box, geom, idx_type);
1344 ReorderParticles(lev, mfi, perm.dataPtr());
1345 }
1346 }
1347}
1348
1349template <typename ParticleType, int NArrayReal, int NArrayInt,
1350 template<class> class Allocator, class CellAssignor>
1352int
1354::hostPartitionTile (ParticleTileType& src_tile,
1355 int lev, int gid, int tid,
1356 int lev_min, int lev_max, int nGrow, int local,
1357 bool remove_negative, int myproc,
1359 Gpu::DeviceVector<int>& levels,
1361 Gpu::DeviceVector<int>& src_indices,
1362 Gpu::DeviceVector<IntVect>& periodic_shift)
1363{
1364 ParticleLocData pld;
1365 auto ptd = src_tile.getParticleTileData();
1366 int num_move = 0;
1367
1368 Long last = src_tile.numParticles() - 1;
1369 Long pindex = 0;
1370 int who = -1;
1371 while (pindex <= last) {
1372 decltype(auto) p = ptd[pindex];
1373
1374 // remove invalid if needed
1375 if (!ptd.id(pindex).is_valid()) {
1376 if (remove_negative) {
1377 copyParticle(ptd, ptd, last, pindex);
1378 correctCellVectors(last, pindex, gid, p);
1379 --last;
1380 } else {
1381 ++pindex;
1382 }
1383 continue;
1384 }
1385
1386 // Fast path: if on the finest level, try the last good location first
1387 // If not on finest level, we need to do the full locate because even if
1388 // the last assignment works, the particle might belong on a finer level.
1389 if (lev == lev_max && pld.m_tile == tid && pld.m_grid == gid && pld.m_lev == lev && who == myproc) {
1390 const auto iv = Index(p, lev);
1391 if (pld.m_tilebox.contains(iv)) {
1392 ++pindex;
1393 continue; // don't need correctCellVectors or particlePostLocate if it stays where it is
1394 }
1395 }
1396
1397 // full locate
1398 locateParticle(p, pld, lev_min, lev_max, nGrow, local ? gid : -1);
1399 particlePostLocate(p, pld, lev);
1400
1401 // particle may have been invalidated, so check if we need to remove again.
1402 if (!ptd.id(pindex).is_valid()) {
1403 copyParticle(ptd, ptd, last, pindex);
1404 correctCellVectors(last, pindex, gid, p);
1405 --last;
1406 continue;
1407 }
1408
1409 // flag particle for move if needed
1410 who = BufferMap().procID(pld.m_grid, pld.m_tile, pld.m_lev);
1411 if (pld.m_lev != lev || pld.m_grid != gid || pld.m_tile != tid || who != myproc) {
1412 // the particle is valid but needs to be sent somewhere else
1413 swapParticle(ptd, ptd, last, pindex);
1414 correctCellVectors(last, pindex, gid, p);
1415 boxes[num_move] = pld.m_grid;
1416 levels[num_move] = pld.m_lev;
1417 tiles[num_move] = pld.m_tile;
1418 src_indices[num_move] = static_cast<int>(last);
1419 ++num_move;
1420 --last;
1421 continue;
1422 }
1423
1424 ++pindex; // if here, particle stays
1425 }
1426
1427 boxes.resize(num_move);
1428 levels.resize(num_move);
1429 tiles.resize(num_move);
1430 src_indices.resize(num_move);
1431 periodic_shift.resize(num_move);
1432
1433 return static_cast<int>(pindex);
1434}
1435
1436//
1437// Shared implementation of Redistribute
1438//
1439template <typename ParticleType, int NArrayReal, int NArrayInt,
1440 template<class> class Allocator, class CellAssignor>
1441void
1442ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1443::Redistribute_impl (int lev_min, int lev_max, int nGrow, int local, bool remove_negative)
1444{
1445 if (local) { AMREX_ASSERT(numParticlesOutOfRange(*this, lev_min, lev_max, local) == 0); }
1446
1447 BL_PROFILE("ParticleContainer::Redistribute_impl()");
1448 BL_PROFILE_VAR_NS("Redistribute_partition", blp_partition);
1449
1450 int theEffectiveFinestLevel = m_gdb->finestLevel();
1451 while (!m_gdb->LevelDefined(theEffectiveFinestLevel)) { theEffectiveFinestLevel--; }
1452
1453 if (std::ssize(m_particles) < theEffectiveFinestLevel+1) {
1454 if (Verbose()) {
1455 amrex::Print() << "ParticleContainer::Redistribute() resizing containers from "
1456 << m_particles.size() << " to "
1457 << theEffectiveFinestLevel + 1 << '\n';
1458 }
1459 m_particles.resize(theEffectiveFinestLevel+1);
1460 m_dummy_mf.resize(theEffectiveFinestLevel+1);
1461 }
1462
1463 for (int lev = 0; lev < theEffectiveFinestLevel+1; ++lev) { RedefineDummyMF(lev); }
1464
1465 int finest_lev_particles;
1466 if (lev_max == -1) {
1467 lev_max = theEffectiveFinestLevel;
1468 finest_lev_particles = m_particles.size() - 1;
1469 } else {
1470 finest_lev_particles = lev_max;
1471 }
1472 AMREX_ASSERT(lev_max <= finestLevel());
1473
1474 this->defineBufferMap();
1475
1476#ifndef AMREX_USE_GPU
1477 if (local > 0) { BuildRedistributeMask(0, local); }
1478#else
1479 if (! m_particle_locator.isValid(GetParGDB(), do_tiling, tile_size)) { m_particle_locator.build(GetParGDB(), do_tiling, tile_size); }
1480 m_particle_locator.setGeometry(GetParGDB());
1481#endif
1482
1483 BL_PROFILE_VAR_START(blp_partition);
1484 ParticleCopyOp op;
1485 int num_levels = finest_lev_particles + 1;
1486 op.setNumLevels(num_levels);
1487 Vector<std::map<std::pair<int, int>, int> > new_sizes(num_levels);
1488#ifndef AMREX_USE_GPU
1489 const int myproc = ParallelContext::MyProcSub();
1490#endif
1491#ifdef AMREX_USE_GPU
1492 auto assign_grid = m_particle_locator.getGridAssignor();
1493 const auto plo = Geom(0).ProbLoArray();
1494 const auto phi = Geom(0).ProbHiArray();
1495 const auto rlo = Geom(0).ProbLoArrayInParticleReal();
1496 const auto rhi = Geom(0).ProbHiArrayInParticleReal();
1497 const auto is_per = Geom(0).isPeriodicArray();
1498#endif
1499
1500#if defined(AMREX_USE_OMP) || defined(AMREX_USE_GPU)
1501 Vector<std::pair<int, int> > grid_tile_ids;
1503 Vector<int> plevs;
1504 Vector<int*> new_size_ptrs;
1505#ifndef AMREX_USE_GPU
1507 Vector<Gpu::DeviceVector<int>*> level_ptrs;
1510 Vector<Gpu::DeviceVector<IntVect>*> periodic_shift_ptrs;
1511#endif
1512 std::size_t num_tiles = 0;
1513 for (int lev = lev_min; lev <= finest_lev_particles; ++lev) {
1514 for (auto const& kv : m_particles[lev]) {
1515 if (kv.second.numParticles() != 0) {
1516 ++num_tiles;
1517 }
1519 }
1520 grid_tile_ids.reserve(num_tiles);
1521 ptile_ptrs.reserve(num_tiles);
1522 plevs.reserve(num_tiles);
1523 new_size_ptrs.reserve(num_tiles);
1524#ifndef AMREX_USE_GPU
1525 box_ptrs.reserve(num_tiles);
1526 level_ptrs.reserve(num_tiles);
1527 tile_ptrs.reserve(num_tiles);
1528 src_index_ptrs.reserve(num_tiles);
1529 periodic_shift_ptrs.reserve(num_tiles);
1530#endif
1531 for (int lev = lev_min; lev <= finest_lev_particles; lev++) {
1532 auto& pmap = m_particles[lev];
1533 for (auto& kv : pmap)
1535 const auto np = kv.second.numParticles();
1536 if (np == 0) { continue; }
1537
1538 grid_tile_ids.push_back(kv.first);
1539 ptile_ptrs.push_back(&(kv.second));
1540 plevs.push_back(lev);
1541 auto index = std::make_pair(kv.first.first, kv.first.second);
1542 auto& new_size = new_sizes[lev][index];
1543 new_size = 0;
1544 new_size_ptrs.push_back(&new_size);
1545#ifndef AMREX_USE_GPU
1546 op.resize(kv.first.first, kv.first.second, lev, static_cast<int>(np));
1547 auto& boxes = op.m_boxes[lev][index];
1548 auto& levels = op.m_levels[lev][index];
1549 auto& tiles = op.m_tiles[lev][index];
1550 auto& src_indices = op.m_src_indices[lev][index];
1551 auto& periodic_shift = op.m_periodic_shift[lev][index];
1552 box_ptrs.push_back(&boxes);
1553 level_ptrs.push_back(&levels);
1554 tile_ptrs.push_back(&tiles);
1555 src_index_ptrs.push_back(&src_indices);
1556 periodic_shift_ptrs.push_back(&periodic_shift);
1557#endif
1558 }
1559 }
1560#endif
1561
1562#if defined(AMREX_USE_OMP) || defined(AMREX_USE_GPU)
1563#ifdef AMREX_USE_OMP
1564#pragma omp parallel for
1565#endif
1566 for (int pmap_it = 0; pmap_it < static_cast<int>(ptile_ptrs.size()); ++pmap_it)
1567 {
1568 int lev = plevs[pmap_it];
1569 int gid = grid_tile_ids[pmap_it].first;
1570 int tid = grid_tile_ids[pmap_it].second;
1571 auto& src_tile = *ptile_ptrs[pmap_it];
1572 const size_t np = src_tile.numParticles();
1573 if (np == 0) {continue;}
1574
1575#ifndef AMREX_USE_GPU
1576 auto& boxes = *box_ptrs[pmap_it];
1577 auto& levels = *level_ptrs[pmap_it];
1578 auto& tiles = *tile_ptrs[pmap_it];
1579 auto& src_indices = *src_index_ptrs[pmap_it];
1580 auto& periodic_shift = *periodic_shift_ptrs[pmap_it];
1581 *new_size_ptrs[pmap_it] = hostPartitionTile(src_tile,
1582 lev, gid, tid,
1583 lev_min, lev_max, nGrow, local,
1584 remove_negative, myproc,
1585 boxes, levels, tiles,
1586 src_indices, periodic_shift);
1587#else // GPU algorithm
1588 int num_stay = partitionParticlesByDest(src_tile, assign_grid,
1589 std::forward<CellAssignor>(CellAssignor{}),
1590 BufferMap(),
1591 plo, phi, rlo, rhi, is_per, lev, gid, tid,
1592 lev_min, lev_max, nGrow, remove_negative);
1593
1594 int num_move = np - num_stay;
1595 *new_size_ptrs[pmap_it] = num_stay;
1596 op.resize(gid, tid, lev, num_move);
1597
1598 auto index = std::make_pair(gid, tid);
1599 auto* p_boxes = op.m_boxes[lev][index].dataPtr();
1600 auto* p_levs = op.m_levels[lev][index].dataPtr();
1601 auto* p_tiles = op.m_tiles[lev][index].dataPtr();
1602 auto* p_src_indices = op.m_src_indices[lev][index].dataPtr();
1603 auto* p_periodic_shift = op.m_periodic_shift[lev][index].dataPtr();
1604 auto ptd = src_tile.getParticleTileData();
1605
1606 amrex::ParallelFor(num_move, [=] AMREX_GPU_DEVICE (int i)
1607 {
1608 const auto p = ptd[i + num_stay];
1609
1610 if (!p.id().is_valid())
1611 {
1612 p_boxes[i] = -1;
1613 p_tiles[i] = -1;
1614 p_levs[i] = -1;
1615 }
1616 else
1617 {
1618 const auto tup = assign_grid(p, lev_min, lev_max, nGrow,
1619 std::forward<CellAssignor>(CellAssignor{}));
1620 p_boxes[i] = amrex::get<0>(tup);
1621 p_tiles[i] = amrex::get<1>(tup);
1622 p_levs[i] = amrex::get<2>(tup);
1623 }
1624 p_periodic_shift[i] = IntVect(AMREX_D_DECL(0,0,0));
1625 p_src_indices[i] = i+num_stay;
1626 });
1627#endif
1628 }
1629#else
1630 for (int lev = lev_min; lev <= finest_lev_particles; ++lev) {
1631 auto& pmap = m_particles[lev];
1632 for (auto& kv : pmap)
1633 {
1634 auto& src_tile = kv.second;
1635 const auto np = src_tile.numParticles();
1636 if (np == 0) { continue; }
1637
1638 int gid = kv.first.first;
1639 int tid = kv.first.second;
1640 auto index = std::make_pair(gid, tid);
1641 auto& new_size = new_sizes[lev][index];
1642 new_size = 0;
1643
1644 op.resize(gid, tid, lev, static_cast<int>(np));
1645 auto& boxes = op.m_boxes[lev][index];
1646 auto& levels = op.m_levels[lev][index];
1647 auto& tiles = op.m_tiles[lev][index];
1648 auto& src_indices = op.m_src_indices[lev][index];
1649 auto& periodic_shift = op.m_periodic_shift[lev][index];
1650 new_size = hostPartitionTile(src_tile,
1651 lev, gid, tid,
1652 lev_min, lev_max, nGrow, local,
1653 remove_negative, myproc,
1654 boxes, levels, tiles,
1655 src_indices, periodic_shift);
1656 }
1657 }
1658#endif
1659 BL_PROFILE_VAR_STOP(blp_partition);
1660
1661 ParticleCopyPlan plan;
1662
1663 plan.build(*this, op, h_redistribute_int_comp,
1664 h_redistribute_real_comp, local);
1665
1666 // by default, this uses The_Arena();
1669
1670 if (use_comms_arena) {
1671 snd_buffer.setArena(The_Comms_Arena());
1672 rcv_buffer.setArena(The_Comms_Arena());
1673 }
1674
1675 packBuffer(*this, op, plan, snd_buffer);
1676
1677 // clear particles from container
1678 for (int lev = lev_min; lev <= lev_max; ++lev)
1679 {
1680 auto& plev = m_particles[lev];
1681 for (auto& kv : plev)
1682 {
1683 int gid = kv.first.first;
1684 int tid = kv.first.second;
1685 auto index = std::make_pair(gid, tid);
1686 auto& tile = plev[index];
1687 tile.resize(new_sizes[lev][index]);
1688 }
1689 }
1690
1691 for (int lev = lev_min; lev <= lev_max; lev++)
1692 {
1693 particle_detail::clearEmptyEntries(m_particles[lev]);
1694 }
1695
1696 if (std::ssize(m_particles) > theEffectiveFinestLevel+1) {
1697 if (m_verbose > 0) {
1698 amrex::Print() << "ParticleContainer::Redistribute() resizing m_particles from "
1699 << m_particles.size() << " to " << theEffectiveFinestLevel+1 << '\n';
1700 }
1701 AMREX_ASSERT(std::ssize(m_particles) >= 2);
1702
1703 m_particles.resize(theEffectiveFinestLevel + 1);
1704 m_dummy_mf.resize(theEffectiveFinestLevel + 1);
1705 }
1706
1708 {
1709 plan.buildMPIFinish(BufferMap());
1710 communicateParticlesStart(*this, plan, snd_buffer, rcv_buffer);
1711 this->ReserveForRedistribute(plan);
1712 unpackBuffer(*this, plan, snd_buffer, RedistributeUnpackPolicy());
1714 unpackRemotes(*this, plan, rcv_buffer, RedistributeUnpackPolicy());
1715 }
1716 else
1717 {
1719 Gpu::PinnedVector<char> pinned_snd_buffer;
1720 Gpu::PinnedVector<char> pinned_rcv_buffer;
1721
1722 if (snd_buffer.arena()->isPinned()) {
1723 plan.buildMPIFinish(BufferMap());
1725 communicateParticlesStart(*this, plan, snd_buffer, pinned_rcv_buffer);
1726 } else {
1727 pinned_snd_buffer.resize(snd_buffer.size());
1728 Gpu::dtoh_memcpy_async(pinned_snd_buffer.dataPtr(), snd_buffer.dataPtr(), snd_buffer.size());
1729 plan.buildMPIFinish(BufferMap());
1731 communicateParticlesStart(*this, plan, pinned_snd_buffer, pinned_rcv_buffer);
1732 }
1733
1734 this->ReserveForRedistribute(plan);
1735
1736 rcv_buffer.resize(pinned_rcv_buffer.size());
1737 unpackBuffer(*this, plan, snd_buffer, RedistributeUnpackPolicy());
1739 Gpu::htod_memcpy_async(rcv_buffer.dataPtr(), pinned_rcv_buffer.dataPtr(), pinned_rcv_buffer.size());
1740 unpackRemotes(*this, plan, rcv_buffer, RedistributeUnpackPolicy());
1741 }
1742
1744 AMREX_ASSERT(numParticlesOutOfRange(*this, lev_min, lev_max, nGrow) == 0);
1745}
1746
1747template <typename ParticleType, int NArrayReal, int NArrayInt,
1748 template<class> class Allocator, class CellAssignor>
1749void
1750ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1751::ReserveForRedistribute (ParticleCopyPlan const& plan)
1752{
1753 BL_PROFILE("ParticleContainer::ReserveForRedistribute()");
1754
1755 std::map<ParticleTileType*, int> addsizes;
1756
1757 for (int lev = 0; lev < this->BufferMap().numLevels(); ++lev) {
1758 for (MFIter mfi = this->MakeMFIter(lev); mfi.isValid(); ++mfi) {
1759 int gid = mfi.index();
1760 int tid = mfi.LocalTileIndex();
1761 auto& tile = this->DefineAndReturnParticleTile(lev, gid, tid);
1762 int num_copies = plan.m_box_counts_h[this->BufferMap().gridAndTileAndLevToBucket(gid, tid, lev)];
1763 if (num_copies > 0) {
1764 addsizes[&tile] += num_copies;
1765 }
1766 }
1767 }
1768
1769 if (plan.m_nrcvs > 0) {
1770 for (int i = 0; i < std::ssize(plan.m_rcv_box_counts); ++i) {
1771 int copy_size = plan.m_rcv_box_counts[i];
1772 int lev = plan.m_rcv_box_levs[i];
1773 int gid = plan.m_rcv_box_ids[i];
1774 int tid = plan.m_rcv_box_tids[i];
1775 auto& tile = this->DefineAndReturnParticleTile(lev, gid, tid);
1776 addsizes[&tile] += copy_size;
1777 }
1778 }
1779
1780 ParticleTileType::reserve(addsizes);
1781}
1782
1783template <typename ParticleType, int NArrayReal, int NArrayInt,
1784 template<class> class Allocator, class CellAssignor>
1785bool
1787{
1788 BL_PROFILE("ParticleContainer::OK()");
1789
1790 if (lev_max == -1) {
1791 lev_max = finestLevel();
1792 }
1793
1794 return (m_dummy_mf.size() >= lev_max+1 && numParticlesOutOfRange(*this, lev_min, lev_max, nGrow) == 0);
1795}
1796
1797template <typename ParticleType, int NArrayReal, int NArrayInt,
1798 template<class> class Allocator, class CellAssignor>
1799void
1801::AddParticlesAtLevel (AoS& particles, int level, int nGrow)
1802{
1803 ParticleTileType ptile;
1804 ptile.GetArrayOfStructs().swap(particles);
1805 AddParticlesAtLevel(ptile, level, nGrow);
1806}
1807
1808template <typename ParticleType, int NArrayReal, int NArrayInt,
1809 template<class> class Allocator, class CellAssignor>
1810void
1812::AddParticlesAtLevel (ParticleTileType& particles, int level, int nGrow)
1813{
1814 BL_PROFILE("ParticleContainer::AddParticlesAtLevel()");
1815
1816 if (std::ssize(m_particles) < level+1)
1817 {
1818 if (Verbose())
1819 {
1820 amrex::Print() << "ParticleContainer::AddParticlesAtLevel resizing m_particles from "
1821 << m_particles.size()
1822 << " to "
1823 << level+1 << '\n';
1824 }
1825 m_particles.resize(level+1);
1826 m_dummy_mf.resize(level+1);
1827 for (int lev = 0; lev < level+1; ++lev) {
1828 RedefineDummyMF(lev);
1829 }
1830 }
1831
1832 auto& ptile = DefineAndReturnParticleTile(level, 0, 0);
1833 int old_np = ptile.size();
1834 int num_to_add = particles.size();
1835 int new_np = old_np + num_to_add;
1836 ptile.resize(new_np);
1837 amrex::copyParticles(ptile, particles, 0, old_np, num_to_add);
1838 Redistribute(level, level, nGrow);
1839 particles.resize(0);
1840}
1841
1842// This is the single-level version for cell-centered density
1843template <typename ParticleType, int NArrayReal, int NArrayInt,
1844 template<class> class Allocator, class CellAssignor>
1845void
1847AssignCellDensitySingleLevel (int rho_index,
1848 MultiFab& mf_to_be_filled,
1849 int lev,
1850 int ncomp,
1851 int particle_lvl_offset) const
1852{
1853 BL_PROFILE("ParticleContainer::AssignCellDensitySingleLevel()");
1854
1855 if (rho_index != 0) { amrex::Abort("AssignCellDensitySingleLevel only works if rho_index = 0"); }
1856
1857 MultiFab* mf_pointer;
1858
1859 if (OnSameGrids(lev, mf_to_be_filled)) {
1860 // If we are already working with the internal mf defined on the
1861 // particle_box_array, then we just work with this.
1862 mf_pointer = &mf_to_be_filled;
1863 }
1864 else {
1865 // If mf_to_be_filled is not defined on the particle_box_array, then we need
1866 // to make a temporary here and copy into mf_to_be_filled at the end.
1867 mf_pointer = new MultiFab(ParticleBoxArray(lev),
1868 ParticleDistributionMap(lev),
1869 ncomp, mf_to_be_filled.nGrow());
1870 }
1871
1872 // We must have ghost cells for each FAB so that a particle in one grid can spread
1873 // its effect to an adjacent grid by first putting the value into ghost cells of its
1874 // own grid. The mf->SumBoundary call then adds the value from one grid's ghost cell
1875 // to another grid's valid region.
1876 if (mf_pointer->nGrow() < 1) {
1877 amrex::Error("Must have at least one ghost cell when in AssignCellDensitySingleLevel");
1878 }
1879
1880 const auto strttime = amrex::second();
1881
1882 const auto dxi = Geom(lev).InvCellSizeArray();
1883 const auto plo = Geom(lev).ProbLoArray();
1884 const auto pdxi = Geom(lev + particle_lvl_offset).InvCellSizeArray();
1885
1886 if (Geom(lev).isAnyPeriodic() && ! Geom(lev).isAllPeriodic())
1887 {
1888 amrex::Error("AssignCellDensitySingleLevel: problem must be periodic in no or all directions");
1889 }
1890
1891 mf_pointer->setVal(0);
1892
1894#ifdef AMREX_USE_OMP
1895#pragma omp parallel if (Gpu::notInLaunchRegion())
1896#endif
1897 {
1898 FArrayBox local_rho;
1899 for (ParConstIter pti(*this, lev); pti.isValid(); ++pti) {
1900 const Long np = pti.numParticles();
1901 auto ptd = pti.GetParticleTile().getConstParticleTileData();
1902 FArrayBox& fab = (*mf_pointer)[pti];
1903 auto rhoarr = fab.array();
1904#ifdef AMREX_USE_OMP
1905 Box tile_box;
1906 if (Gpu::notInLaunchRegion())
1907 {
1908 tile_box = pti.tilebox();
1909 tile_box.grow(mf_pointer->nGrow());
1910 local_rho.resize(tile_box,ncomp);
1911 local_rho.setVal<RunOn::Host>(0.0);
1912 rhoarr = local_rho.array();
1913 }
1914#endif
1915
1916 if (particle_lvl_offset == 0)
1917 {
1919 {
1920 auto p = ptd[i];
1921 amrex_deposit_cic(p, ncomp, rhoarr, plo, dxi);
1922 });
1923 }
1924 else
1925 {
1927 {
1928 auto p = ptd[i];
1929 amrex_deposit_particle_dx_cic(p, ncomp, rhoarr, plo, dxi, pdxi);
1930 });
1931 }
1932
1933#ifdef AMREX_USE_OMP
1934 if (Gpu::notInLaunchRegion())
1935 {
1936 fab.atomicAdd<RunOn::Host>(local_rho, tile_box, tile_box, 0, 0, ncomp);
1937 }
1938#endif
1939 }
1940 }
1941
1942 mf_pointer->SumBoundary(Geom(lev).periodicity());
1943
1944 // If ncomp > 1, first divide the momenta (component n)
1945 // by the mass (component 0) in order to get velocities.
1946 // Be careful not to divide by zero.
1947 for (int n = 1; n < ncomp; n++)
1948 {
1949 for (MFIter mfi(*mf_pointer); mfi.isValid(); ++mfi)
1950 {
1951 (*mf_pointer)[mfi].protected_divide<RunOn::Device>((*mf_pointer)[mfi],0,n,1);
1952 }
1953 }
1954
1955 // Only multiply the first component by (1/vol) because this converts mass
1956 // to density. If there are additional components (like velocity), we don't
1957 // want to divide those by volume.
1958 const Real* dx = Geom(lev).CellSize();
1959 const Real vol = AMREX_D_TERM(dx[0], *dx[1], *dx[2]);
1960
1961 mf_pointer->mult(Real(1.0)/vol, 0, 1, mf_pointer->nGrow());
1962
1963 // If mf_to_be_filled is not defined on the particle_box_array, then we need
1964 // to copy here from mf_pointer into mf_to_be_filled.
1965 if (mf_pointer != &mf_to_be_filled)
1966 {
1967 mf_to_be_filled.ParallelCopy(*mf_pointer,0,0,ncomp,0,0);
1968 delete mf_pointer;
1969 }
1970
1971 if (m_verbose > 1)
1972 {
1973 auto stoptime = amrex::second() - strttime;
1974
1975 ParallelReduce::Max(stoptime, ParallelContext::IOProcessorNumberSub(),
1976 ParallelContext::CommunicatorSub());
1977
1978 amrex::Print() << "ParticleContainer::AssignCellDensitySingleLevel) time: "
1979 << stoptime << '\n';
1980 }
1981}
1982
1983template <typename ParticleType, int NArrayReal, int NArrayInt,
1984 template<class> class Allocator, class CellAssignor>
1985void
1987ResizeRuntimeRealComp (int new_size, bool communicate)
1988{
1989 int old_size = m_num_runtime_real;
1990
1991 m_runtime_comps_defined = (new_size > 0);
1992 m_num_runtime_real = new_size;
1993 int cur_size = h_redistribute_real_comp.size();
1994 h_redistribute_real_comp.resize(cur_size-old_size+new_size, communicate);
1995 SetParticleSize();
1996
1997 for (int lev = 0; lev < numLevels(); ++lev) {
1998 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
1999 auto& tile = DefineAndReturnParticleTile(lev, pti);
2000 auto np = tile.numParticles();
2001 if (np > 0 && new_size > old_size) {
2002 auto& soa = tile.GetStructOfArrays();
2003 soa.resize(np);
2004 }
2005 }
2006 }
2007}
2008
2009template <typename ParticleType, int NArrayReal, int NArrayInt,
2010 template<class> class Allocator, class CellAssignor>
2011void
2013ResizeRuntimeIntComp (int new_size, bool communicate)
2014{
2015 int old_size = m_num_runtime_int;
2016
2017 m_runtime_comps_defined = (new_size > 0);
2018 m_num_runtime_int = new_size;
2019 int cur_size = h_redistribute_int_comp.size();
2020 h_redistribute_int_comp.resize(cur_size-old_size+new_size, communicate);
2021 SetParticleSize();
2022
2023 for (int lev = 0; lev < numLevels(); ++lev) {
2024 for (ParIterType pti(*this,lev); pti.isValid(); ++pti) {
2025 auto& tile = DefineAndReturnParticleTile(lev, pti);
2026 auto np = tile.numParticles();
2027 if (np > 0 && new_size > old_size) {
2028 auto& soa = tile.GetStructOfArrays();
2029 soa.resize(np);
2030 }
2031 }
2032 }
2033}
2034
2035}
#define BL_PROFILE_VAR_START(vname)
Definition AMReX_BLProfiler.H:562
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_PROFILE_VAR_STOP(vname)
Definition AMReX_BLProfiler.H:563
#define BL_PROFILE_SYNC_STOP()
Definition AMReX_BLProfiler.H:645
#define BL_PROFILE_SYNC_START_TIMED(fname)
Definition AMReX_BLProfiler.H:644
#define BL_PROFILE_VAR_NS(fname, vname)
Definition AMReX_BLProfiler.H:561
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_HOST_DEVICE_FOR_1D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:105
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Print on all processors of the default communicator.
Definition AMReX_Print.H:113
BaseFab< T > & atomicAdd(const BaseFab< T > &x) noexcept
Atomic FAB addition (a[i] <- a[i] + b[i]).
Definition AMReX_BaseFab.H:2487
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:387
void setVal(T const &x, const Box &bx, int dcomp, int ncomp) noexcept
The setVal functions set sub-regions in the BaseFab to a constant value. This most general form speci...
Definition AMReX_BaseFab.H:1404
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:854
BoxArray & grow(int n)
Grow each Box in the BoxArray by the specified amount.
Definition AMReX_BoxArray.cpp:706
std::vector< std::pair< int, Box > > intersections(const Box &bx) const
Return intersections of Box and BoxArray.
Definition AMReX_BoxArray.cpp:1186
Box getCellCenteredBox(int index) const noexcept
Return cell-centered box at element index of this BoxArray.
Definition AMReX_BoxArray.H:744
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:837
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
BoxList & complementIn(const Box &b, const BoxList &bl)
Definition AMReX_BoxList.cpp:307
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:649
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:216
GpuArray< Real, 3 > InvCellSizeArray() const noexcept
Definition AMReX_CoordSys.H:87
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
void resize(const Box &b, int N=1, Arena *ar=nullptr)
For debugging purposes we hide BaseFab version and do some extra work.
Definition AMReX_FArrayBox.cpp:178
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:110
int nGrow(int direction=0) const noexcept
Return the grow factor that defines the region of definition.
Definition AMReX_FabArrayBase.H:78
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:131
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:95
void ParallelCopy(const FabArray< FAB > &src, const Periodicity &period=Periodicity::NonPeriodic(), CpOp op=FabArrayBase::COPY)
Definition AMReX_FabArray.H:873
void setVal(value_type val)
Set all components in the entire region of each FAB to val.
Definition AMReX_FabArray.H:2705
Array4< typename FabArray< FAB >::value_type const > array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:561
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:216
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:192
static void streamSynchronize() noexcept
Definition AMReX_GpuDevice.cpp:856
__host__ __device__ bool cellCentered() const noexcept
True if the IndexTypeND is CELL based in all directions.
Definition AMReX_IndexType.H:102
__host__ static __device__ constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:771
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
void mult(Real val, int comp, int num_comp, int nghost=0)
Scales the value of each cell in the specified subregion of the MultiFab by the scalar val (a[i] <- a...
Definition AMReX_MultiFab.cpp:1417
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
size_type size() const noexcept
Definition AMReX_PODVector.H:648
void resize(size_type a_new_size, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_PODVector.H:728
T * dataPtr() noexcept
Definition AMReX_PODVector.H:670
Definition AMReX_ParIter.H:147
Definition AMReX_ParIter.H:118
int queryAdd(std::string_view name, T &ref)
If name is found, the value in the ParmParse database will be stored in the ref argument....
Definition AMReX_ParmParse.H:1047
int queryarr(std::string_view name, std::vector< int > &ref, int start_ix=FIRST, int num_val=ALL) const
Same as queryktharr() but searches for last occurrence of name.
Definition AMReX_ParmParse.cpp:2034
int query(std::string_view name, bool &ref, int ival=FIRST) const
Same as querykth() but searches for the last occurrence of name.
Definition AMReX_ParmParse.cpp:1948
virtual void reserveData()
Definition AMReX_ParticleContainerBase.cpp:41
virtual void resizeData()
Definition AMReX_ParticleContainerBase.cpp:46
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:149
IntVect Index(const P &p, int lev) const
Definition AMReX_ParticleContainerI.H:202
std::map< std::pair< int, int >, ParticleTileType > ParticleLevel
Definition AMReX_ParticleContainer.H:196
void SetSoACompileTimeNames(std::vector< std::string > const &rdata_name, std::vector< std::string > const &idata_name)
Definition AMReX_ParticleContainerI.H:110
typename ParticleTileType::AoS AoS
Definition AMReX_ParticleContainer.H:199
void SetParticleSize()
Definition AMReX_ParticleContainerI.H:17
void RemoveParticlesAtLevel(int level)
The Following methods are for managing Virtual and Ghost Particles.
Definition AMReX_ParticleContainerI.H:759
void clearParticles()
Clear all the particles in this container. This does not free memory.
Definition AMReX_ParticleContainerI.H:1130
void Increment(MultiFab &mf, int level)
Definition AMReX_ParticleContainerI.H:722
bool HasIntComp(std::string const &name)
Definition AMReX_ParticleContainerI.H:149
void ShrinkToFit()
Definition AMReX_ParticleContainerI.H:702
bool HasRealComp(std::string const &name)
Definition AMReX_ParticleContainerI.H:141
void resizeData() override
This resizes the vector of dummy MultiFabs used by the ParticleContainer for the current number of le...
Definition AMReX_ParticleContainerI.H:436
Long IncrementWithTotal(MultiFab &mf, int level, bool local=false)
Definition AMReX_ParticleContainerI.H:749
Long NumberOfParticlesAtLevel(int level, bool only_valid=true, bool only_local=false) const
Returns # of particles at specified the level.
Definition AMReX_ParticleContainerI.H:552
void SortParticlesByCell()
Sort the particles on each tile by cell, using Fortran ordering.
Definition AMReX_ParticleContainerI.H:1283
int GetRealCompIndex(std::string const &name)
Definition AMReX_ParticleContainerI.H:162
int GetIntCompIndex(std::string const &name)
Definition AMReX_ParticleContainerI.H:184
void RemoveParticlesNotAtFinestLevel()
Definition AMReX_ParticleContainerI.H:773
Vector< Long > NumberOfParticlesInGrid(int level, bool only_valid=true, bool only_local=false) const
Definition AMReX_ParticleContainerI.H:497
ParticleLocData Reset(ParticleType &prt, bool update, bool verbose=true, ParticleLocData pld=ParticleLocData()) const
Updates a particle's location (Where), tries to periodic shift any particles that have left the domai...
Definition AMReX_ParticleContainerI.H:393
std::conditional_t< is_rtsoa_pc, ParticleTileRT< typename ParticleType::RealType, typename ParticleType::IntType >, ParticleTile< ParticleType, NArrayReal, NArrayInt, Allocator > > ParticleTileType
Definition AMReX_ParticleContainer.H:191
Long TotalNumberOfParticles(bool only_valid=true, bool only_local=false) const
Returns # of particles at all levels.
Definition AMReX_ParticleContainerI.H:482
void reserveData() override
This reserves data in the vector of dummy MultiFabs used by the ParticleContainer for the maximum num...
Definition AMReX_ParticleContainerI.H:427
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:151
Definition AMReX_ParticleLocator.H:123
void build(const BoxArray &ba, const Geometry &geom, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000))
Definition AMReX_ParticleLocator.H:130
AssignGrid< BinIteratorFactory > getGridAssignor() const noexcept
Definition AMReX_ParticleLocator.H:206
This class provides the user with a few print options.
Definition AMReX_Print.H:35
Definition AMReX_Reduce.H:438
Type value()
Definition AMReX_Reduce.H:473
Definition AMReX_Reduce.H:597
void eval(MF const &mf, IntVect const &nghost, D &reduce_data, F &&f)
Definition AMReX_Reduce.H:731
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:29
Long size() const noexcept
Definition AMReX_Vector.H:54
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
Long sum(int comp, int nghost=0, bool local=false) const
Returns the sum in component comp.
Definition AMReX_iMultiFab.cpp:454
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
T ExclusiveSum(N n, T const *in, T *out, RetSum a_ret_sum=retSum)
Exclusive sum.
Definition AMReX_Scan.H:1040
void SumBoundary(const Periodicity &period=Periodicity::NonPeriodic(), bool deterministic=false)
Sum values in overlapped cells.
Definition AMReX_FabArray.H:3655
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1289
Arena * The_Comms_Arena()
Definition AMReX_Arena.cpp:880
void Min(KeyValuePair< K, V > &vi, int root, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:287
void Max(KeyValuePair< K, V > &vi, int root, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:254
void Sum(T &v, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:221
void Sum(T &v, int root, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:352
__host__ __device__ AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
Definition AMReX_GpuAtomic.H:283
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void dtoh_memcpy_async(void *p_h, const void *p_d, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:435
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:421
void QueueReduction(Func f)
Definition AMReX_Lazy.cpp:7
constexpr Long GhostParticleID
Definition AMReX_Particle.H:19
constexpr Long VirtualParticleID
Definition AMReX_Particle.H:20
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
int MyProcSub() noexcept
my sub-rank in current frame
Definition AMReX_ParallelContext.H:76
int IOProcessorNumberSub() noexcept
IO sub-rank in current frame.
Definition AMReX_ParallelContext.H:78
bool UseGpuAwareMpi()
Definition AMReX_ParallelDescriptor.H:113
void Bcast(void *, int, MPI_Datatype, int, MPI_Comm)
Definition AMReX_ParallelDescriptor.cpp:1295
void GatherLayoutDataToVector(const LayoutData< T > &sendbuf, Vector< T > &recvbuf, int root)
Gather LayoutData values to a vector on root.
Definition AMReX_ParallelDescriptor.H:1277
static constexpr RetSum retSum
Definition AMReX_Scan.H:33
Definition AMReX_Amr.cpp:50
void gatherParticles(PTile &dst, const PTile &src, N np, const Index *inds)
Gather particles copies particles into contiguous order from an arbitrary order. Specifically,...
Definition AMReX_ParticleTransformation.H:730
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2852
__host__ __device__ void swapParticle(const ParticleTileData< T_ParticleType, NAR, NAI > &dst, const ParticleTileData< T_ParticleType, NAR, NAI > &src, int src_i, int dst_i) noexcept
A general single particle swapping routine that can run on the GPU.
Definition AMReX_ParticleTransformation.H:120
__host__ __device__ int getTileIndex(const IntVect &iv, const Box &box, const bool a_do_tiling, const IntVect &a_tile_size, Box &tbx)
Definition AMReX_ParticleUtil.H:191
__host__ __device__ void copyParticle(const ParticleTileData< T_ParticleType, NAR, NAI > &dst, const ConstParticleTileData< T_ParticleType, NAR, NAI > &src, int src_i, int dst_i) noexcept
A general single particle copying routine that can run on the GPU.
Definition AMReX_ParticleTransformation.H:32
void copyParticles(DstTile &dst, const SrcTile &src) noexcept
Copy particles from src to dst. This version copies all the particles, writing them to the beginning ...
Definition AMReX_ParticleTransformation.H:222
void communicateParticlesStart(const PC &pc, ParticleCopyPlan &plan, const SndBuffer &snd_buffer, RcvBuffer &rcv_buffer)
Definition AMReX_ParticleCommunication.H:909
int partitionParticlesByDest(PTile &ptile, const PLocator &ploc, CellAssignor const &assignor, const ParticleBufferMap &pmap, const GpuArray< Real, 3 > &plo, const GpuArray< Real, 3 > &phi, const GpuArray< ParticleReal, 3 > &rlo, const GpuArray< ParticleReal, 3 > &rhi, const GpuArray< int, 3 > &is_per, int lev, int gid, int tid, int lev_min, int lev_max, int nGrow, bool remove_negative)
Definition AMReX_ParticleUtil.H:654
void ReorderParticles(PTile &ptile, const index_type *permutations)
Reorder particles on the tile ptile using a the permutations array.
Definition AMReX_ParticleUtil.H:953
__host__ __device__ bool enforcePeriodic(P &p, amrex::GpuArray< amrex::Real, 3 > const &plo, amrex::GpuArray< amrex::Real, 3 > const &phi, amrex::GpuArray< amrex::ParticleReal, 3 > const &rlo, amrex::GpuArray< amrex::ParticleReal, 3 > const &rhi, amrex::GpuArray< int, 3 > const &is_per) noexcept
Definition AMReX_ParticleUtil.H:423
void ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:202
__host__ __device__ int numTilesInBox(const Box &box, const bool a_do_tiling, const IntVect &a_tile_size)
Definition AMReX_ParticleUtil.H:239
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
void ParticleToMesh(PC const &pc, const Vector< MultiFab * > &mf, int lev_min, int lev_max, F &&f, bool zero_out_input=true, bool vol_weight=true)
Deposit particles onto a hierarchy of MultiFabs.
Definition AMReX_AmrParticles.H:189
void unpackBuffer(PC &pc, const ParticleCopyPlan &plan, const Buffer &snd_buffer, UnpackPolicy const &policy)
Definition AMReX_ParticleCommunication.H:782
double second() noexcept
Definition AMReX_Utility.cpp:940
void packBuffer(const PC &pc, const ParticleCopyOp &op, const ParticleCopyPlan &plan, Buffer &snd_buffer)
Definition AMReX_ParticleCommunication.H:582
void communicateParticlesFinish(const ParticleCopyPlan &plan)
Definition AMReX_ParticleCommunication.cpp:445
Index filterAndTransformParticles(DstTile &dst, const SrcTile &src, Index *mask, F const &f, Index src_start, Index dst_start) noexcept
Conditionally copy particles from src to dst based on the value of mask. A transformation will also b...
Definition AMReX_ParticleTransformation.H:513
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void unpackRemotes(PC &pc, const ParticleCopyPlan &plan, Buffer &rcv_buffer, UnpackPolicy const &policy)
Definition AMReX_ParticleCommunication.H:1009
Index filterParticles(DstTile &dst, const SrcTile &src, const Index *mask) noexcept
Conditionally copy particles from src to dst based on the value of mask.
Definition AMReX_ParticleTransformation.H:389
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:235
int numParticlesOutOfRange(Iterator const &pti, int nGrow)
Returns the number of particles that are more than nGrow cells from the box correspond to the input i...
Definition AMReX_ParticleUtil.H:36
int Verbose() noexcept
Definition AMReX.cpp:181
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
const int[]
Definition AMReX_BLProfiler.cpp:1664
void transformParticles(DstTile &dst, const SrcTile &src, F &&f) noexcept
Apply the function f to all the particles in src, writing the result to dst. This version does all th...
Definition AMReX_ParticleTransformation.H:273
__host__ __device__ IntVect getParticleCell(P const &p, amrex::GpuArray< amrex::Real, 3 > const &plo, amrex::GpuArray< amrex::Real, 3 > const &dxi) noexcept
Returns the cell index for a given particle using the provided lower bounds and cell sizes.
Definition AMReX_ParticleUtil.H:343
Definition AMReX_ParticleLocator.H:248
A multidimensional array accessor.
Definition AMReX_Array4.H:285
Definition AMReX_ParticleContainerI.H:1040
amrex::AmrAssignGrid< amrex::DenseBinIteratorFactory< amrex::Box > > m_assign_grid
Definition AMReX_ParticleContainerI.H:1043
int m_lev_max
Definition AMReX_ParticleContainerI.H:1042
AssignGridFilter(amrex::AmrAssignGrid< amrex::DenseBinIteratorFactory< amrex::Box > > assign_grid, int gid, int level, int nGrow)
This filters based on matching grids.
Definition AMReX_ParticleContainerI.H:1049
int m_nGrow
Definition AMReX_ParticleContainerI.H:1042
int m_lev_min
Definition AMReX_ParticleContainerI.H:1042
AMREX_GPU_HOST_DEVICE int operator()(const SrcData &src, int src_i) const noexcept
Definition AMReX_ParticleContainerI.H:1055
int m_gid
Definition AMReX_ParticleContainerI.H:1042
Definition AMReX_ParticleLocator.H:17
Definition AMReX_ParticleUtil.H:396
Definition AMReX_DenseBins.H:32
Definition AMReX_ParticleContainerI.H:801
Box m_domain
Definition AMReX_ParticleContainerI.H:805
GpuArray< Real, 3 > m_dxi
Definition AMReX_ParticleContainerI.H:804
FilterVirt(const amrex::AssignGrid< amrex::DenseBinIteratorFactory< amrex::Box > > &assign_buffer_grid, const GpuArray< Real, 3 > &plo, const GpuArray< Real, 3 > &dxi, const Box &domain)
Definition AMReX_ParticleContainerI.H:807
AMREX_GPU_HOST_DEVICE int operator()(const SrcData &src, int src_i) const noexcept
Definition AMReX_ParticleContainerI.H:814
amrex::AssignGrid< amrex::DenseBinIteratorFactory< amrex::Box > > m_assign_buffer_grid
Definition AMReX_ParticleContainerI.H:803
GpuArray< Real, 3 > m_plo
Definition AMReX_ParticleContainerI.H:804
Definition AMReX_ParticleUtil.H:310
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_GpuControl.H:180
Definition AMReX_TypeTraits.H:90
Definition AMReX_ParticleCommunication.H:92
Vector< int > m_rcv_box_ids
Definition AMReX_ParticleCommunication.H:366
Vector< int > m_rcv_box_counts
Definition AMReX_ParticleCommunication.H:364
Vector< int > m_rcv_box_levs
Definition AMReX_ParticleCommunication.H:369
Vector< int > m_rcv_box_tids
Definition AMReX_ParticleCommunication.H:367
int m_nrcvs
Definition AMReX_ParticleCommunication.H:372
Gpu::HostVector< unsigned int > m_box_counts_h
Definition AMReX_ParticleCommunication.H:361
A struct used for storing a particle's position in the AMR hierarchy.
Definition AMReX_ParticleContainer.H:93
Box m_grown_gridbox
Definition AMReX_ParticleContainer.H:100
IntVect m_cell
Definition AMReX_ParticleContainer.H:97
int m_grid
Definition AMReX_ParticleContainer.H:95
int m_tile
Definition AMReX_ParticleContainer.H:96
int m_lev
Definition AMReX_ParticleContainer.H:94
Box m_tilebox
Definition AMReX_ParticleContainer.H:99
Box m_gridbox
Definition AMReX_ParticleContainer.H:98
The struct used to store particles.
Definition AMReX_Particle.H:405
__host__ __device__ RealVect pos() const &
Definition AMReX_Particle.H:456
Definition AMReX_ParticleContainerI.H:1067
AMREX_GPU_HOST_DEVICE void operator()(DstData &dst, const SrcData &src, int src_i, int dst_i) const noexcept
Definition AMReX_ParticleContainerI.H:1071
Definition AMReX_ParticleContainerI.H:822
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(DstData &dst, const SrcData &src, int src_i, int dst_i) const noexcept
Definition AMReX_ParticleContainerI.H:825