1#ifndef AMREX_PARTICLEUTIL_H_
2#define AMREX_PARTICLEUTIL_H_
3#include <AMReX_Config.H>
32template <class Iterator, std::enable_if_t<IsParticleIterator<Iterator>::value,
int> foo = 0>
50template <class Iterator, std::enable_if_t<IsParticleIterator<Iterator>::value,
int> foo = 0>
54 const auto& tile = pti.GetParticleTile();
55 const auto np = tile.numParticles();
56 const auto& ptd = tile.getConstParticleTileData();
57 const auto& geom = pti.Geom(pti.GetLevel());
59 const auto& domain = geom.Domain();
60 const auto& plo = geom.ProbLoArray();
61 const auto& dxi = geom.InvCellSizeArray();
63 Box box = pti.tilebox();
68 using ReduceTuple =
typename decltype(reduce_data)::Type;
70 reduce_op.
eval(np, reduce_data,
74 if (!p.id().is_valid()) {
return false; }
75 using AssignorType =
typename Iterator::CellAssignor;
76 AssignorType assignor;
77 IntVect iv = assignor(p, plo, dxi, domain);
80 int hv = amrex::get<0>(reduce_data.
value(reduce_op));
96template <class PC, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
115template <class PC, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
136template <class PC, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
160template <class PC, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
166 using ParIter =
typename PC::ParConstIterType;
168 for (
int lev = lev_min; lev <= lev_max; ++lev)
171#pragma omp parallel if (Gpu::notInLaunchRegion() && !system::regtest_reduction) reduction(+:num_wrong)
187 if (a_do_tiling ==
false) {
194 auto tiling_1d = [](
int i,
int lo,
int hi,
int tilesize,
195 int& ntile,
int& tileidx,
int& tlo,
int& thi) {
196 int ncells = hi-lo+1;
198 int ts_right = ncells/ntile;
199 int ts_left = ts_right+1;
200 int nleft = ncells - ntile*ts_right;
202 int nbndry = nleft*ts_left;
204 tileidx = ii / ts_left;
205 tlo = lo + tileidx * ts_left;
206 thi = tlo + ts_left - 1;
208 tileidx = nleft + (ii-nbndry) / ts_right;
209 tlo = lo + tileidx * ts_right + nleft;
210 thi = tlo + ts_right - 1;
215 IntVect ntiles, ivIndex, tilelo, tilehi;
221 AMREX_D_TERM(tiling_1d(iv0, sml[0], big[0], a_tile_size[0], ntiles[0], ivIndex[0], tilelo[0], tilehi[0]);,
222 tiling_1d(iv1, sml[1], big[1], a_tile_size[1], ntiles[1], ivIndex[1], tilelo[1], tilehi[1]);,
223 tiling_1d(iv2, sml[2], big[2], a_tile_size[2], ntiles[2], ivIndex[2], tilelo[2], tilehi[2]););
225 tbx =
Box(tilelo, tilehi);
227 return AMREX_D_TERM(ivIndex[0], + ntiles[0]*ivIndex[1], + ntiles[0]*ntiles[1]*ivIndex[2]);
234 if (a_do_tiling ==
false) {
240 auto tiling_1d = [](
int lo,
int hi,
int tilesize,
int& ntile) {
241 int ncells = hi-lo+1;
249 AMREX_D_TERM(tiling_1d(sml[0], big[0], a_tile_size[0], ntiles[0]);,
250 tiling_1d(sml[1], big[1], a_tile_size[1], ntiles[1]);,
251 tiling_1d(sml[2], big[2], a_tile_size[2], ntiles[2]););
253 return AMREX_D_TERM(ntiles[0], *=ntiles[1], *=ntiles[2]);
264 int* bin_type_array=
nullptr)
268 template <
typename T>
282 static_cast<int>(amrex::Math::floor((p.pos(2)-
m_plo_p[type][2])*
m_dxi_p[type][2])) -
m_lo_p[type].
z));
283 auto iv3 = iv.dim3();
290 return static_cast<unsigned int>( (uix * ny + uiy) * nz + uiz +
offset );
310 template <
typename ParticleType>
312 unsigned int operator() (
const ParticleType& p)
const noexcept
317 return static_cast<unsigned int>(tid);
341 AMREX_D_DECL(
int(amrex::Math::floor((p.pos(0)-plo[0])*dxi[0])),
342 int(amrex::Math::floor((p.pos(1)-plo[1])*dxi[1])),
343 int(amrex::Math::floor((p.pos(2)-plo[2])*dxi[2]))));
366 const Box& domain)
noexcept
369 iv += domain.smallEnd();
373template <
typename PTD>
378 const Box& domain)
noexcept
380 if constexpr (PTD::ParticleType::is_soa_particle)
383 AMREX_D_DECL(
int(amrex::Math::floor((ptd.m_rdata[0][i]-plo[0])*dxi[0])),
384 int(amrex::Math::floor((ptd.m_rdata[1][i]-plo[1])*dxi[1])),
385 int(amrex::Math::floor((ptd.m_rdata[2][i]-plo[2])*dxi[2]))));
386 iv += domain.smallEnd();
396 template <
typename P>
401 const Box& domain)
const noexcept
412 const Box& domain)
noexcept
414 if (!p.id().is_valid()) {
return -1; }
428 bool shifted =
false;
429 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim)
431 if (! is_per[idim]) {
continue; }
432 if (p.pos(idim) > rhi[idim]) {
433 while (p.pos(idim) > rhi[idim]) {
434 p.pos(idim) -=
static_cast<ParticleReal
>(phi[idim] - plo[idim]);
437 if (p.pos(idim) < rlo[idim]) {
438 p.pos(idim) = rlo[idim];
442 else if (p.pos(idim) < rlo[idim]) {
443 while (p.pos(idim) < rlo[idim]) {
444 p.pos(idim) +=
static_cast<ParticleReal
>(phi[idim] - plo[idim]);
447 if (p.pos(idim) > rhi[idim]) {
448 p.pos(idim) = rhi[idim];
452 AMREX_ASSERT( (p.pos(idim) >= rlo[idim] ) && ( p.pos(idim) <= rhi[idim] ));
470template <
typename PTile,
typename ParFunc>
474 const int np = ptile.numParticles();
475 if (np == 0) {
return 0; }
477 auto ptd = ptile.getParticleTileData();
479 const int num_left = Reduce::Sum<int>(np,
482 return int(is_left(ptd, i));
496 const int max_num_swaps = std::min(num_left, np - num_left);
497 if (max_num_swaps == 0) {
return num_left; }
501 int *
const p_index_left = index_left.
dataPtr();
502 int *
const p_index_right = index_right.
dataPtr();
518 Scan::PrefixSum<int>(np,
521 return int(!is_left(ptd, i));
525 if (!is_left(ptd, i)) {
527 if (dst < max_num_swaps) {
528 p_index_right[dst] = i;
531 int dst = num_left-1-(i-s);
532 if (dst < max_num_swaps) {
533 p_index_left[dst] = i;
560 int left_i = p_index_left[i];
561 int right_i = p_index_right[i];
562 if (right_i < left_i) {
587template <
typename PTile,
typename ParFunc>
591 const int np = ptile.numParticles();
592 if (np == 0) {
return; }
594 auto ptd = ptile.getParticleTileData();
596 const int max_num_swaps = std::min(num_left, np - num_left);
597 if (max_num_swaps == 0) {
return; }
601 int *
const p_index_left = index_left.
dataPtr();
602 int *
const p_index_right = index_right.
dataPtr();
604 Scan::PrefixSum<int>(np,
607 return int(!is_left(ptd, i));
611 if (!is_left(ptd, i)) {
613 if (dst < max_num_swaps) {
614 p_index_right[dst] = i;
617 int dst = num_left-1-(i-s);
618 if (dst < max_num_swaps) {
619 p_index_left[dst] = i;
628 int left_i = p_index_left[i];
629 int right_i = p_index_right[i];
630 if (right_i < left_i) {
639template <
typename PTile>
645 return ptd.id(i).is_valid();
647 ptile.resize(new_size);
650#if defined(AMREX_USE_GPU)
652template <
typename PTile,
typename PLocator,
typename CellAssignor>
661 int lev,
int gid,
int ,
662 int lev_min,
int lev_max,
int nGrow,
bool remove_negative)
668 uint8_t *
const p_particle_stays = particle_stays.
dataPtr();
669 auto ptd = ptile.getParticleTileData();
679 if (!ptd.id(i).is_valid())
686 auto p_prime = ptd.getSuperParticle(i);
687 enforcePeriodic(p_prime, plo, phi, rlo, rhi, is_per);
688 auto tup_prime = ploc(p_prime, lev_min, lev_max, nGrow, assignor);
689 assigned_grid = amrex::get<0>(tup_prime);
690 assigned_lev = amrex::get<1>(tup_prime);
691 if (assigned_grid >= 0)
693 AMREX_D_TERM(ptd.pos(0, i) = p_prime.pos(0);,
694 ptd.pos(1, i) = p_prime.pos(1);,
695 ptd.pos(2, i) = p_prime.pos(2););
697 else if (lev_min > 0)
700 p_prime.pos(1) = ptd.pos(1, i);,
701 p_prime.pos(2) = ptd.pos(2, i););
702 auto tup = ploc(p_prime, lev_min, lev_max, nGrow, assignor);
703 assigned_grid = amrex::get<0>(tup);
704 assigned_lev = amrex::get<1>(tup);
708 p_particle_stays[i] = uint8_t(
709 ((remove_negative ==
false) && (!ptd.id(i).is_valid())) ||
710 ((assigned_grid == gid) && (assigned_lev == lev) && (getPID(lev, gid) == pid)));
715 return p_particle_stays[i];
721template <
class PC1,
class PC2>
723 if (pc1.numLevels() != pc2.numLevels()) {
return false;}
724 if (pc1.do_tiling != pc2.do_tiling) {
return false;}
725 if (pc1.tile_size != pc2.tile_size) {
return false;}
726 for (
int lev = 0; lev < pc1.numLevels(); ++lev) {
727 if (pc1.ParticleBoxArray(lev) != pc2.ParticleBoxArray(lev)) {
return false;}
728 if (pc1.ParticleDistributionMap(lev) != pc2.ParticleDistributionMap(lev)) {
return false;}
735 using Iter =
typename PC::ParIterType;
736 for (
int lev = 0; lev < pc.numLevels(); ++lev) {
737 for (Iter pti(pc, lev); pti.isValid(); ++pti) {
738 pc.DefineAndReturnParticleTile(lev, pti);
743IntVect computeRefFac (
const ParGDBBase* a_gdb,
int src_lev,
int lev);
745Vector<int> computeNeighborProcs (
const ParGDBBase* a_gdb,
int ngrow);
747namespace particle_detail
752 for (
auto c_it = c.begin(); c_it != c.end(); )
754 if (c_it->second.empty()) { c.erase(c_it++); }
760template <
class index_type,
typename F>
762 index_type nbins, F
const& f)
766#if defined(AMREX_USE_HIP)
769 static constexpr index_type gpu_block_size = 64;
770 static constexpr bool compressed_layout =
true;
774 static constexpr index_type gpu_block_size = 1024;
775 static constexpr bool compressed_layout =
false;
778 static constexpr index_type gpu_block_size_m1 = gpu_block_size - 1;
779 static constexpr index_type llist_guard = std::numeric_limits<index_type>::max();
782 nbins = (nbins + gpu_block_size_m1) / gpu_block_size * gpu_block_size;
789 index_type* pllist_start = llist_start.
dataPtr();
790 index_type* pllist_next = llist_next.
dataPtr();
791 index_type* pperm = perm.
dataPtr();
792 index_type* pglobal_idx = global_idx.
dataPtr();
797 pllist_next[i] = Gpu::Atomic::Exch(pllist_start + f(i), i);
800#if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
801 amrex::launch<gpu_block_size>(nbins / gpu_block_size, Gpu::gpuStream(),
803 __shared__ index_type sdata[gpu_block_size];
804 __shared__ index_type global_idx_start;
805 __shared__ index_type idx_start;
807 index_type current_idx = 0;
809 if constexpr (compressed_layout) {
813 current_idx = pllist_start[threadIdx.x + gpu_block_size * blockIdx.x];
815 index_type num_particles_thread = 0;
816 while (current_idx != llist_guard) {
817 ++num_particles_thread;
818 current_idx = pllist_next[current_idx];
821 index_type num_particles_block =
822 Gpu::blockReduceSum<gpu_block_size>(num_particles_thread);
824 if (threadIdx.x == 0) {
825 global_idx_start = Gpu::Atomic::Add(pglobal_idx, num_particles_block);
829 current_idx = pllist_start[threadIdx.x + gpu_block_size * blockIdx.x];
832 sdata[threadIdx.x] = index_type(current_idx != llist_guard);
836 for (index_type i = 1; i<gpu_block_size; i*=2) {
838 if (threadIdx.x >= i) {
839 x = sdata[threadIdx.x - i];
842 if (threadIdx.x >= i) {
843 sdata[threadIdx.x] +=
x;
847 if (sdata[gpu_block_size_m1] == 0) {
850 if (threadIdx.x == gpu_block_size_m1) {
851 if constexpr (compressed_layout) {
852 idx_start = global_idx_start;
853 global_idx_start += sdata[gpu_block_size_m1];
855 idx_start = Gpu::Atomic::Add(pglobal_idx, sdata[gpu_block_size_m1]);
859 sdata[threadIdx.x] += idx_start;
860 if (current_idx != llist_guard) {
861 pperm[sdata[threadIdx.x] - 1] = current_idx;
862 current_idx = pllist_next[current_idx];
868 Abort(
"PermutationForDeposition only implemented for CUDA and HIP");
871 Gpu::Device::streamSynchronize();
874template <
class index_type,
class PTile>
881 const IntVect type_vect = idx_type - idx_type / 2 * 2;
894 const int ref_product =
AMREX_D_TERM(refine_vect[0], * refine_vect[1], * refine_vect[2]);
895 const IntVect ref_offset(
AMREX_D_DECL(1, refine_vect[0], refine_vect[0] * refine_vect[1]));
897 auto ptd = ptile.getConstParticleTileData();
898 PermutationForDeposition<index_type>(perm, nitems, bx.
numPts() * ref_product,
901 const auto p = ptd[idx];
903 IntVect iv = ((p.pos() - pos_offset) * dxi).round();
905 IntVect iv_coarse = iv / refine_vect;
906 IntVect iv_remainder = iv - iv_coarse * refine_vect;
910 return bx.
index(iv_coarse) + bx.
numPts() * (iv_remainder * ref_offset).sum();
916 int first_r_name = 0;
917 if constexpr (P::is_soa_particle) {
918 if (i < AMREX_SPACEDIM) {
919 constexpr int x_in_ascii = 120;
920 std::string
const name{char(x_in_ascii+i)};
923 first_r_name = AMREX_SPACEDIM;
925 std::string
const name{(
"real_comp" + std::to_string(i-first_r_name))};
931 std::string
const name{(
"int_comp" + std::to_string(i))};
935#ifdef AMREX_USE_HDF5_ASYNC
936void async_vol_es_wait_particle();
937void async_vol_es_wait_close_particle();
#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_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:1089
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:630
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Get the bigend.
Definition AMReX_Box.H:119
__host__ __device__ Long numPts() const noexcept
Returns the number of points contained in the BoxND.
Definition AMReX_Box.H:349
__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:930
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition AMReX_Box.H:207
__host__ __device__ Long index(const IntVectND< dim > &v) const noexcept
Returns offset of point from smallend; i.e. index(smallend) -> 0, bigend would return numPts()-1....
Definition AMReX_Box.H:1011
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition AMReX_Box.H:108
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:73
void refine(IntVect const &rr)
Definition AMReX_Geometry.H:411
const Real * ProbHi() const noexcept
Returns the hi end of the problem domain in each dimension.
Definition AMReX_Geometry.H:180
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:210
const Real * ProbLo() const noexcept
Returns the lo end of the problem domain in each dimension.
Definition AMReX_Geometry.H:178
__host__ __device__ 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:448
__host__ __device__ int min() const noexcept
minimum (no absolute values) value
Definition AMReX_IntVect.H:230
__host__ __device__ int max() const noexcept
maximum (no absolute values) value
Definition AMReX_IntVect.H:219
__host__ __device__ 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:398
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
Definition AMReX_PODVector.H:297
void resize(size_type a_new_size, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_PODVector.H:717
void free_async() noexcept
Definition AMReX_PODVector.H:842
T * dataPtr() noexcept
Definition AMReX_PODVector.H:659
Definition AMReX_ParIter.H:115
Definition AMReX_ParticleBufferMap.H:53
GetPID getPIDFunctor() const noexcept
Definition AMReX_ParticleBufferMap.H:156
Definition AMReX_Reduce.H:249
Type value()
Definition AMReX_Reduce.H:281
Definition AMReX_Reduce.H:364
std::enable_if_t< IsFabArray< MF >::value > eval(MF const &mf, IntVect const &nghost, D &reduce_data, F &&f)
Definition AMReX_Reduce.H:433
void Sum(T &v, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:204
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:33
void clearEmptyEntries(C &c)
Definition AMReX_ParticleUtil.H:750
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:138
__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:119
__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:184
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:191
void EnsureThreadSafeTiles(PC &pc)
Definition AMReX_ParticleUtil.H:734
std::string getDefaultCompNameInt(const int i)
Definition AMReX_ParticleUtil.H:930
__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:421
void removeInvalidParticles(PTile &ptile)
Definition AMReX_ParticleUtil.H:641
__host__ __device__ int numTilesInBox(const Box &box, const bool a_do_tiling, const IntVect &a_tile_size)
Definition AMReX_ParticleUtil.H:232
BoxND< 3 > Box
Definition AMReX_BaseFwd.H:27
std::string getDefaultCompNameReal(const int i)
Definition AMReX_ParticleUtil.H:915
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:21
IntVectND< 3 > IntVect
Definition AMReX_BaseFwd.H:30
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:35
__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:409
void PermutationForDeposition(Gpu::DeviceVector< index_type > &perm, index_type nitems, index_type nbins, F const &f)
Definition AMReX_ParticleUtil.H:761
bool SameIteratorsOK(const PC1 &pc1, const PC2 &pc2)
Definition AMReX_ParticleUtil.H:722
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:472
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:34
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, int lev_min, int lev_max, int nGrow, bool remove_negative)
Definition AMReX_ParticleUtil.H:654
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
__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:336
Definition AMReX_Array4.H:61
Definition AMReX_ParticleUtil.H:258
const Dim3 * m_hi_p
Definition AMReX_ParticleUtil.H:298
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:259
const Dim3 * m_lo_p
Definition AMReX_ParticleUtil.H:297
const GpuArray< Real, 3 > * m_plo_p
Definition AMReX_ParticleUtil.H:296
int * m_bin_type_array
Definition AMReX_ParticleUtil.H:299
const int * m_off_bins_p
Definition AMReX_ParticleUtil.H:294
__host__ __device__ unsigned int operator()(const T &ptd, int i) const
Definition AMReX_ParticleUtil.H:270
const GpuArray< Real, 3 > * m_dxi_p
Definition AMReX_ParticleUtil.H:295
Definition AMReX_ParticleUtil.H:394
__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:398
Definition AMReX_Dim3.H:12
int x
Definition AMReX_Dim3.H:12
int z
Definition AMReX_Dim3.H:12
int y
Definition AMReX_Dim3.H:12
Definition AMReX_ParticleUtil.H:303
GpuArray< Real, 3 > plo
Definition AMReX_ParticleUtil.H:304
IntVect bin_size
Definition AMReX_ParticleUtil.H:307
__host__ __device__ unsigned int operator()(const ParticleType &p) const noexcept
Definition AMReX_ParticleUtil.H:312
GpuArray< Real, 3 > dxi
Definition AMReX_ParticleUtil.H:305
Box domain
Definition AMReX_ParticleUtil.H:306
Box box
Definition AMReX_ParticleUtil.H:308
Definition AMReX_Array.H:34
Definition AMReX_GpuMemory.H:56
T * dataPtr()
Definition AMReX_GpuMemory.H:90