WarpX
SortingUtils.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Andrew Myers, Maxence Thevenet, Remi Lehe
2  * Weiqun Zhang
3  *
4  * This file is part of WarpX.
5  *
6  * License: BSD-3-Clause-LBNL
7  */
8 #ifndef WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
9 #define WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
10 
12 
13 #include <AMReX_Gpu.H>
14 #include <AMReX_Partition.H>
15 
16 
23 
34 template< typename ForwardIterator >
35 ForwardIterator stablePartition(ForwardIterator const index_begin,
36  ForwardIterator const index_end,
37  amrex::Gpu::DeviceVector<int> const& predicate)
38 {
39 #ifdef AMREX_USE_GPU
40  // On GPU: Use amrex
41  int const* AMREX_RESTRICT predicate_ptr = predicate.dataPtr();
42  int N = static_cast<int>(std::distance(index_begin, index_end));
43  auto num_true = amrex::StablePartition(&(*index_begin), N,
44  [predicate_ptr] AMREX_GPU_DEVICE (long i) { return predicate_ptr[i]; });
45 
46  ForwardIterator sep = index_begin;
47  std::advance(sep, num_true);
48 #else
49  // On CPU: Use std library
50  ForwardIterator const sep = std::stable_partition(
51  index_begin, index_end,
52  [&predicate](long i) { return predicate[i]; }
53  );
54 #endif
55  return sep;
56 }
57 
65 template< typename ForwardIterator >
66 int iteratorDistance(ForwardIterator const first,
67  ForwardIterator const last)
68 {
69  return std::distance( first, last );
70 }
71 
83 {
84  public:
85  fillBufferFlag( WarpXParIter const& pti, amrex::iMultiFab const* bmasks,
87  amrex::Geometry const& geom ):
88  // Extract simple structure that can be used directly on the GPU
89  m_domain{geom.Domain()},
90  m_inexflag_ptr{inexflag.dataPtr()},
91  m_particles{pti.GetArrayOfStructs().data()},
92  m_buffer_mask{(*bmasks)[pti].array()}
93  {
94  for (int idim=0; idim<AMREX_SPACEDIM; idim++) {
95  m_prob_lo[idim] = geom.ProbLo(idim);
96  m_inv_cell_size[idim] = geom.InvCellSize(idim);
97  }
98  }
99 
100 
102  void operator()( const long i ) const {
103  // Select a particle
104  auto const& p = m_particles[i];
105  // Find the index of the cell where this particle is located
108  // Find the value of the buffer flag in this cell and
109  // store it at the corresponding particle position in the array `inexflag`
111  }
112 
113  private:
120 };
121 
137 {
138  public:
140  WarpXParIter const& pti,
141  amrex::iMultiFab const* bmasks,
143  amrex::Geometry const& geom,
144  amrex::Gpu::DeviceVector<long> const& particle_indices,
145  long const start_index ) :
146  m_domain{geom.Domain()},
147  // Extract simple structure that can be used directly on the GPU
148  m_inexflag_ptr{inexflag.dataPtr()},
149  m_particles{pti.GetArrayOfStructs().data()},
150  m_buffer_mask{(*bmasks)[pti].array()},
151  m_start_index{start_index},
152  m_indices_ptr{particle_indices.dataPtr()}
153  {
154  for (int idim=0; idim<AMREX_SPACEDIM; idim++) {
155  m_prob_lo[idim] = geom.ProbLo(idim);
156  m_inv_cell_size[idim] = geom.InvCellSize(idim);
157  }
158  }
159 
160 
162  void operator()( const long i ) const {
163  // Select a particle
164  auto const& p = m_particles[m_indices_ptr[i+m_start_index]];
165  // Find the index of the cell where this particle is located
168  // Find the value of the buffer flag in this cell and
169  // store it at the corresponding particle position in the array `inexflag`
171  }
172 
173  private:
180  long const m_start_index;
181  long const* m_indices_ptr;
182 };
183 
191 template <typename T>
193 {
194  public:
196  amrex::Gpu::DeviceVector<T> const& src,
198  amrex::Gpu::DeviceVector<long> const& indices ):
199  // Extract simple structure that can be used directly on the GPU
200  m_src_ptr{src.dataPtr()},
201  m_dst_ptr{dst.dataPtr()},
202  m_indices_ptr{indices.dataPtr()}
203  {}
204 
206  void operator()( const long ip ) const {
207  m_dst_ptr[ip] = m_src_ptr[ m_indices_ptr[ip] ];
208  }
209 
210  private:
211  T const* m_src_ptr;
213  long const* m_indices_ptr;
214 };
215 
216 #endif // WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_GPU_DEVICE
#define AMREX_GPU_HOST_DEVICE
void fillWithConsecutiveIntegers(amrex::Gpu::DeviceVector< long > &v)
Fill the elements of the input vector with consecutive integer, starting from 0.
Definition: SortingUtils.cpp:11
ForwardIterator stablePartition(ForwardIterator const index_begin, ForwardIterator const index_end, amrex::Gpu::DeviceVector< int > const &predicate)
Find the indices that would reorder the elements of predicate so that the elements with non-zero valu...
Definition: SortingUtils.H:35
int iteratorDistance(ForwardIterator const first, ForwardIterator const last)
Return the number of elements between first and last
Definition: SortingUtils.H:66
Definition: WarpXParticleContainer.H:52
const Real * InvCellSize() const noexcept
const Real * ProbLo() const noexcept
T * dataPtr() noexcept
Functor that copies the elements of src into dst, while reordering them according to indices
Definition: SortingUtils.H:193
T const * m_src_ptr
Definition: SortingUtils.H:211
copyAndReorder(amrex::Gpu::DeviceVector< T > const &src, amrex::Gpu::DeviceVector< T > &dst, amrex::Gpu::DeviceVector< long > const &indices)
Definition: SortingUtils.H:195
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(const long ip) const
Definition: SortingUtils.H:206
long const * m_indices_ptr
Definition: SortingUtils.H:213
T * m_dst_ptr
Definition: SortingUtils.H:212
Functor that fills the elements of the particle array inexflag with the value of the spatial array bm...
Definition: SortingUtils.H:83
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_inv_cell_size
Definition: SortingUtils.H:119
fillBufferFlag(WarpXParIter const &pti, amrex::iMultiFab const *bmasks, amrex::Gpu::DeviceVector< int > &inexflag, amrex::Geometry const &geom)
Definition: SortingUtils.H:85
amrex::Array4< int const > m_buffer_mask
Definition: SortingUtils.H:117
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_prob_lo
Definition: SortingUtils.H:118
WarpXParticleContainer::ParticleType const * m_particles
Definition: SortingUtils.H:116
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i) const
Definition: SortingUtils.H:102
int * m_inexflag_ptr
Definition: SortingUtils.H:115
amrex::Box m_domain
Definition: SortingUtils.H:114
Functor that fills the elements of the particle array inexflag with the value of the spatial array bm...
Definition: SortingUtils.H:137
long const * m_indices_ptr
Definition: SortingUtils.H:181
amrex::Box m_domain
Definition: SortingUtils.H:176
long const m_start_index
Definition: SortingUtils.H:180
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_inv_cell_size
Definition: SortingUtils.H:175
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_prob_lo
Definition: SortingUtils.H:174
fillBufferFlagRemainingParticles(WarpXParIter const &pti, amrex::iMultiFab const *bmasks, amrex::Gpu::DeviceVector< int > &inexflag, amrex::Geometry const &geom, amrex::Gpu::DeviceVector< long > const &particle_indices, long const start_index)
Definition: SortingUtils.H:139
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i) const
Definition: SortingUtils.H:162
amrex::Array4< int const > m_buffer_mask
Definition: SortingUtils.H:179
WarpXParticleContainer::ParticleType const * m_particles
Definition: SortingUtils.H:178
int * m_inexflag_ptr
Definition: SortingUtils.H:177
int StablePartition(T *data, int beg, int end, F &&f)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE IntVect getParticleCell(P const &p, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &plo, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, const Box &domain) noexcept
i
Definition: check_interp_points_and_weights.py:174
data
Definition: run_alltests_1node.py:325