1#ifndef AMREX_NEIGHBORPARTICLESGPUIMPL_H_
2#define AMREX_NEIGHBORPARTICLESGPUIMPL_H_
3#include <AMReX_Config.H>
12 void forEachIntersectingTile (
IntVect const& iv,
int nGrow,
13 Box const& grid_box,
IntVect const& periodic_shift,
14 bool do_tiling,
IntVect const& tile_size,
F&& f)
18 cell_box += periodic_shift;
21 if (!cell_box.ok()) {
return; }
23 auto const& cb_lo = cell_box.smallEnd();
24 auto const& cb_hi = cell_box.bigEnd();
26#if (AMREX_SPACEDIM == 1)
27 for (
int i = cb_lo[0]; i <= cb_hi[0]; ++i) {
30 int tile =
getTileIndex(cell, grid_box, do_tiling, tile_size, tbx);
36#elif (AMREX_SPACEDIM == 2)
37 for (
int j = cb_lo[1]; j <= cb_hi[1]; ++j) {
38 for (
int i = cb_lo[0]; i <= cb_hi[0]; ++i) {
41 int tile =
getTileIndex(cell, grid_box, do_tiling, tile_size, tbx);
51 for (
int k = cb_lo[2]; k <= cb_hi[2]; ++k) {
52 for (
int j = cb_lo[1]; j <= cb_hi[1]; ++j) {
53 for (
int i = cb_lo[0]; i <= cb_hi[0]; ++i) {
56 int tile =
getTileIndex(cell, grid_box, do_tiling, tile_size, tbx);
69 inline Vector<Box> getBoundaryBoxes(
const Box& box,
int ncells)
72 "Too many cells requested in getBoundaryBoxes");
75 "Box must be cell-centered");
78 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
82 face_boxes.push_back(hi_face_box); bl.push_back(hi_face_box);
83 face_boxes.push_back(lo_face_box); bl.push_back(lo_face_box);
84 for (
auto face_box : face_boxes) {
85 for (
int j = 0; j < AMREX_SPACEDIM; ++j) {
86 if (i == j) {
continue; }
90 edge_boxes.push_back(hi_edge_box); bl.push_back(hi_edge_box);
91 edge_boxes.push_back(lo_edge_box); bl.push_back(lo_edge_box);
92 for (
auto edge_box : edge_boxes) {
93 for (
int k = 0; k < AMREX_SPACEDIM; ++k) {
94 if ((j == k) || (i == k)) {
continue; }
97 bl.push_back(hi_corner_box);
98 bl.push_back(lo_corner_box);
111template <
int NStructReal,
int NStructInt,
int NArrayReal,
int NArrayInt>
116 BL_PROFILE(
"NeighborParticleContainer<NStructReal, NStructInt, NArrayReal, NArrayInt>::buildNeighborMask");
117 m_neighbor_mask_initialized =
true;
119 const Geometry& geom = this->Geom(lev);
120 const BoxArray& ba = this->ParticleBoxArray(lev);
125 if (!m_neighbor_mask_initialized ||
130 const std::vector<IntVect>& pshifts = periodicity.
shiftIntVect();
134 int grid = mfi.index();
136 std::set<NeighborTask> neighbor_grids;
137 for (
auto pshift : pshifts)
139 const Box box = ba[mfi] + pshift;
141 const bool first_only =
false;
142 auto isecs = ba.
intersections(box, first_only, m_num_neighbor_cells);
144 for (
auto& isec : isecs)
146 int nbor_grid = isec.first;
147 const Box isec_box = isec.second - pshift;
148 neighbor_grids.insert(
NeighborTask(nbor_grid, isec_box, pshift));
149 const int global_rank = dmap[nbor_grid];
156 for (
auto nbor_grid : neighbor_grids)
159 code.
grid_id = nbor_grid.grid_id;
160 code.
grid_box = ba[nbor_grid.grid_id];
166 m_code_array[grid].resize(h_code_arr.
size());
168 m_code_array[grid].begin());
169 m_isec_boxes[grid].resize(h_isec_boxes.
size());
171 m_isec_boxes[grid].begin());
176 m_cached_neighbor_ba = ba;
177 m_cached_neighbor_dm = dmap;
178 m_neighbor_mask_initialized =
true;
183template <
int NStructReal,
int NStructInt,
int NArrayReal,
int NArrayInt>
188 BL_PROFILE(
"NeighborParticleContainer<NStructReal, NStructInt, NArrayReal, NArrayInt>::buildNeighborCopyOp()");
193 const auto& geom = this->Geom(lev);
194 const auto dxi = this->Geom(lev).InvCellSizeArray();
195 const auto plo = this->Geom(lev).ProbLoArray();
196 const auto domain = this->Geom(lev).Domain();
197 auto& plev = this->GetParticles(lev);
198 auto& ba = this->ParticleBoxArray(lev);
200 if (ba.size() == 1 && (! geom.isAnyPeriodic()) ) {
return; }
202 for(
MFIter mfi = this->MakeMFIter(lev); mfi.
isValid(); ++mfi)
204 int gid = mfi.index();
205 int tid = mfi.LocalTileIndex();
206 auto index = std::make_pair(gid, tid);
208 auto& src_tile = plev[index];
209 auto& aos = src_tile.GetArrayOfStructs();
210 const size_t np_real = aos.numParticles();
213 if (use_boundary_neighbor) {
214 np = m_boundary_particle_ids[lev][index].size();
217 m_boundary_particle_ids.resize(1);
218 m_boundary_particle_ids[lev][index];
221 const auto* p_bndry_pid = m_boundary_particle_ids[lev][index].dataPtr();
225 auto p_counts = counts.
dataPtr();
226 auto p_offsets = offsets.
dataPtr();
229 auto p_code_array = m_code_array[gid].dataPtr();
230 auto p_isec_boxes = m_isec_boxes[gid].dataPtr();
231 const int nisec_box = m_isec_boxes[gid].size();
232 const bool do_tiling = this->do_tiling;
233 const IntVect tile_size = this->tile_size;
234 const int nGrow = m_num_neighbor_cells;
241 if (use_boundary_neighbor) {
242 pid = p_bndry_pid[i];
245 for (
int j=0; j<nisec_box; ++j) {
246 if (p_isec_boxes[j].contains(iv)) {
247 detail::forEachIntersectingTile(iv, nGrow,
248 p_code_array[j].grid_box,
249 p_code_array[j].periodic_shift,
250 do_tiling, tile_size,
253 bool is_self = (p_code_array[j].grid_id == gid) && (dst_tile == tid)
254 AMREX_D_TERM( && (p_code_array[j].periodic_shift[0] == 0),
255 && (p_code_array[j].periodic_shift[1] == 0),
256 && (p_code_array[j].periodic_shift[2] == 0));
271 neighbor_copy_op.resize(gid, tid, lev, num_copies);
273 auto tile_index = std::make_pair(gid, tid);
274 auto p_boxes = neighbor_copy_op.m_boxes[lev][tile_index].dataPtr();
275 auto p_levs = neighbor_copy_op.m_levels[lev][tile_index].dataPtr();
276 auto p_tiles = neighbor_copy_op.m_tiles[lev][tile_index].dataPtr();
277 auto p_src_indices = neighbor_copy_op.m_src_indices[lev][tile_index].dataPtr();
278 auto p_periodic_shift = neighbor_copy_op.m_periodic_shift[lev][tile_index].dataPtr();
284 if (use_boundary_neighbor) {
285 pid = p_bndry_pid[i];
288 int k = p_offsets[i];
289 for (
int j=0; j<nisec_box; ++j) {
290 if (p_isec_boxes[j].contains(iv)) {
291 detail::forEachIntersectingTile(iv, nGrow,
292 p_code_array[j].grid_box,
293 p_code_array[j].periodic_shift,
294 do_tiling, tile_size,
297 bool is_self = (p_code_array[j].grid_id == gid) && (dst_tile == tid)
298 AMREX_D_TERM( && (p_code_array[j].periodic_shift[0] == 0),
299 && (p_code_array[j].periodic_shift[1] == 0),
300 && (p_code_array[j].periodic_shift[2] == 0));
302 p_boxes[k] = p_code_array[j].grid_id;
303 p_tiles[k] = dst_tile;
305 p_periodic_shift[k] = p_code_array[j].periodic_shift;
306 p_src_indices[k] = pid;
319template <
int NStructReal,
int NStructInt,
int NArrayReal,
int NArrayInt>
324 BL_PROFILE(
"NeighborParticleContainer::fillNeighbors");
329 this->defineBufferMap();
331 neighbor_copy_op.clear();
332 neighbor_copy_plan.clear();
333 buildNeighborCopyOp();
334 neighbor_copy_plan.build(*
this, neighbor_copy_op, ghost_int_comp,
335 ghost_real_comp,
true);
336 updateNeighborsGPU(
false);
339template <
int NStructReal,
int NStructInt,
int NArrayReal,
int NArrayInt>
344 BL_PROFILE(
"NeighborParticleContainer::updateNeighborsGPU");
346 if (boundary_neighbors_only) {
347 neighbor_copy_op.clear();
348 neighbor_copy_plan.clear();
349 buildNeighborCopyOp(
true);
350 neighbor_copy_plan.build(*
this, neighbor_copy_op, ghost_int_comp,
351 ghost_real_comp,
true);
356 if (this->use_comms_arena) {
361 packBuffer(*
this, neighbor_copy_op, neighbor_copy_plan, snd_buffer);
364 neighbor_copy_plan.buildMPIFinish(this->BufferMap());
373 if (snd_buffer.arena()->isPinned()) {
374 neighbor_copy_plan.buildMPIFinish(this->BufferMap());
378 pinned_snd_buffer.resize(snd_buffer.size());
380 neighbor_copy_plan.buildMPIFinish(this->BufferMap());
385 rcv_buffer.resize(pinned_rcv_buffer.size());
395template <
int NStructReal,
int NStructInt,
int NArrayReal,
int NArrayInt>
400 BL_PROFILE(
"NeighborParticleContainer::clearNeighborsGPU");
404 for (
int lev = 0; lev < this->numLevels(); ++lev)
406 for(
MFIter mfi = this->MakeMFIter(lev); mfi.
isValid(); ++mfi)
408 auto& ptile = this->DefineAndReturnParticleTile(lev, mfi);
409 ptile.setNumNeighbors(0);
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:37
#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_FOR_1D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:97
#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
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
std::vector< std::pair< int, Box > > intersections(const Box &bx) const
Return intersections of Box and BoxArray.
Definition AMReX_BoxArray.cpp:1186
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
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static bool SameRefs(const DistributionMapping &lhs, const DistributionMapping &rhs)
Definition AMReX_DistributionMapping.H:189
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
bool isAnyPeriodic() const noexcept
Is domain periodic in any direction?
Definition AMReX_Geometry.H:334
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:356
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
void clearNeighborsGPU()
Definition AMReX_NeighborParticlesGPUImpl.H:398
typename ParticleContainerType::ParticleType ParticleType
Definition AMReX_NeighborParticles.H:40
void buildNeighborCopyOp(bool use_boundary_neighbor=false)
Definition AMReX_NeighborParticlesGPUImpl.H:186
void updateNeighborsGPU(bool boundary_neighbors_only=false)
Definition AMReX_NeighborParticlesGPUImpl.H:342
void buildNeighborMask()
Definition AMReX_NeighborParticlesGPUImpl.H:114
void fillNeighborsGPU()
Definition AMReX_NeighborParticlesGPUImpl.H:322
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
size_type size() const noexcept
Definition AMReX_PODVector.H:648
iterator begin() noexcept
Definition AMReX_PODVector.H:674
iterator end() noexcept
Definition AMReX_PODVector.H:678
T * dataPtr() noexcept
Definition AMReX_PODVector.H:670
T * data() noexcept
Definition AMReX_PODVector.H:666
void push_back(const T &a_value)
Definition AMReX_PODVector.H:629
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
std::vector< IntVect > shiftIntVect(IntVect const &nghost=IntVect(0)) const
Definition AMReX_Periodicity.cpp:8
OutIter exclusive_scan(InIter begin, InIter end, OutIter result)
Definition AMReX_Scan.H:1442
Arena * The_Comms_Arena()
Definition AMReX_Arena.cpp:880
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 dtoh_memcpy_async(void *p_h, const void *p_d, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:435
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:421
int global_to_local_rank(int rank) noexcept
Definition AMReX_ParallelContext.H:98
bool UseGpuAwareMpi()
Definition AMReX_ParallelDescriptor.H:113
Definition AMReX_Amr.cpp:49
__host__ __device__ BoxND< dim > adjCellHi(const BoxND< dim > &b, int dir, int len=1) noexcept
Similar to adjCellLo but builds an adjacent BoxND on the high end.
Definition AMReX_Box.H:1735
__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
void communicateParticlesStart(const PC &pc, ParticleCopyPlan &plan, const SndBuffer &snd_buffer, RcvBuffer &rcv_buffer)
Definition AMReX_ParticleCommunication.H:908
__host__ __device__ BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the cell centered BoxND of length len adjacent to b on the low end along the coordinate direct...
Definition AMReX_Box.H:1714
void unpackRemotes(PC &pc, const ParticleCopyPlan &plan, Buffer &rcv_buffer, UnpackPolicy const &policy)
Definition AMReX_ParticleCommunication.H:1008
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
void communicateParticlesFinish(const ParticleCopyPlan &plan)
Definition AMReX_ParticleCommunication.cpp:443
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
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:35
void RemoveDuplicates(Vector< T > &vec)
Definition AMReX_Vector.H:211
void unpackBuffer(PC &pc, const ParticleCopyPlan &plan, const Buffer &snd_buffer, UnpackPolicy const &policy)
Definition AMReX_ParticleCommunication.H:781
void packBuffer(const PC &pc, const ParticleCopyOp &op, const ParticleCopyPlan &plan, Buffer &snd_buffer)
Definition AMReX_ParticleCommunication.H:581
__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:337
Definition AMReX_NeighborParticles.H:17
Box grid_box
Definition AMReX_NeighborParticles.H:19
IntVect periodic_shift
Definition AMReX_NeighborParticles.H:20
int grid_id
Definition AMReX_NeighborParticles.H:18
Definition AMReX_NeighborParticles.H:434
Definition AMReX_ParticleCommunication.H:26