Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleUtil.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLEUTIL_H_
2#define AMREX_PARTICLEUTIL_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_IntVect.H>
6#include <AMReX_Box.H>
7#include <AMReX_Gpu.H>
8#include <AMReX_Print.H>
10#include <AMReX_Math.H>
11#include <AMReX_MFIter.H>
12#include <AMReX_ParGDB.H>
13#include <AMReX_ParticleTile.H>
16#include <AMReX_TypeTraits.H>
17#include <AMReX_Scan.H>
18
19#include <limits>
20
21namespace amrex {
22
33template <class Iterator>
34requires (IsParticleIterator<Iterator>::value)
35int
36numParticlesOutOfRange (Iterator const& pti, int nGrow)
37{
38 return numParticlesOutOfRange(pti,
39 IntVect(AMREX_D_DECL(nGrow, nGrow, nGrow)));
40}
41
52template <class Iterator>
53requires (IsParticleIterator<Iterator>::value)
54int
55numParticlesOutOfRange (Iterator const& pti, IntVect nGrow)
56{
57 const auto& tile = pti.GetParticleTile();
58 const auto np = tile.numParticles();
59 const auto& ptd = tile.getConstParticleTileData();
60 const auto& geom = pti.Geom(pti.GetLevel());
61
62 const auto& domain = geom.Domain();
63 const auto& plo = geom.ProbLoArray();
64 const auto& dxi = geom.InvCellSizeArray();
65
66 Box box = pti.tilebox();
67 box.grow(nGrow);
68
69 ReduceOps<ReduceOpSum> reduce_op;
70 ReduceData<int> reduce_data(reduce_op);
71 using ReduceTuple = typename decltype(reduce_data)::Type;
72
73 reduce_op.eval(np, reduce_data,
74 [=] AMREX_GPU_DEVICE (int i) -> ReduceTuple
75 {
76 auto p = ptd[i];
77 if (!p.id().is_valid()) { return false; }
78 using AssignorType = typename Iterator::CellAssignor;
79 AssignorType assignor;
80 IntVect iv = assignor(p, plo, dxi, domain);
81 return !box.contains(iv);
82 });
83 int hv = amrex::get<0>(reduce_data.value(reduce_op));
84 return hv;
85}
86
99template <class PC>
100requires (IsParticleContainer<PC>::value)
101Long
102numParticlesOutOfRange (PC const& pc, int nGrow)
103{
104 return numParticlesOutOfRange(pc, 0, pc.finestLevel(), nGrow);
105}
106
119template <class PC>
120requires (IsParticleContainer<PC>::value)
121Long
122numParticlesOutOfRange (PC const& pc, IntVect nGrow)
123{
124 return numParticlesOutOfRange(pc, 0, pc.finestLevel(), nGrow);
125}
126
141template <class PC>
142requires (IsParticleContainer<PC>::value)
143Long
144numParticlesOutOfRange (PC const& pc, int lev_min, int lev_max, int nGrow)
145{
146 BL_PROFILE("numParticlesOutOfRange()");
147
148 return numParticlesOutOfRange(pc, lev_min, lev_max,
149 IntVect(AMREX_D_DECL(nGrow, nGrow, nGrow)));
150}
151
166template <class PC>
167requires (IsParticleContainer<PC>::value)
168Long
169numParticlesOutOfRange (PC const& pc, int lev_min, int lev_max, IntVect nGrow)
170{
171 BL_PROFILE("numParticlesOutOfRange()");
172
173 using ParIter = typename PC::ParConstIterType;
174 Long num_wrong = 0;
175 for (int lev = lev_min; lev <= lev_max; ++lev)
176 {
177#ifdef AMREX_USE_OMP
178#pragma omp parallel if (Gpu::notInLaunchRegion() && !system::regtest_reduction) reduction(+:num_wrong)
179#endif
180 for(ParIter pti(pc, lev); pti.isValid(); ++pti)
181 {
182 num_wrong += numParticlesOutOfRange(pti, nGrow);
183 }
184 }
186
187 return num_wrong;
188}
189
191int getTileIndex (const IntVect& iv, const Box& box, const bool a_do_tiling,
192 const IntVect& a_tile_size, Box& tbx)
193{
194 if (a_do_tiling == false) {
195 tbx = box;
196 return 0;
197 } else {
198 //
199 // This function must be consistent with FabArrayBase::buildTileArray function!!!
200 //
201 auto tiling_1d = [](int i, int lo, int hi, int tilesize,
202 int& ntile, int& tileidx, int& tlo, int& thi) {
203 int ncells = hi-lo+1;
204 ntile = amrex::max(ncells/tilesize, 1);
205 int ts_right = ncells/ntile;
206 int ts_left = ts_right+1;
207 int nleft = ncells - ntile*ts_right;
208 int ii = i - lo;
209 int nbndry = nleft*ts_left;
210 if (ii < nbndry) {
211 tileidx = ii / ts_left; // tiles on the left of nbndry have size of ts_left
212 tlo = lo + tileidx * ts_left;
213 thi = tlo + ts_left - 1;
214 } else {
215 tileidx = nleft + (ii-nbndry) / ts_right; // tiles on the right: ts_right
216 tlo = lo + tileidx * ts_right + nleft;
217 thi = tlo + ts_right - 1;
218 }
219 };
220 const IntVect& sml = box.smallEnd();
221 const IntVect& big = box.bigEnd();
222 IntVect ntiles, ivIndex, tilelo, tilehi;
223
224 AMREX_D_TERM(int iv0 = amrex::min(amrex::max(iv[0], sml[0]), big[0]);,
225 int iv1 = amrex::min(amrex::max(iv[1], sml[1]), big[1]);,
226 int iv2 = amrex::min(amrex::max(iv[2], sml[2]), big[2]););
227
228 AMREX_D_TERM(tiling_1d(iv0, sml[0], big[0], a_tile_size[0], ntiles[0], ivIndex[0], tilelo[0], tilehi[0]);,
229 tiling_1d(iv1, sml[1], big[1], a_tile_size[1], ntiles[1], ivIndex[1], tilelo[1], tilehi[1]);,
230 tiling_1d(iv2, sml[2], big[2], a_tile_size[2], ntiles[2], ivIndex[2], tilelo[2], tilehi[2]););
231
232 tbx = Box(tilelo, tilehi);
233
234 return AMREX_D_TERM(ivIndex[0], + ntiles[0]*ivIndex[1], + ntiles[0]*ntiles[1]*ivIndex[2]);
235 }
236}
237
239int numTilesInBox (const Box& box, const bool a_do_tiling, const IntVect& a_tile_size)
240{
241 if (a_do_tiling == false) {
242 return 1;
243 } else {
244 //
245 // This function must be consistent with FabArrayBase::buildTileArray function!!!
246 //
247 auto tiling_1d = [](int lo, int hi, int tilesize, int& ntile) {
248 int ncells = hi-lo+1;
249 ntile = amrex::max(ncells/tilesize, 1);
250 };
251
252 const IntVect& sml = box.smallEnd();
253 const IntVect& big = box.bigEnd();
254 IntVect ntiles;
255
256 AMREX_D_TERM(tiling_1d(sml[0], big[0], a_tile_size[0], ntiles[0]);,
257 tiling_1d(sml[1], big[1], a_tile_size[1], ntiles[1]);,
258 tiling_1d(sml[2], big[2], a_tile_size[2], ntiles[2]););
259
260 return AMREX_D_TERM(ntiles[0], *=ntiles[1], *=ntiles[2]);
261 }
262}
263
265{
266 BinMapper(const int* off_bins_p,
269 const Dim3* lo_p,
270 const Dim3* hi_p,
271 int* bin_type_array=nullptr)
272 : m_off_bins_p(off_bins_p), m_dxi_p(dxi_p), m_plo_p(plo_p) ,
273 m_lo_p(lo_p) , m_hi_p(hi_p) , m_bin_type_array(bin_type_array) {}
274
275 template <typename T>
277 unsigned int operator() (const T& ptd, int i) const
278 {
279 auto p = ptd[i];
280 int type = (m_bin_type_array) ? m_bin_type_array[i] : 0;
281 int offset = m_off_bins_p[type];
282
283 AMREX_D_TERM(AMREX_ASSERT((p.pos(0)-m_plo_p[type][0])*m_dxi_p[type][0] - m_lo_p[type].x >= 0.0);,
284 AMREX_ASSERT((p.pos(1)-m_plo_p[type][1])*m_dxi_p[type][1] - m_lo_p[type].y >= 0.0);,
285 AMREX_ASSERT((p.pos(2)-m_plo_p[type][2])*m_dxi_p[type][2] - m_lo_p[type].z >= 0.0));
286
287 auto iv = IntVect(AMREX_D_DECL(static_cast<int>(amrex::Math::floor((p.pos(0)-m_plo_p[type][0])*m_dxi_p[type][0])) - m_lo_p[type].x,
288 static_cast<int>(amrex::Math::floor((p.pos(1)-m_plo_p[type][1])*m_dxi_p[type][1])) - m_lo_p[type].y,
289 static_cast<int>(amrex::Math::floor((p.pos(2)-m_plo_p[type][2])*m_dxi_p[type][2])) - m_lo_p[type].z));
290 auto iv3 = iv.dim3();
291 int nx = m_hi_p[type].x-m_lo_p[type].x+1;
292 int ny = m_hi_p[type].y-m_lo_p[type].y+1;
293 int nz = m_hi_p[type].z-m_lo_p[type].z+1;
294 int uix = amrex::min(nx-1,amrex::max(0,iv3.x));
295 int uiy = amrex::min(ny-1,amrex::max(0,iv3.y));
296 int uiz = amrex::min(nz-1,amrex::max(0,iv3.z));
297 return static_cast<unsigned int>( (uiz * ny + uiy) * nx + uix + offset );
298 }
299
300private:
301 const int* m_off_bins_p;
302 const GpuArray<Real,AMREX_SPACEDIM>* m_dxi_p;
303 const GpuArray<Real,AMREX_SPACEDIM>* m_plo_p;
304 const Dim3* m_lo_p;
305 const Dim3* m_hi_p;
306 int* m_bin_type_array;
307};
308
310{
316
317 template <typename ParticleType>
319 unsigned int operator() (const ParticleType& p) const noexcept
320 {
321 Box tbx;
322 auto iv = getParticleCell(p, plo, dxi, domain);
323 auto tid = getTileIndex(iv, box, true, bin_size, tbx);
324 return static_cast<unsigned int>(tid);
325 }
326};
327
341template <typename P>
346{
347 IntVect iv(
348 AMREX_D_DECL(int(amrex::Math::floor((p.pos(0)-plo[0])*dxi[0])),
349 int(amrex::Math::floor((p.pos(1)-plo[1])*dxi[1])),
350 int(amrex::Math::floor((p.pos(2)-plo[2])*dxi[2]))));
351 return iv;
352}
353
368template <typename P>
373 const Box& domain) noexcept
374{
375 IntVect iv = getParticleCell(p, plo, dxi);
376 iv += domain.smallEnd();
377 return iv;
378}
379
380template <typename PTD>
382IntVect getParticleCell (PTD const& ptd, int i,
385 const Box& domain) noexcept
386{
387 IntVect iv(
388 AMREX_D_DECL(int(amrex::Math::floor((ptd.pos(0, i)-plo[0])*dxi[0])),
389 int(amrex::Math::floor((ptd.pos(1, i)-plo[1])*dxi[1])),
390 int(amrex::Math::floor((ptd.pos(2, i)-plo[2])*dxi[2]))));
391 iv += domain.smallEnd();
392 return iv;
393}
394
396{
397
398 template <typename P>
403 const Box& domain) const noexcept
404 {
405 return getParticleCell(p, plo, dxi, domain);
406 }
407};
408
409template <typename P>
414 const Box& domain) noexcept
415{
416 if (!p.id().is_valid()) { return -1; }
417 IntVect iv = getParticleCell(p, plo, dxi, domain);
418 return mask(iv);
419}
420
421template <typename P>
428 amrex::GpuArray<int,AMREX_SPACEDIM> const& is_per) noexcept
429{
430 bool shifted = false;
431 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim)
432 {
433 if (! is_per[idim]) { continue; }
434 if (p.pos(idim) > rhi[idim]) {
435 while (p.pos(idim) > rhi[idim]) {
436 p.pos(idim) -= static_cast<ParticleReal>(phi[idim] - plo[idim]);
437 }
438 // clamp to avoid precision issues;
439 if (p.pos(idim) < rlo[idim]) {
440 p.pos(idim) = rlo[idim];
441 }
442 shifted = true;
443 }
444 else if (p.pos(idim) < rlo[idim]) {
445 while (p.pos(idim) < rlo[idim]) {
446 p.pos(idim) += static_cast<ParticleReal>(phi[idim] - plo[idim]);
447 }
448 // clamp to avoid precision issues;
449 if (p.pos(idim) > rhi[idim]) {
450 p.pos(idim) = rhi[idim];
451 }
452 shifted = true;
453 }
454 AMREX_ASSERT( (p.pos(idim) >= rlo[idim] ) && ( p.pos(idim) <= rhi[idim] ));
455 }
456
457 return shifted;
458}
459
472template <typename PTile, typename ParFunc>
473int
474partitionParticles (PTile& ptile, ParFunc const& is_left)
475{
476 const int np = ptile.numParticles();
477 if (np == 0) { return 0; }
478
479 auto ptd = ptile.getParticleTileData();
480
481 const int num_left = Reduce::Sum<int>(np,
482 [=] AMREX_GPU_DEVICE (int i) -> int
483 {
484 return int(is_left(ptd, i));
485 });
486
487 // The ptile will be partitioned into left [0, num_left-1] and right [num_left, np-1].
488 //
489 // Note that currently the number of particles in [0, num_left-1] that should belong to the
490 // right partition is equal to the number of particles in [num_left, np-1] that should belong
491 // in the left partition. We will define num_swaps to be this number. This is the minimum
492 // number of swaps that need to be performed to partition the ptile in place for any algorithm.
493 //
494 // From this it is easy to see that
495 // max_num_swaps = min(size([0, num_left-1]), size([num_left, np-1]))
496 // is an upper bound for num_swaps.
497
498 const int max_num_swaps = std::min(num_left, np - num_left);
499 if (max_num_swaps == 0) { return num_left; }
500
501 Gpu::DeviceVector<int> index_left(max_num_swaps);
502 Gpu::DeviceVector<int> index_right(max_num_swaps);
503 int * const p_index_left = index_left.dataPtr();
504 int * const p_index_right = index_right.dataPtr();
505
506 // The num_swaps particles that are in [0, num_left-1] but should be moved to the right
507 // partition are at the same time the first num_swaps particles for which is_left is false
508 // in all the ptile.
509 // Similarly, the num_swaps particles in [num_left, np-1] that should be moved to the left
510 // partition are the last num_swaps particles of the ptile for which is_left is true.
511 //
512 // The PrefixSum is used to find exactly these particles and store their indexes in
513 // index_left and index_right. Since num_swaps is not known, the first max_num_swaps
514 // particles are stored instead. Here, dst = num_left-1-(i-s) is used to effectively reverse
515 // the PrefixSum to store the last particles for which is_left is true.
516 //
517 // This way, all indexes in index_right are in ascending order, and all indexes in
518 // index_left are in descending order.
519
520 Scan::PrefixSum<int>(np,
521 [=] AMREX_GPU_DEVICE (int i) -> int
522 {
523 return int(!is_left(ptd, i));
524 },
525 [=] AMREX_GPU_DEVICE (int i, int const& s)
526 {
527 if (!is_left(ptd, i)) {
528 int dst = s;
529 if (dst < max_num_swaps) {
530 p_index_right[dst] = i;
531 }
532 } else {
533 int dst = num_left-1-(i-s); // avoid integer overflow
534 if (dst < max_num_swaps) {
535 p_index_left[dst] = i;
536 }
537 }
538 },
540
541 // Finally, the particles are swapped. Since max_num_swaps is only an upper bound for num_swaps,
542 // some swaps should not be performed and need to be skipped. This is the case if the index
543 // in index_left[i] is already in the left partition or the index in index_right[i] is already
544 // in the right partition. These two cases coincide for the same i because index_right is in
545 // ascending order and index_left in descending order. This means for both index_left and
546 // index_right the first num_swaps particles need to be swapped, and the particles after that
547 // should be skipped.
548 //
549 // The check right_i < left_i makes sure that the particle going to the right partition is
550 // actually coming from the left partition, which has a lower index than the other particle and
551 // visa-versa.
552 //
553 // Since exactly num_swaps swap operations are performed in the end, which is the smallest
554 // number possible, this algorithm is optimal in the number of swap operations.
555 // This results in good performance in practice if the size of a particle is large enough that
556 // it compensates for the extra kernel launches and evaluations of is_left which this
557 // algorithm needs.
558
559 ParallelForOMP(max_num_swaps,
560 [=] AMREX_GPU_DEVICE (int i)
561 {
562 int left_i = p_index_left[i];
563 int right_i = p_index_right[i];
564 if (right_i < left_i) {
565 swapParticle(ptd, ptd, left_i, right_i);
566 }
567 });
568
569 index_left.free_async();
570 index_right.free_async();
571
572 return num_left;
573}
574
589template <typename PTile, typename ParFunc>
590void
591partitionParticles (PTile& ptile, int num_left, ParFunc const& is_left)
592{
593 const int np = ptile.numParticles();
594 if (np == 0) { return; }
595
596 auto ptd = ptile.getParticleTileData();
597
598 const int max_num_swaps = std::min(num_left, np - num_left);
599 if (max_num_swaps == 0) { return; }
600
601 Gpu::DeviceVector<int> index_left(max_num_swaps);
602 Gpu::DeviceVector<int> index_right(max_num_swaps);
603 int * const p_index_left = index_left.dataPtr();
604 int * const p_index_right = index_right.dataPtr();
605
606 Scan::PrefixSum<int>(np,
607 [=] AMREX_GPU_DEVICE (int i) -> int
608 {
609 return int(!is_left(ptd, i));
610 },
611 [=] AMREX_GPU_DEVICE (int i, int const& s)
612 {
613 if (!is_left(ptd, i)) {
614 int dst = s;
615 if (dst < max_num_swaps) {
616 p_index_right[dst] = i;
617 }
618 } else {
619 int dst = num_left-1-(i-s); // avoid integer overflow
620 if (dst < max_num_swaps) {
621 p_index_left[dst] = i;
622 }
623 }
624 },
626
627 ParallelForOMP(max_num_swaps,
628 [=] AMREX_GPU_DEVICE (int i)
629 {
630 int left_i = p_index_left[i];
631 int right_i = p_index_right[i];
632 if (right_i < left_i) {
633 swapParticle(ptd, ptd, left_i, right_i);
634 }
635 });
636
637 index_left.free_async();
638 index_right.free_async();
639}
640
641template <typename PTile>
642void
644{
645 const int new_size = partitionParticles(ptile,
646 [=] AMREX_GPU_DEVICE (auto& ptd, int i) {
647 return ptd.id(i).is_valid();
648 });
649 ptile.resize(new_size);
650}
651
652template <typename PTile, typename PLocator, typename CellAssignor>
653int
654partitionParticlesByDest (PTile& ptile, const PLocator& ploc, CellAssignor const& assignor,
655 const ParticleBufferMap& pmap,
660 const GpuArray<int ,AMREX_SPACEDIM>& is_per,
661 int lev, int gid, int tid,
662 int lev_min, int lev_max, int nGrow, bool remove_negative)
663{
664 auto getPID = pmap.getPIDFunctor();
665 int pid = ParallelContext::MyProcSub();
666
667 Gpu::DeviceVector<uint8_t> particle_stays(ptile.numParticles());
668 uint8_t * const p_particle_stays = particle_stays.dataPtr();
669 auto ptd = ptile.getParticleTileData();
670
671 // the function for determining if a particle stays on this grid is very slow,
672 // so we cache it in particle_stays to avoid evaluating it multiple times.
673 ParallelForOMP(ptile.numParticles(),
674 [=] AMREX_GPU_DEVICE (int i)
675 {
676 int assigned_grid;
677 int assigned_tile;
678 int assigned_lev;
679
680 if (!ptd.id(i).is_valid())
681 {
682 assigned_grid = -1;
683 assigned_tile = -1;
684 assigned_lev = -1;
685 }
686 else
687 {
688 auto p_prime = ptd.getSuperParticle(i);
689 enforcePeriodic(p_prime, plo, phi, rlo, rhi, is_per);
690 auto tup_prime = ploc(p_prime, lev_min, lev_max, nGrow, assignor);
691 assigned_grid = amrex::get<0>(tup_prime);
692 assigned_tile = amrex::get<1>(tup_prime);
693 assigned_lev = amrex::get<2>(tup_prime);
694 if (assigned_grid >= 0)
695 {
696 AMREX_D_TERM(ptd.pos(0, i) = p_prime.pos(0);,
697 ptd.pos(1, i) = p_prime.pos(1);,
698 ptd.pos(2, i) = p_prime.pos(2););
699 }
700 else if (lev_min > 0)
701 {
702 AMREX_D_TERM(p_prime.pos(0) = ptd.pos(0, i);,
703 p_prime.pos(1) = ptd.pos(1, i);,
704 p_prime.pos(2) = ptd.pos(2, i););
705 auto tup = ploc(p_prime, lev_min, lev_max, nGrow, assignor);
706 assigned_grid = amrex::get<0>(tup);
707 assigned_tile = amrex::get<1>(tup);
708 assigned_lev = amrex::get<2>(tup);
709 }
710 }
711
712 p_particle_stays[i] = uint8_t(
713 ((remove_negative == false) && (!ptd.id(i).is_valid())) ||
714 ((assigned_grid == gid) && (assigned_tile == tid) &&
715 (assigned_lev == lev) && (getPID(lev, gid, tid) == pid)));
716 });
717
718 return partitionParticles(ptile,
719 [=] AMREX_GPU_DEVICE (auto& /* ptd */, int i) -> bool {
720 return p_particle_stays[i];
721 });
722}
723
724template <class PC1, class PC2>
725bool SameIteratorsOK (const PC1& pc1, const PC2& pc2) {
726 if (pc1.numLevels() != pc2.numLevels()) {return false;}
727 if (pc1.do_tiling != pc2.do_tiling) {return false;}
728 if (pc1.tile_size != pc2.tile_size) {return false;}
729 for (int lev = 0; lev < pc1.numLevels(); ++lev) {
730 if (pc1.ParticleBoxArray(lev) != pc2.ParticleBoxArray(lev)) {return false;}
731 if (pc1.ParticleDistributionMap(lev) != pc2.ParticleDistributionMap(lev)) {return false;}
732 }
733 return true;
734}
735
736template <class PC>
738 using Iter = typename PC::ParIterType;
739 for (int lev = 0; lev < pc.numLevels(); ++lev) {
740 for (Iter pti(pc, lev); pti.isValid(); ++pti) {
741 pc.DefineAndReturnParticleTile(lev, pti);
742 }
743 }
744}
745
746IntVect computeRefFac (const ParGDBBase* a_gdb, int src_lev, int lev);
747
748Vector<int> computeNeighborProcs (const ParGDBBase* a_gdb, int ngrow);
749
751namespace particle_detail
752{
753template <typename C>
754void clearEmptyEntries (C& c)
755{
756 for (auto c_it = c.begin(); c_it != c.end(); /* no ++ */)
757 {
758 if (c_it->second.empty()) { c.erase(c_it++); }
759 else { ++c_it; }
760 }
761}
762}
764
765template <class index_type, typename F>
767 index_type nbins, F const& f)
768{
769 BL_PROFILE("PermutationForDeposition()");
770
771#if defined(AMREX_USE_HIP)
772 // MI250X has a small L2 cache and is more tolerant of atomic add contention,
773 // so we use a small block size of 64 and the compressed layout.
774 static constexpr index_type gpu_block_size = 64;
775 static constexpr bool compressed_layout = true;
776#else
777 // A100 has a larger L2 cache and is very sensitive to atomic add contention,
778 // so we use a large bock size of 1024 and not the compressed layout.
779 static constexpr index_type gpu_block_size = 1024;
780 static constexpr bool compressed_layout = false;
781#endif
782
783 static constexpr index_type gpu_block_size_m1 = gpu_block_size - 1;
784 static constexpr index_type llist_guard = std::numeric_limits<index_type>::max();
785
786 // round up to gpu_block_size
787 nbins = (nbins + gpu_block_size_m1) / gpu_block_size * gpu_block_size;
788
789 Gpu::DeviceVector<index_type> llist_start(nbins, llist_guard);
790 Gpu::DeviceVector<index_type> llist_next(nitems);
791 perm.resize(nitems);
792 Gpu::DeviceScalar<index_type> global_idx(0);
793
794 index_type* pllist_start = llist_start.dataPtr();
795 index_type* pllist_next = llist_next.dataPtr();
796 index_type* pperm = perm.dataPtr();
797 index_type* pglobal_idx = global_idx.dataPtr();
798
799 amrex::ParallelFor(nitems, [=] AMREX_GPU_DEVICE (index_type i) noexcept
800 {
801 i = nitems - i - 1;
802 pllist_next[i] = Gpu::Atomic::Exch(pllist_start + f(i), i);
803 });
804
805#if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
806 amrex::launch<gpu_block_size>(nbins / gpu_block_size, Gpu::gpuStream(),
807 [pllist_start,pllist_next,pperm,pglobal_idx] AMREX_GPU_DEVICE () {
808 __shared__ index_type sdata[gpu_block_size];
809 __shared__ index_type global_idx_start;
810 __shared__ index_type idx_start;
811
812 index_type current_idx = 0;
813
814 if constexpr (compressed_layout) {
815 // Compressed layout: subsequent sweeps of up to gpu_block_size contiguous particles
816 // are put right next to each other, while without the compressed layout,
817 // there can be other particle sweeps from different locations between them.
818 current_idx = pllist_start[threadIdx.x + gpu_block_size * blockIdx.x];
819
820 index_type num_particles_thread = 0;
821 while (current_idx != llist_guard) {
822 ++num_particles_thread;
823 current_idx = pllist_next[current_idx];
824 }
825
826 index_type num_particles_block =
827 Gpu::blockReduceSum<gpu_block_size>(num_particles_thread);
828
829 if (threadIdx.x == 0) {
830 global_idx_start = Gpu::Atomic::Add(pglobal_idx, num_particles_block);
831 }
832 }
833
834 current_idx = pllist_start[threadIdx.x + gpu_block_size * blockIdx.x];
835
836 while (true) {
837 sdata[threadIdx.x] = index_type(current_idx != llist_guard);
838 index_type x = 0;
839
840 // simple block wide prefix sum
841 for (index_type i = 1; i<gpu_block_size; i*=2) {
842 __syncthreads();
843 if (threadIdx.x >= i) {
844 x = sdata[threadIdx.x - i];
845 }
846 __syncthreads();
847 if (threadIdx.x >= i) {
848 sdata[threadIdx.x] += x;
849 }
850 }
851 __syncthreads();
852 if (sdata[gpu_block_size_m1] == 0) {
853 break;
854 }
855 if (threadIdx.x == gpu_block_size_m1) {
856 if constexpr (compressed_layout) {
857 idx_start = global_idx_start;
858 global_idx_start += sdata[gpu_block_size_m1];
859 } else {
860 idx_start = Gpu::Atomic::Add(pglobal_idx, sdata[gpu_block_size_m1]);
861 }
862 }
863 __syncthreads();
864 sdata[threadIdx.x] += idx_start;
865 if (current_idx != llist_guard) {
866 pperm[sdata[threadIdx.x] - 1] = current_idx;
867 current_idx = pllist_next[current_idx];
868 }
869 }
870 });
871#else
872 amrex::ignore_unused(pperm, pglobal_idx, compressed_layout);
873 Abort("PermutationForDeposition only implemented for CUDA and HIP");
874#endif
875
876 Gpu::Device::streamSynchronize();
877}
878
879template <class index_type, class PTile>
881 const PTile& ptile, Box bx, Geometry geom, const IntVect idx_type)
882{
883 AMREX_ALWAYS_ASSERT(idx_type.allGE(0) && idx_type.allLE(2));
884
885 const IntVect refine_vect = max(idx_type, IntVect(1)).min(IntVect(2));
886 const IntVect type_vect = idx_type - idx_type / 2 * 2;
887
888 geom.refine(refine_vect);
889
890 Box domain = geom.Domain();
891
892 bx.convert(type_vect);
893 domain.convert(type_vect);
894
895 const RealVect dxi(geom.InvCellSize());
896 const RealVect pos_offset = Real(0.5) * (RealVect(geom.ProbLo()) + RealVect(geom.ProbHi())
897 - RealVect(geom.CellSize()) * RealVect(domain.smallEnd() + domain.bigEnd()));
898
899 const int ref_product = AMREX_D_TERM(refine_vect[0], * refine_vect[1], * refine_vect[2]);
900 const IntVect ref_offset(AMREX_D_DECL(1, refine_vect[0], refine_vect[0] * refine_vect[1]));
901
902 auto ptd = ptile.getConstParticleTileData();
903 PermutationForDeposition<index_type>(perm, nitems, bx.numPts() * ref_product,
904 [=] AMREX_GPU_DEVICE (index_type idx) noexcept
905 {
906 const auto p = ptd[idx];
907
908 IntVect iv = ((p.pos() - pos_offset) * dxi).round();
909
910 IntVect iv_coarse = iv / refine_vect;
911 IntVect iv_remainder = iv - iv_coarse * refine_vect;
912
913 iv_coarse = iv_coarse.max(bx.smallEnd());
914 iv_coarse = iv_coarse.min(bx.bigEnd());
915 return bx.index(iv_coarse) + bx.numPts() * (iv_remainder * ref_offset).sum();
916 });
917}
918
919template <typename P>
920std::string getDefaultCompNameReal (const int i) {
921 int first_r_name = 0;
922 if constexpr (P::is_soa_particle) {
923 if (i < AMREX_SPACEDIM) {
924 constexpr int x_in_ascii = 120;
925 std::string const name{char(x_in_ascii+i)};
926 return name;
927 }
928 first_r_name = AMREX_SPACEDIM;
929 }
930 std::string const name{("real_comp" + std::to_string(i-first_r_name))};
931 return name;
932}
933
934template <typename P>
935std::string getDefaultCompNameInt (const int i) {
936 std::string const name{("int_comp" + std::to_string(i))};
937 return name;
938}
939
940
950template <class PTile, class index_type>
951requires (!PTile::ParticleType::is_rtsoa_particle)
952void
953ReorderParticles (PTile& ptile, const index_type* permutations)
954{
955 const size_t np = ptile.numParticles();
956 const size_t np_total = np + ptile.numNeighborParticles();
957
958#if defined(AMREX_USE_CUDA) && defined(_WIN32)
959 if (!PTile::ParticleType::is_soa_particle)
960#else
961 if constexpr (!PTile::ParticleType::is_soa_particle)
962#endif
963 {
964 static_assert(sizeof(typename PTile::ParticleType)%4 == 0 && sizeof(uint32_t) == 4);
965 using tmp_t = std::conditional_t<sizeof(typename PTile::ParticleType)%8 == 0,
966 uint64_t, uint32_t>;
967 constexpr std::size_t nchunks = sizeof(typename PTile::ParticleType) / sizeof(tmp_t);
969 auto* ptmp = tmp.data();
970 auto* paos = (tmp_t*)(ptile.getParticleTileData().m_aos);
971 for (std::size_t ichunk = 0; ichunk < nchunks; ++ichunk) {
972 // Do not need to reorder neighbor particles
974 {
975 ptmp[i] = paos[permutations[i]*nchunks+ichunk];
976 });
978 {
979 paos[i*nchunks+ichunk] = ptmp[i];
980 });
981 }
982 Gpu::streamSynchronize();
983 } else {
984 typename PTile::SoA::IdCPU tmp_idcpu;
985 if constexpr (PTile::has_polymorphic_allocator) {
986 tmp_idcpu.setArena(ptile.GetStructOfArrays().GetIdCPUData().arena());
987 }
988 tmp_idcpu.resize(np_total);
989 auto src = ptile.GetStructOfArrays().GetIdCPUData().data();
990 uint64_t* dst = tmp_idcpu.data();
991 AMREX_HOST_DEVICE_FOR_1D( np_total, i,
992 {
993 dst[i] = i < np ? src[permutations[i]] : src[i];
994 });
995
996 Gpu::streamSynchronize();
997
998 ptile.GetStructOfArrays().GetIdCPUData().swap(tmp_idcpu);
999 }
1000
1001 { // Create a scope for the temporary vector below
1002 typename PTile::RealVector tmp_real;
1003 if (ptile.NumRealComps() > 0) {
1004 if constexpr (PTile::has_polymorphic_allocator) {
1005 tmp_real.setArena(ptile.GetStructOfArrays().GetRealData(0).arena());
1006 }
1007 tmp_real.resize(np_total);
1008 }
1009 for (int comp = 0; comp < ptile.NumRealComps(); ++comp) {
1010 auto src = ptile.GetStructOfArrays().GetRealData(comp).data();
1011 ParticleReal* dst = tmp_real.data();
1012 AMREX_HOST_DEVICE_FOR_1D( np_total, i,
1013 {
1014 dst[i] = i < np ? src[permutations[i]] : src[i];
1015 });
1016
1017 Gpu::streamSynchronize();
1018
1019 ptile.GetStructOfArrays().GetRealData(comp).swap(tmp_real);
1020 }
1021 }
1022
1023 typename PTile::IntVector tmp_int;
1024 if (ptile.NumIntComps() > 0) {
1025 if constexpr (PTile::has_polymorphic_allocator) {
1026 tmp_int.setArena(ptile.GetStructOfArrays().GetIntData(0).arena());
1027 }
1028 tmp_int.resize(np_total);
1029 }
1030
1031 for (int comp = 0; comp < ptile.NumIntComps(); ++comp) {
1032 auto src = ptile.GetStructOfArrays().GetIntData(comp).data();
1033 int* dst = tmp_int.data();
1034 AMREX_HOST_DEVICE_FOR_1D( np_total , i,
1035 {
1036 dst[i] = i < np ? src[permutations[i]] : src[i];
1037 });
1038
1039 Gpu::streamSynchronize();
1040
1041 ptile.GetStructOfArrays().GetIntData(comp).swap(tmp_int);
1042 }
1043}
1044
1054template <class PTile, class index_type>
1055requires (PTile::ParticleType::is_rtsoa_particle)
1056void
1057ReorderParticles (PTile& ptile, const index_type* permutations)
1058{
1059 const size_t np = ptile.numParticles();
1060 {
1062
1063 auto src = ptile.GetIdCPUData().data();
1064 uint64_t* dst = tmp_idcpu.data();
1066 [=] AMREX_GPU_DEVICE (size_t i) {
1067 dst[i] = src[permutations[i]];
1068 });
1070 [=] AMREX_GPU_DEVICE (size_t i) {
1071 src[i] = dst[i];
1072 });
1073 }
1074 {
1076
1077 for (int comp = 0; comp < ptile.NumRealComps(); ++comp) {
1078 auto src = ptile.GetRealData(comp).data();
1079 auto dst = tmp_real.data();
1081 [=] AMREX_GPU_DEVICE (size_t i) {
1082 dst[i] = src[permutations[i]];
1083 });
1085 [=] AMREX_GPU_DEVICE (size_t i) {
1086 src[i] = dst[i];
1087 });
1088 }
1089 }
1090 {
1092
1093 for (int comp = 0; comp < ptile.NumIntComps(); ++comp) {
1094 auto src = ptile.GetIntData(comp).data();
1095 auto dst = tmp_int.data();
1097 [=] AMREX_GPU_DEVICE (size_t i) {
1098 dst[i] = src[permutations[i]];
1099 });
1101 [=] AMREX_GPU_DEVICE (size_t i) {
1102 src[i] = dst[i];
1103 });
1104 }
1105 }
1106 // stream sync for deallocation of permutations variable in user code
1107 Gpu::streamSynchronize();
1108}
1109
1110#ifdef AMREX_USE_HDF5_ASYNC
1111void async_vol_es_wait_particle();
1112void async_vol_es_wait_close_particle();
1113#endif
1114}
1115
1116#endif // include guard
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#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
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:649
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:124
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:364
__host__ __device__ BoxND & convert(IndexTypeND< dim > typ) noexcept
Convert the BoxND from the current type into the argument type. This may change the BoxND coordinates...
Definition AMReX_Box.H:982
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:216
__host__ __device__ Long index(const IntVectND< dim > &v) const noexcept
Return offset of point from smallend; i.e. index(smallend) -> 0, bigend would return numPts()-1....
Definition AMReX_Box.H:1063
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:112
const Real * CellSize() const noexcept
Returns the cellsize for each coordinate direction.
Definition AMReX_CoordSys.H:71
const Real * InvCellSize() const noexcept
Returns the inverse cellsize for each coordinate direction.
Definition AMReX_CoordSys.H:82
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
void refine(IntVect const &rr)
Refine the Geometry by rr.
Definition AMReX_Geometry.H:429
const Real * ProbHi() const noexcept
Returns the hi end of the problem domain in each dimension.
Definition AMReX_Geometry.H:186
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:216
const Real * ProbLo() const noexcept
Returns the lo end of the problem domain in each dimension.
Definition AMReX_Geometry.H:184
__host__ __device__ constexpr int min() const noexcept
minimum (no absolute values) value
Definition AMReX_IntVect.H:324
__host__ __device__ constexpr bool allGE(const IntVectND< dim > &rhs) const noexcept
Returns true if this is greater than or equal to argument for all components. NOTE: This is NOT a str...
Definition AMReX_IntVect.H:542
__host__ __device__ constexpr int max() const noexcept
maximum (no absolute values) value
Definition AMReX_IntVect.H:313
__host__ __device__ constexpr bool allLE(const IntVectND< dim > &rhs) const noexcept
Returns true if this is less than or equal to argument for all components. NOTE: This is NOT a strict...
Definition AMReX_IntVect.H:492
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
void resize(size_type a_new_size, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_PODVector.H:728
void free_async() noexcept
Definition AMReX_PODVector.H:853
T * dataPtr() noexcept
Definition AMReX_PODVector.H:670
T * data() noexcept
Definition AMReX_PODVector.H:666
Definition AMReX_ParIter.H:118
Definition AMReX_ParticleBufferMap.H:59
GetPID getPIDFunctor() const noexcept
Definition AMReX_ParticleBufferMap.H:195
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
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
void ParallelForOMP(T n, L const &f) noexcept
Performance-portable kernel launch function with optional OpenMP threading.
Definition AMReX_GpuLaunch.H:327
void Sum(T &v, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:221
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
static constexpr struct amrex::Scan::Type::Exclusive exclusive
static constexpr RetSum noRetSum
Definition AMReX_Scan.H:34
Definition AMReX_Amr.cpp:50
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__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
void EnsureThreadSafeTiles(PC &pc)
Definition AMReX_ParticleUtil.H:737
std::string getDefaultCompNameInt(const int i)
Definition AMReX_ParticleUtil.H:935
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
__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
void removeInvalidParticles(PTile &ptile)
Definition AMReX_ParticleUtil.H:643
__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
std::string getDefaultCompNameReal(const int i)
Definition AMReX_ParticleUtil.H:920
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:25
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:45
__host__ __device__ int getParticleGrid(P const &p, amrex::Array4< int > const &mask, amrex::GpuArray< amrex::Real, 3 > const &plo, amrex::GpuArray< amrex::Real, 3 > const &dxi, const Box &domain) noexcept
Definition AMReX_ParticleUtil.H:411
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
void PermutationForDeposition(Gpu::DeviceVector< index_type > &perm, index_type nitems, index_type nbins, F const &f)
Definition AMReX_ParticleUtil.H:766
bool SameIteratorsOK(const PC1 &pc1, const PC2 &pc2)
Definition AMReX_ParticleUtil.H:725
int partitionParticles(PTile &ptile, ParFunc const &is_left)
Reorders the ParticleTile into two partitions left [0, num_left-1] and right [num_left,...
Definition AMReX_ParticleUtil.H:474
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
__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
A multidimensional array accessor.
Definition AMReX_Array4.H:285
Definition AMReX_ParticleUtil.H:265
BinMapper(const int *off_bins_p, const GpuArray< Real, 3 > *dxi_p, const GpuArray< Real, 3 > *plo_p, const Dim3 *lo_p, const Dim3 *hi_p, int *bin_type_array=nullptr)
Definition AMReX_ParticleUtil.H:266
__host__ __device__ unsigned int operator()(const T &ptd, int i) const
Definition AMReX_ParticleUtil.H:277
Definition AMReX_ParticleUtil.H:396
__host__ __device__ IntVect operator()(P const &p, amrex::GpuArray< amrex::Real, 3 > const &plo, amrex::GpuArray< amrex::Real, 3 > const &dxi, const Box &domain) const noexcept
Definition AMReX_ParticleUtil.H:400
Definition AMReX_Dim3.H:13
int x
Definition AMReX_Dim3.H:13
int z
Definition AMReX_Dim3.H:13
int y
Definition AMReX_Dim3.H:13
Definition AMReX_ParticleUtil.H:310
GpuArray< Real, 3 > plo
Definition AMReX_ParticleUtil.H:311
IntVect bin_size
Definition AMReX_ParticleUtil.H:314
__host__ __device__ unsigned int operator()(const ParticleType &p) const noexcept
Definition AMReX_ParticleUtil.H:319
GpuArray< Real, 3 > dxi
Definition AMReX_ParticleUtil.H:312
Box domain
Definition AMReX_ParticleUtil.H:313
Box box
Definition AMReX_ParticleUtil.H:315
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_GpuMemory.H:57
T * dataPtr()
Definition AMReX_GpuMemory.H:91