Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleTransformation.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLETRANSFORMATION_H_
2#define AMREX_PARTICLETRANSFORMATION_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>
11#include <AMReX_ParticleUtil.H>
12
13namespace amrex
14{
15
30template <typename T_ParticleType, int NAR, int NAI>
34 int src_i, int dst_i) noexcept
35{
36 AMREX_ASSERT(dst.m_num_runtime_real == src.m_num_runtime_real);
37 AMREX_ASSERT(dst.m_num_runtime_int == src.m_num_runtime_int );
38
39 if constexpr(!T_ParticleType::is_soa_particle) {
40 dst.m_aos[dst_i] = src.m_aos[src_i];
41 } else {
42 dst.m_idcpu[dst_i] = src.m_idcpu[src_i];
43 }
44 if constexpr(NAR > 0) {
45 for (int j = 0; j < NAR; ++j) {
46 dst.m_rdata[j][dst_i] = src.m_rdata[j][src_i];
47 }
48 }
49 for (int j = 0; j < dst.m_num_runtime_real; ++j) {
50 dst.m_runtime_rdata[j][dst_i] = src.m_runtime_rdata[j][src_i];
51 }
52 if constexpr(NAI > 0) {
53 for (int j = 0; j < NAI; ++j) {
54 dst.m_idata[j][dst_i] = src.m_idata[j][src_i];
55 }
56 }
57 for (int j = 0; j < dst.m_num_runtime_int; ++j) {
58 dst.m_runtime_idata[j][dst_i] = src.m_runtime_idata[j][src_i];
59 }
60}
61
76template <typename T_ParticleType, int NAR, int NAI>
80 int src_i, int dst_i) noexcept
81{
82 AMREX_ASSERT(dst.m_num_runtime_real == src.m_num_runtime_real);
83 AMREX_ASSERT(dst.m_num_runtime_int == src.m_num_runtime_int );
84
85 if constexpr(T_ParticleType::is_soa_particle) {
86 dst.m_idcpu[dst_i] = src.m_idcpu[src_i];
87 } else {
88 dst.m_aos[dst_i] = src.m_aos[src_i];
89 }
90 for (int j = 0; j < NAR; ++j) {
91 dst.m_rdata[j][dst_i] = src.m_rdata[j][src_i];
92 }
93 for (int j = 0; j < dst.m_num_runtime_real; ++j) {
94 dst.m_runtime_rdata[j][dst_i] = src.m_runtime_rdata[j][src_i];
95 }
96 for (int j = 0; j < NAI; ++j) {
97 dst.m_idata[j][dst_i] = src.m_idata[j][src_i];
98 }
99 for (int j = 0; j < dst.m_num_runtime_int; ++j) {
100 dst.m_runtime_idata[j][dst_i] = src.m_runtime_idata[j][src_i];
101 }
102}
103
118template <typename T_ParticleType, int NAR, int NAI>
122 int src_i, int dst_i) noexcept
123{
124 AMREX_ASSERT(dst.m_num_runtime_real == src.m_num_runtime_real);
125 AMREX_ASSERT(dst.m_num_runtime_int == src.m_num_runtime_int );
126
127 if constexpr(T_ParticleType::is_soa_particle) {
128 amrex::Swap(src.m_idcpu[src_i], dst.m_idcpu[dst_i]);
129 } else {
130 amrex::Swap(src.m_aos[src_i], dst.m_aos[dst_i]);
131 }
132 if constexpr (NAR > 0) {
133 for (int j = 0; j < NAR; ++j) {
134 amrex::Swap(dst.m_rdata[j][dst_i], src.m_rdata[j][src_i]);
135 }
136 }
137 for (int j = 0; j < dst.m_num_runtime_real; ++j) {
138 amrex::Swap(dst.m_runtime_rdata[j][dst_i], src.m_runtime_rdata[j][src_i]);
139 }
140 if constexpr (NAI > 0) {
141 for (int j = 0; j < NAI; ++j) {
142 amrex::Swap(dst.m_idata[j][dst_i], src.m_idata[j][src_i]);
143 }
144 }
145 for (int j = 0; j < dst.m_num_runtime_int; ++j) {
146 amrex::Swap(dst.m_runtime_idata[j][dst_i], src.m_runtime_idata[j][src_i]);
147 }
148}
149
159template <class DRType, class DIType, class SRType, class SIType>
164 typename ParticleTileDataRT<DRType, DIType>::size_type dst_i) noexcept
165{
166 AMREX_ASSERT(dst.m_n_real == src.m_n_real);
167 AMREX_ASSERT(dst.m_n_int == src.m_n_int);
168
169 dst.idcpu(dst_i) = src.idcpu(src_i);
170
171 for (int j = 0; j < dst.m_n_real; ++j) {
172 dst.rdata(j)[dst_i] = src.rdata(j)[src_i];
173 }
174
175 for (int j = 0; j < dst.m_n_int; ++j) {
176 dst.idata(j)[dst_i] = src.idata(j)[src_i];
177 }
178}
179
189template <class DRType, class DIType, class SRType, class SIType>
194 typename ParticleTileDataRT<DRType, DIType>::size_type dst_i) noexcept
195{
196 AMREX_ASSERT(dst.m_n_real == src.m_n_real);
197 AMREX_ASSERT(dst.m_n_int == src.m_n_int);
198
199 amrex::Swap(dst.idcpu(dst_i), src.idcpu(src_i));
200
201 for (int j = 0; j < dst.m_n_real; ++j) {
202 amrex::Swap(dst.rdata(j)[dst_i], src.rdata(j)[src_i]);
203 }
204
205 for (int j = 0; j < dst.m_n_int; ++j) {
206 amrex::Swap(dst.idata(j)[dst_i], src.idata(j)[src_i]);
207 }
208}
209
221template <typename DstTile, typename SrcTile>
222void copyParticles (DstTile& dst, const SrcTile& src) noexcept
223{
224 auto np = src.numParticles();
225 copyParticles(dst, src, 0, 0, np);
226}
227
244template <typename DstTile, typename SrcTile, typename Index, typename N,
245 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
246void copyParticles (DstTile& dst, const SrcTile& src,
247 Index src_start, Index dst_start, N n) noexcept
248{
249 const auto src_data = src.getConstParticleTileData();
250 auto dst_data = dst.getParticleTileData();
251
253 {
254 copyParticle(dst_data, src_data, src_start+i, dst_start+i);
255 });
256
258}
259
273template <typename DstTile, typename SrcTile, typename F>
274void transformParticles (DstTile& dst, const SrcTile& src, F&& f) noexcept
275{
276 auto np = src.numParticles();
277 using Index = decltype(np);
278 transformParticles(dst, src, Index{0}, Index{0}, np, std::forward<F>(f));
279}
280
300template <typename DstTile, typename SrcTile, typename Index, typename N, typename F,
301 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
302void transformParticles (DstTile& dst, const SrcTile& src,
303 Index src_start, Index dst_start, N n, F const& f) noexcept
304{
305 const auto src_data = src.getConstParticleTileData();
306 auto dst_data = dst.getParticleTileData();
307
309 {
310 f(dst_data, src_data, src_start+i, dst_start+i);
311 });
312
314}
315
331template <typename DstTile1, typename DstTile2, typename SrcTile, typename F>
332void transformParticles (DstTile1& dst1, DstTile2& dst2, const SrcTile& src, F&& f) noexcept
333{
334 auto np = src.numParticles();
335 using Index = decltype(np);
336 transformParticles(dst1, dst2, src, Index{0}, Index{0}, Index{0}, np, std::forward<F>(f));
337}
338
361template <typename DstTile1, typename DstTile2, typename SrcTile,
362 typename Index, typename N, typename F,
363 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
364void transformParticles (DstTile1& dst1, DstTile2& dst2, const SrcTile& src,
365 Index src_start, Index dst1_start, Index dst2_start, N n, F const& f) noexcept
366{
367 const auto src_data = src.getConstParticleTileData();
368 auto dst1_data = dst1.getParticleTileData();
369 auto dst2_data = dst2.getParticleTileData();
370
372 {
373 f(dst1_data, dst2_data, src_data, src_start+i, dst1_start+i, dst2_start+i);
374 });
375
377}
378
391template <typename DstTile, typename SrcTile, typename Index, typename N,
392 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
393Index filterParticles (DstTile& dst, const SrcTile& src, const Index* mask) noexcept
394{
395 return filterParticles(dst, src, mask, Index{0}, Index{0}, src.numParticles());
396}
397
415template <typename DstTile, typename SrcTile, typename Index, typename N,
416 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
417Index filterParticles (DstTile& dst, const SrcTile& src, const Index* mask,
418 Index src_start, Index dst_start, N n) noexcept
419{
420 Gpu::DeviceVector<Index> offsets(n);
421 Gpu::exclusive_scan(mask, mask+n, offsets.begin());
422
423 Index last_mask=0, last_offset=0;
424 Gpu::copyAsync(Gpu::deviceToHost, mask+n-1, mask + n, &last_mask);
425 Gpu::copyAsync(Gpu::deviceToHost, offsets.data()+n-1, offsets.data()+n, &last_offset);
426
427 auto* p_offsets = offsets.dataPtr();
428
429 const auto src_data = src.getConstParticleTileData();
430 auto dst_data = dst.getParticleTileData();
431
433 {
434 if (mask[i]) { copyParticle(dst_data, src_data, src_start+i, dst_start+p_offsets[i]); }
435 });
436
438 return last_mask + last_offset;
439}
440
453template <typename DstTile, typename SrcTile, typename Pred,
454 std::enable_if_t<!std::is_pointer_v<std::decay_t<Pred>>,int> foo = 0>
455int filterParticles (DstTile& dst, const SrcTile& src, Pred&& p) noexcept
456{
457 return filterParticles(dst, src, std::forward<Pred>(p), 0, 0, src.numParticles());
458}
459
477template <typename DstTile, typename SrcTile, typename Pred, typename Index, typename N,
478 std::enable_if_t<!std::is_pointer_v<std::decay_t<Pred>>,Index> nvccfoo = 0>
479Index filterParticles (DstTile& dst, const SrcTile& src, Pred const& p,
480 Index src_start, Index dst_start, N n) noexcept
481{
483
484 auto* p_mask = mask.dataPtr();
485 const auto src_data = src.getConstParticleTileData();
486
488 [p, p_mask, src_data, src_start] AMREX_GPU_DEVICE (int i, amrex::RandomEngine const& engine) noexcept
489 {
490 amrex::ignore_unused(p, p_mask, src_data, src_start, engine);
491 if constexpr (IsCallable<Pred,decltype(src_data),Index,RandomEngine>::value) {
492 p_mask[i] = p(src_data, src_start+i, engine);
493 } else {
494 p_mask[i] = p(src_data, src_start+i);
495 }
496 });
497 return filterParticles(dst, src, mask.dataPtr(), src_start, dst_start, n);
498}
499
517template <typename DstTile, typename SrcTile, typename Index, typename F,
518 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
519Index filterAndTransformParticles (DstTile& dst, const SrcTile& src, Index* mask, F const& f,
520 Index src_start, Index dst_start) noexcept
521{
522 auto np = src.numParticles();
523 Gpu::DeviceVector<Index> offsets(np);
524 Gpu::exclusive_scan(mask, mask+np, offsets.begin());
525
526 Index last_mask=0, last_offset=0;
527 Gpu::copyAsync(Gpu::deviceToHost, mask+np-1, mask + np, &last_mask);
528 Gpu::copyAsync(Gpu::deviceToHost, offsets.data()+np-1, offsets.data()+np, &last_offset);
529
530 auto const* p_offsets = offsets.dataPtr();
531
532 const auto src_data = src.getConstParticleTileData();
533 auto dst_data = dst.getParticleTileData();
534
536 {
537 if (mask[i]) {
538 f(dst_data, src_data, src_start+i,
539 dst_start+p_offsets[i]);
540 }
541 });
542
544 return last_mask + last_offset;
545}
546
562template <typename DstTile, typename SrcTile, typename Index, typename F,
563 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
564Index filterAndTransformParticles (DstTile& dst, const SrcTile& src, Index* mask, F&& f) noexcept
565{
566 return filterAndTransformParticles(dst, src, mask, std::forward<F>(f), Index{0}, Index{0});
567}
568
584template <typename DstTile, typename SrcTile, typename Pred, typename F,
585 std::enable_if_t<!std::is_pointer_v<std::decay_t<Pred>>,int> foo = 0>
586int filterAndTransformParticles (DstTile& dst, const SrcTile& src, Pred&& p, F&& f) noexcept
587{
588 using Index = decltype(src.numParticles());
589 return filterAndTransformParticles(dst, src, std::forward<Pred>(p), std::forward<F>(f),
590 Index{0}, Index{0});
591}
592
610template <typename DstTile1, typename DstTile2, typename SrcTile, typename Index, typename F,
611 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
612Index filterAndTransformParticles (DstTile1& dst1, DstTile2& dst2,
613 const SrcTile& src, Index* mask, F const& f) noexcept
614{
615 auto np = src.numParticles();
616 Gpu::DeviceVector<Index> offsets(np);
617 Gpu::exclusive_scan(mask, mask+np, offsets.begin());
618
619 Index last_mask=0, last_offset=0;
620 Gpu::copyAsync(Gpu::deviceToHost, mask+np-1, mask + np, &last_mask);
621 Gpu::copyAsync(Gpu::deviceToHost, offsets.data()+np-1, offsets.data()+np, &last_offset);
622
623 auto* p_offsets = offsets.dataPtr();
624
625 const auto src_data = src.getConstParticleTileData();
626 auto dst_data1 = dst1.getParticleTileData();
627 auto dst_data2 = dst2.getParticleTileData();
628
630 {
631 if (mask[i]) { f(dst_data1, dst_data2, src_data, i, p_offsets[i], p_offsets[i]); }
632 });
633
635 return last_mask + last_offset;
636}
637
655template <typename DstTile1, typename DstTile2, typename SrcTile, typename Pred, typename F,
656 std::enable_if_t<!std::is_pointer_v<std::decay_t<Pred>>, int> foo = 0>
657int filterAndTransformParticles (DstTile1& dst1, DstTile2& dst2, const SrcTile& src,
658 Pred const& p, F&& f) noexcept
659{
660 auto np = src.numParticles();
662
663 auto* p_mask = mask.dataPtr();
664 const auto src_data = src.getConstParticleTileData();
665
667 [p, p_mask, src_data] AMREX_GPU_DEVICE (int i, amrex::RandomEngine const& engine) noexcept
668 {
669 amrex::ignore_unused(p, p_mask, src_data, engine);
670 if constexpr (IsCallable<Pred,decltype(src_data),int,RandomEngine>::value) {
671 p_mask[i] = p(src_data, i, engine);
672 } else {
673 p_mask[i] = p(src_data, i);
674 }
675 });
676 return filterAndTransformParticles(dst1, dst2, src, mask.dataPtr(), std::forward<F>(f));
677}
678
679
697template <typename DstTile, typename SrcTile, typename Pred, typename F, typename Index,
698 std::enable_if_t<!std::is_pointer_v<std::decay_t<Pred>>,Index> nvccfoo = 0>
699Index filterAndTransformParticles (DstTile& dst, const SrcTile& src, Pred const& p, F&& f,
700 Index src_start, Index dst_start) noexcept
701{
702 auto np = src.numParticles();
704
705 auto* p_mask = mask.dataPtr();
706 const auto src_data = src.getConstParticleTileData();
707
709 [p, p_mask, src_data, src_start] AMREX_GPU_DEVICE (int i, amrex::RandomEngine const& engine) noexcept
710 {
711 amrex::ignore_unused(p, p_mask, src_data, src_start, engine);
712 if constexpr (IsCallable<Pred,decltype(src_data),Index,RandomEngine>::value) {
713 p_mask[i] = p(src_data, src_start+i, engine);
714 } else {
715 p_mask[i] = p(src_data, src_start+i);
716 }
717 });
718 return filterAndTransformParticles(dst, src, mask.dataPtr(), std::forward<F>(f), src_start, dst_start);
719}
720
721
737template <typename PTile, typename N, typename Index,
738 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
739void gatherParticles (PTile& dst, const PTile& src, N np, const Index* inds)
740{
741 const auto src_data = src.getConstParticleTileData();
742 auto dst_data = dst.getParticleTileData();
743
745 {
746 copyParticle(dst_data, src_data, inds[i], i);
747 });
748
750}
751
767template <typename PTile, typename N, typename Index,
768 std::enable_if_t<std::is_integral_v<Index>, int> foo = 0>
769void scatterParticles (PTile& dst, const PTile& src, N np, const Index* inds)
770{
771 const auto src_data = src.getConstParticleTileData();
772 auto dst_data = dst.getParticleTileData();
773
775 {
776 copyParticle(dst_data, src_data, i, inds[i]);
777 });
778
780}
781
782}
783
784#endif // include guard
#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
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
iterator begin() noexcept
Definition AMReX_PODVector.H:674
T * dataPtr() noexcept
Definition AMReX_PODVector.H:670
T * data() noexcept
Definition AMReX_PODVector.H:666
OutIter exclusive_scan(InIter begin, InIter end, OutIter result)
Definition AMReX_Scan.H:1440
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
A host-to-device copy routine. Note this is just a wrapper around memcpy, so it assumes contiguous st...
Definition AMReX_GpuContainers.H:228
static constexpr DeviceToHost deviceToHost
Definition AMReX_GpuContainers.H:106
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
Definition AMReX_Amr.cpp:49
__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__ void Swap(T &t1, T &t2) noexcept
Definition AMReX_Algorithm.H:93
__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
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:519
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:393
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:739
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:274
void scatterParticles(PTile &dst, const PTile &src, N np, const Index *inds)
Scatter particles copies particles from contiguous order into an arbitrary order. Specifically,...
Definition AMReX_ParticleTransformation.H:769
AMREX_ATTRIBUTE_FLATTEN_FOR void ParallelForRNG(T n, L const &f) noexcept
Definition AMReX_GpuLaunchFunctsC.H:1231
Definition AMReX_ParticleTile.H:516
Test if a given type T is callable with arguments of type Args...
Definition AMReX_TypeTraits.H:213
Definition AMReX_ParticleTileRT.H:71
Long size_type
Definition AMReX_ParticleTileRT.H:73
Definition AMReX_ParticleTile.H:33
Definition AMReX_RandomEngine.H:72