Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleLocator.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLE_LOCATOR_H_
2#define AMREX_PARTICLE_LOCATOR_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_ParGDB.H>
8#include <AMReX_Tuple.H>
9
10#include <utility>
11
12namespace amrex
13{
14
15template <class BinIteratorFactory>
17{
18 BinIteratorFactory m_bif;
19
24
29 bool m_do_tiling = false;
30
32 AssignGrid () = default;
33
34 AssignGrid (BinIteratorFactory a_bif,
35 const IntVect& a_bins_lo, const IntVect& a_bins_hi, const IntVect& a_bin_size,
36 const IntVect& a_num_bins, const Geometry& a_geom,
37 bool a_do_tiling, const IntVect& a_tile_size)
38 : m_bif(a_bif),
39 m_lo(a_bins_lo.dim3()), m_hi(a_bins_hi.dim3()), m_bin_size(a_bin_size.dim3()),
40 m_num_bins(a_num_bins.dim3()), m_domain(a_geom.Domain()),
41 m_plo(a_geom.ProbLoArray()), m_dxi(a_geom.InvCellSizeArray()),
42 m_tile_size(a_tile_size.dim3()), m_do_tiling(a_do_tiling)
43 {
44 // clamp bin size and num_bins to 1 for AMREX_SPACEDIM < 3
45 if (m_bin_size.x >= 0) {m_bin_size.x = amrex::max(m_bin_size.x, 1);}
46 if (m_bin_size.y >= 0) {m_bin_size.y = amrex::max(m_bin_size.y, 1);}
47 if (m_bin_size.z >= 0) {m_bin_size.z = amrex::max(m_bin_size.z, 1);}
48
49 if (m_bin_size.x >= 0) {m_num_bins.x = amrex::max(m_num_bins.x, 1);}
50 if (m_bin_size.y >= 0) {m_num_bins.y = amrex::max(m_num_bins.y, 1);}
51 if (m_bin_size.z >= 0) {m_num_bins.z = amrex::max(m_num_bins.z, 1);}
52 }
53
55 int getTile (const IntVect& iv, const Box& bx) const noexcept
56 {
57 Box tbx;
59 return getTileIndex(iv, bx, m_do_tiling, tile_size, tbx);
60 }
61
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
65 {
66 const auto iv = assignor(p, m_plo, m_dxi, m_domain);
67 return this->operator()(iv, nGrow);
68 }
69
71 std::pair<int, int> operator() (const IntVect& iv, int nGrow=0) const noexcept
72 {
73 if (AMREX_D_TERM((m_num_bins.x == 0), && (m_num_bins.y == 0), && (m_num_bins.z == 0))) {
74 return {-1, -1};
75 }
76 const auto lo = iv.dim3();
77 int ix_lo = amrex::max((lo.x - nGrow - m_lo.x) / m_bin_size.x - 1, 0);
78 int iy_lo = amrex::max((lo.y - nGrow - m_lo.y) / m_bin_size.y - 1, 0);
79 int iz_lo = amrex::max((lo.z - nGrow - m_lo.z) / m_bin_size.z - 1, 0);
80
81 int ix_hi = amrex::min((lo.x + nGrow - m_lo.x) / m_bin_size.x, m_num_bins.x-1);
82 int iy_hi = amrex::min((lo.y + nGrow - m_lo.y) / m_bin_size.y, m_num_bins.y-1);
83 int iz_hi = amrex::min((lo.z + nGrow - m_lo.z) / m_bin_size.z, m_num_bins.z-1);
84 int loc_gid = -1;
85 int loc_tid = -1;
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) {
89 int index = (ii * m_num_bins.y + jj) * m_num_bins.z + kk;
90 for (const auto& nbor : m_bif.getBinIterator(index)) {
91 Box bx = nbor.second;
92 if (bx.contains(iv)) {
93 return {nbor.first, getTile(iv, bx)};
94 }
95 Box gbx = bx;
96 gbx.grow(nGrow);
97 if (gbx.contains(iv)) {
98 int tile = getTile(iv, bx);
99 if (loc_gid < 0) {
100 loc_gid = nbor.first;
101 loc_tid = tile;
102 }
103 // Prefer particle not in corner ghost cells
104 for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
105 Box gdbx = bx;
106 gdbx.grow(dir, nGrow);
107 if (gdbx.contains(iv)) {
108 loc_gid = nbor.first;
109 loc_tid = tile;
110 }
111 }
112 }
113 }
114 }
115 }
116 }
117 return {loc_gid, loc_tid};
118 }
119};
120
121template <class Bins>
123{
124public:
125
126 using BinIteratorFactory = typename Bins::BinIteratorFactory;
127
128 ParticleLocator () = default;
129
130 void build (const BoxArray& ba, const Geometry& geom,
131 bool a_do_tiling = false,
132 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000)))
133 {
134 m_defined = true;
135 m_do_tiling = a_do_tiling;
136 m_tile_size = a_tile_size;
137 m_ba = ba;
138 m_geom = geom;
139 int num_boxes = static_cast<int>(ba.size());
140 m_host_boxes.resize(0);
141 for (int i = 0; i < num_boxes; ++i) { m_host_boxes.push_back(ba[i]); }
142
143 m_device_boxes.resize(num_boxes);
145
146 if (num_boxes == 0) {
147 m_bins_lo = IntVect(AMREX_D_DECL( 0, 0, 0));
148 m_bins_hi = IntVect(AMREX_D_DECL(-1, -1, -1));
149 m_bin_size = IntVect(AMREX_D_DECL(-1, -1, -1));
150 m_num_bins = IntVect(AMREX_D_DECL( 0, 0, 0));
151 return;
152 }
153
154 // compute the lo, hi and the max box size in each direction
158 ReduceData<AMREX_D_DECL(int, int, int),
159 AMREX_D_DECL(int, int, int),
160 AMREX_D_DECL(int, int, int)> reduce_data(reduce_op);
161 using ReduceTuple = typename decltype(reduce_data)::Type;
162
163 auto *const boxes_ptr = m_device_boxes.dataPtr();
164 reduce_op.eval(num_boxes, reduce_data,
165 [=] AMREX_GPU_DEVICE (int i) -> ReduceTuple
166 {
167 const Box& box = boxes_ptr[i];
168 IntVect lo = box.smallEnd();
169 IntVect hi = box.bigEnd();
170 IntVect si = box.length();
171 return {AMREX_D_DECL(lo[0], lo[1], lo[2]),
172 AMREX_D_DECL(hi[0], hi[1], hi[2]),
173 AMREX_D_DECL(si[0], si[1], si[2])};
174 });
175
176 ReduceTuple hv = reduce_data.value(reduce_op);
177
178 m_bins_lo = IntVect(AMREX_D_DECL(amrex::get<0>(hv),
179 amrex::get<1>(hv),
180 amrex::get<2>(hv)));
181 m_bins_hi = IntVect(AMREX_D_DECL(amrex::get< AMREX_SPACEDIM >(hv),
182 amrex::get< AMREX_SPACEDIM+1>(hv),
183 amrex::get< AMREX_SPACEDIM+2>(hv)));
184 m_bin_size = IntVect(AMREX_D_DECL(amrex::get<2*AMREX_SPACEDIM>(hv),
185 amrex::get<2*AMREX_SPACEDIM+1>(hv),
186 amrex::get<2*AMREX_SPACEDIM+2>(hv)));
187
189
191 IntVect bin_size = m_bin_size;
192 IntVect bins_lo = m_bins_lo;
193 m_bins.build(num_boxes, boxes_ptr, bins_box,
194 [=] AMREX_GPU_DEVICE (const Box& box) noexcept -> IntVect
195 {
196 return (box.smallEnd() - bins_lo) / bin_size;
197 });
198 }
199
200 void setGeometry (const Geometry& a_geom) noexcept
201 {
203 m_geom = a_geom;
204 }
205
213
214 bool isValid (const BoxArray& ba,
215 bool a_do_tiling = false,
216 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000))) const noexcept
217 {
218 if (m_defined) {
219 return BoxArray::SameRefs(m_ba, ba) &&
220 (m_do_tiling == a_do_tiling) &&
221 (m_tile_size == a_tile_size);
222 }
223 return false;
224 }
225
226protected:
227
228 bool m_defined{false};
229 bool m_do_tiling{false};
230
233 IntVect m_tile_size{AMREX_D_DECL(1024000, 1024000, 1024000)};
234
239
240 Bins m_bins;
241
244};
245
246template <class BinIteratorFactory>
248{
250 std::size_t m_size;
251
252 AmrAssignGrid(const AssignGrid<BinIteratorFactory>* a_funcs, std::size_t a_size)
253 : m_funcs(a_funcs), m_size(a_size)
254 {}
255
256 template <typename P, typename Assignor = DefaultAssignor>
258 GpuTuple<int, int, int> operator() (const P& p, int lev_min=-1, int lev_max=-1, int nGrow=0,
259 Assignor const& assignor = {}) const noexcept
260 {
261 lev_min = (lev_min == -1) ? 0 : lev_min;
262 lev_max = (lev_max == -1) ? m_size - 1 : lev_max;
263
264 for (int lev = lev_max; lev >= lev_min; --lev)
265 {
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); }
268 }
269
270 {
271 int lev = lev_min;
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); }
274 }
275
276 return makeTuple(-1, -1, -1);
277 }
278};
279
280template <class Bins>
282{
283public:
284 using BinIteratorFactory = typename Bins::BinIteratorFactory;
285
286private:
287 Vector<ParticleLocator<Bins> > m_locators;
289 bool m_defined = false;
290
291public:
292
294
296 const Vector<Geometry>& a_geom,
297 bool a_do_tiling = false,
298 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000)))
299 {
300 build(a_ba, a_geom, a_do_tiling, a_tile_size);
301 }
302
304 bool a_do_tiling = false,
305 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000)))
306 {
307 build(a_gdb, a_do_tiling, a_tile_size);
308 }
309
310 void build (const Vector<BoxArray>& a_ba,
311 const Vector<Geometry>& a_geom,
312 bool a_do_tiling = false,
313 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000)))
314 {
315 m_defined = true;
316 int num_levels = static_cast<int>(a_ba.size());
317 m_locators.resize(num_levels);
318 m_grid_assignors.resize(num_levels);
319#ifdef AMREX_USE_GPU
320 Gpu::HostVector<AssignGrid<BinIteratorFactory> > h_grid_assignors(num_levels);
321 for (int lev = 0; lev < num_levels; ++lev)
322 {
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();
325 }
326 Gpu::htod_memcpy_async(m_grid_assignors.data(), h_grid_assignors.data(),
327 sizeof(AssignGrid<BinIteratorFactory>)*num_levels);
329#else
330 for (int lev = 0; lev < num_levels; ++lev)
331 {
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();
334 }
335#endif
336 }
337
338 void build (const ParGDBBase* a_gdb,
339 bool a_do_tiling = false,
340 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000)))
341 {
343 Vector<Geometry> geom;
344 int num_levels = a_gdb->finestLevel()+1;
345 for (int lev = 0; lev < num_levels; ++lev)
346 {
347 ba.push_back(a_gdb->ParticleBoxArray(lev));
348 geom.push_back(a_gdb->Geom(lev));
349 }
350 build(ba, geom, a_do_tiling, a_tile_size);
351 }
352
353 [[nodiscard]] bool isValid (const Vector<BoxArray>& a_ba,
354 bool a_do_tiling = false,
355 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000))) const
356 {
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);
363 }
364 return all_valid;
365 }
366
367 bool isValid (const ParGDBBase* a_gdb,
368 bool a_do_tiling = false,
369 const IntVect& a_tile_size = IntVect(AMREX_D_DECL(1024000, 1024000, 1024000))) const
370 {
372 int num_levels = a_gdb->finestLevel()+1;
373 for (int lev = 0; lev < num_levels; ++lev) {
374 ba.push_back(a_gdb->ParticleBoxArray(lev));
375 }
376 return this->isValid(ba, a_do_tiling, a_tile_size);
377 }
378
379 void setGeometry (const ParGDBBase* a_gdb)
380 {
381 int num_levels = a_gdb->finestLevel()+1;
382#ifdef AMREX_USE_GPU
383 Gpu::HostVector<AssignGrid<BinIteratorFactory> > h_grid_assignors(num_levels);
384 for (int lev = 0; lev < num_levels; ++lev)
385 {
386 m_locators[lev].setGeometry(a_gdb->Geom(lev));
387 h_grid_assignors[lev] = m_locators[lev].getGridAssignor();
388 }
389 Gpu::htod_memcpy_async(m_grid_assignors.data(), h_grid_assignors.data(),
390 sizeof(AssignGrid<BinIteratorFactory>)*num_levels);
392#else
393 for (int lev = 0; lev < num_levels; ++lev)
394 {
395 m_locators[lev].setGeometry(a_gdb->Geom(lev));
396 m_grid_assignors[lev] = m_locators[lev].getGridAssignor();
397 }
398#endif
399 }
400
402 {
403 AMREX_ASSERT(m_defined);
404 return AmrAssignGrid<BinIteratorFactory>(m_grid_assignors.dataPtr(), m_locators.size());
405 }
406};
407
408}
409
410#endif
#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
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
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