1#ifndef AMREX_PARTICLE_LOCATOR_H_
2#define AMREX_PARTICLE_LOCATOR_H_
3#include <AMReX_Config.H>
15template <
class BinIteratorFactory>
37 bool a_do_tiling,
const IntVect& a_tile_size)
41 m_plo(a_geom.ProbLoArray()),
m_dxi(a_geom.InvCellSizeArray()),
62 template <
typename P,
typename Assignor = DefaultAssignor>
64 std::pair<int, int>
operator() (
const P& p,
int nGrow=0, Assignor
const& assignor = Assignor{})
const noexcept
76 const auto lo = iv.
dim3();
86 for (
int ii = ix_lo; ii <= ix_hi; ++ii) {
87 for (
int jj = iy_lo; jj <= iy_hi; ++jj) {
88 for (
int kk = iz_lo; kk <= iz_hi; ++kk) {
90 for (
const auto& nbor :
m_bif.getBinIterator(index)) {
93 return {nbor.first,
getTile(iv, bx)};
100 loc_gid = nbor.first;
104 for (
int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
106 gdbx.
grow(dir, nGrow);
108 loc_gid = nbor.first;
117 return {loc_gid, loc_tid};
131 bool a_do_tiling =
false,
139 int num_boxes =
static_cast<int>(ba.
size());
141 for (
int i = 0; i < num_boxes; ++i) {
m_host_boxes.push_back(ba[i]); }
146 if (num_boxes == 0) {
161 using ReduceTuple =
typename decltype(reduce_data)::Type;
164 reduce_op.
eval(num_boxes, reduce_data,
167 const Box& box = boxes_ptr[i];
176 ReduceTuple hv = reduce_data.
value(reduce_op);
182 amrex::get< AMREX_SPACEDIM+1>(hv),
183 amrex::get< AMREX_SPACEDIM+2>(hv)));
185 amrex::get<2*AMREX_SPACEDIM+1>(hv),
186 amrex::get<2*AMREX_SPACEDIM+2>(hv)));
193 m_bins.build(num_boxes, boxes_ptr, bins_box,
196 return (box.
smallEnd() - bins_lo) / bin_size;
215 bool a_do_tiling =
false,
246template <
class BinIteratorFactory>
256 template <
typename P,
typename Assignor = DefaultAssignor>
259 Assignor
const& assignor = {})
const noexcept
261 lev_min = (lev_min == -1) ? 0 : lev_min;
262 lev_max = (lev_max == -1) ?
m_size - 1 : lev_max;
264 for (
int lev = lev_max; lev >= lev_min; --lev)
266 auto grid_tile =
m_funcs[lev](p, 0, assignor);
267 if (grid_tile.first >= 0) {
return makeTuple(grid_tile.first, grid_tile.second, lev); }
272 auto grid_tile =
m_funcs[lev](p, nGrow, assignor);
273 if (grid_tile.first >= 0) {
return makeTuple(grid_tile.first, grid_tile.second, lev); }
289 bool m_defined =
false;
297 bool a_do_tiling =
false,
300 build(a_ba, a_geom, a_do_tiling, a_tile_size);
304 bool a_do_tiling =
false,
307 build(a_gdb, a_do_tiling, a_tile_size);
312 bool a_do_tiling =
false,
316 int num_levels =
static_cast<int>(a_ba.
size());
317 m_locators.resize(num_levels);
318 m_grid_assignors.
resize(num_levels);
321 for (
int lev = 0; lev < num_levels; ++lev)
323 m_locators[lev].build(a_ba[lev], a_geom[lev], a_do_tiling, a_tile_size);
324 h_grid_assignors[lev] = m_locators[lev].getGridAssignor();
330 for (
int lev = 0; lev < num_levels; ++lev)
332 m_locators[lev].build(a_ba[lev], a_geom[lev], a_do_tiling, a_tile_size);
333 m_grid_assignors[lev] = m_locators[lev].getGridAssignor();
339 bool a_do_tiling =
false,
345 for (
int lev = 0; lev < num_levels; ++lev)
348 geom.push_back(a_gdb->
Geom(lev));
350 build(ba, geom, a_do_tiling, a_tile_size);
354 bool a_do_tiling =
false,
357 if ( !m_defined || (m_locators.empty()) ||
358 (m_locators.
size() != a_ba.
size()) ) {
return false; }
359 bool all_valid =
true;
360 int num_levels = m_locators.
size();
361 for (
int lev = 0; lev < num_levels; ++lev) {
362 all_valid = all_valid && m_locators[lev].isValid(a_ba[lev], a_do_tiling, a_tile_size);
368 bool a_do_tiling =
false,
373 for (
int lev = 0; lev < num_levels; ++lev) {
376 return this->
isValid(ba, a_do_tiling, a_tile_size);
384 for (
int lev = 0; lev < num_levels; ++lev)
386 m_locators[lev].setGeometry(a_gdb->
Geom(lev));
387 h_grid_assignors[lev] = m_locators[lev].getGridAssignor();
393 for (
int lev = 0; lev < num_levels; ++lev)
395 m_locators[lev].setGeometry(a_gdb->
Geom(lev));
396 m_grid_assignors[lev] = m_locators[lev].getGridAssignor();
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#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
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Definition AMReX_ParticleLocator.H:282
void build(const ParGDBBase *a_gdb, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000))
Definition AMReX_ParticleLocator.H:338
AmrParticleLocator(const ParGDBBase *a_gdb, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000))
Definition AMReX_ParticleLocator.H:303
AmrAssignGrid< BinIteratorFactory > getGridAssignor() const noexcept
Definition AMReX_ParticleLocator.H:401
typename Bins::BinIteratorFactory BinIteratorFactory
Definition AMReX_ParticleLocator.H:284
bool isValid(const ParGDBBase *a_gdb, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000)) const
Definition AMReX_ParticleLocator.H:367
void build(const Vector< BoxArray > &a_ba, const Vector< Geometry > &a_geom, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000))
Definition AMReX_ParticleLocator.H:310
AmrParticleLocator(const Vector< BoxArray > &a_ba, const Vector< Geometry > &a_geom, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000))
Definition AMReX_ParticleLocator.H:295
bool isValid(const Vector< BoxArray > &a_ba, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000)) const
Definition AMReX_ParticleLocator.H:353
void setGeometry(const ParGDBBase *a_gdb)
Definition AMReX_ParticleLocator.H:379
AmrParticleLocator()=default
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:841
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:615
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:123
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:212
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:111
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
GPU-compatible tuple.
Definition AMReX_Tuple.H:98
__host__ static __device__ constexpr IntVectND< dim > TheUnitVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:781
__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
__host__ __device__ constexpr Dim3 dim3() const noexcept
Definition AMReX_IntVect.H:265
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
T * dataPtr() noexcept
Definition AMReX_PODVector.H:670
T * data() noexcept
Definition AMReX_PODVector.H:666
Definition AMReX_ParGDB.H:13
virtual const BoxArray & ParticleBoxArray(int level) const =0
virtual const Geometry & Geom(int level) const =0
virtual int finestLevel() const =0
Definition AMReX_ParticleLocator.H:123
ParticleLocator()=default
IntVect m_bin_size
Definition AMReX_ParticleLocator.H:237
IntVect m_bins_hi
Definition AMReX_ParticleLocator.H:236
bool isValid(const BoxArray &ba, bool a_do_tiling=false, const IntVect &a_tile_size=IntVect(1024000, 1024000, 1024000)) const noexcept
Definition AMReX_ParticleLocator.H:214
BoxArray m_ba
Definition AMReX_ParticleLocator.H:231
bool m_do_tiling
Definition AMReX_ParticleLocator.H:229
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
Geometry m_geom
Definition AMReX_ParticleLocator.H:232
Gpu::HostVector< Box > m_host_boxes
Definition AMReX_ParticleLocator.H:242
IntVect m_bins_lo
Definition AMReX_ParticleLocator.H:235
IntVect m_tile_size
Definition AMReX_ParticleLocator.H:233
Bins m_bins
Definition AMReX_ParticleLocator.H:240
IntVect m_num_bins
Definition AMReX_ParticleLocator.H:238
AssignGrid< BinIteratorFactory > getGridAssignor() const noexcept
Definition AMReX_ParticleLocator.H:206
typename Bins::BinIteratorFactory BinIteratorFactory
Definition AMReX_ParticleLocator.H:126
Gpu::DeviceVector< Box > m_device_boxes
Definition AMReX_ParticleLocator.H:243
void setGeometry(const Geometry &a_geom) noexcept
Definition AMReX_ParticleLocator.H:200
bool m_defined
Definition AMReX_ParticleLocator.H:228
Definition AMReX_Reduce.H:453
Type value()
Definition AMReX_Reduce.H:488
Definition AMReX_Reduce.H:612
std::enable_if_t< IsFabArray< MF >::value > eval(MF const &mf, IntVect const &nghost, D &reduce_data, F &&f)
Definition AMReX_Reduce.H:746
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Long size() const noexcept
Definition AMReX_Vector.H:53
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 HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:105
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:421
Definition AMReX_Amr.cpp:49
__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:185
__host__ __device__ constexpr GpuTuple< detail::tuple_decay_t< Ts >... > makeTuple(Ts &&... args)
Definition AMReX_Tuple.H:263
__host__ __device__ constexpr const T & min(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:24
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:44
Definition AMReX_ParticleLocator.H:248
std::size_t m_size
Definition AMReX_ParticleLocator.H:250
AmrAssignGrid(const AssignGrid< BinIteratorFactory > *a_funcs, std::size_t a_size)
Definition AMReX_ParticleLocator.H:252
__host__ __device__ GpuTuple< int, int, int > operator()(const P &p, int lev_min=-1, int lev_max=-1, int nGrow=0, Assignor const &assignor={}) const noexcept
Definition AMReX_ParticleLocator.H:258
const AssignGrid< BinIteratorFactory > * m_funcs
Definition AMReX_ParticleLocator.H:249
Definition AMReX_ParticleLocator.H:17
Dim3 m_hi
Definition AMReX_ParticleLocator.H:21
__host__ __device__ int getTile(const IntVect &iv, const Box &bx) const noexcept
Definition AMReX_ParticleLocator.H:55
__host__ __device__ std::pair< int, int > operator()(const P &p, int nGrow=0, Assignor const &assignor=Assignor{}) const noexcept
Definition AMReX_ParticleLocator.H:64
Dim3 m_bin_size
Definition AMReX_ParticleLocator.H:22
bool m_do_tiling
Definition AMReX_ParticleLocator.H:29
Dim3 m_tile_size
Definition AMReX_ParticleLocator.H:28
Box m_domain
Definition AMReX_ParticleLocator.H:25
Dim3 m_num_bins
Definition AMReX_ParticleLocator.H:23
BinIteratorFactory m_bif
Definition AMReX_ParticleLocator.H:18
GpuArray< Real, 3 > m_dxi
Definition AMReX_ParticleLocator.H:27
__host__ __device__ AssignGrid()=default
Dim3 m_lo
Definition AMReX_ParticleLocator.H:20
GpuArray< Real, 3 > m_plo
Definition AMReX_ParticleLocator.H:26
AssignGrid(BinIteratorFactory a_bif, const IntVect &a_bins_lo, const IntVect &a_bins_hi, const IntVect &a_bin_size, const IntVect &a_num_bins, const Geometry &a_geom, bool a_do_tiling, const IntVect &a_tile_size)
Definition AMReX_ParticleLocator.H:34
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
Fixed-size array that can be used on GPU.
Definition AMReX_Array.H:43
Definition AMReX_Reduce.H:348
Definition AMReX_Reduce.H:314