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 
22 void fillWithConsecutiveIntegers( amrex::Gpu::DeviceVector<long>& v )
23 {
24 #ifdef AMREX_USE_GPU
25  // On GPU: Use amrex
26  auto data = v.data();
27  auto N = v.size();
28  AMREX_FOR_1D( N, i, {data[i] = i;});
29 #else
30  // On CPU: Use std library
31  std::iota( v.begin(), v.end(), 0L );
32 #endif
33 }
34 
44 template< typename ForwardIterator >
45 ForwardIterator stablePartition(ForwardIterator const index_begin,
46  ForwardIterator const index_end,
47  amrex::Gpu::DeviceVector<int> const& predicate)
48 {
49 #ifdef AMREX_USE_GPU
50  // On GPU: Use amrex
51  int const* AMREX_RESTRICT predicate_ptr = predicate.dataPtr();
52  int N = static_cast<int>(std::distance(index_begin, index_end));
53  auto num_true = amrex::StablePartition(&(*index_begin), N,
54  [predicate_ptr] AMREX_GPU_DEVICE (long i) { return predicate_ptr[i]; });
55 
56  ForwardIterator sep = index_begin;
57  std::advance(sep, num_true);
58 #else
59  // On CPU: Use std library
60  ForwardIterator const sep = std::stable_partition(
61  index_begin, index_end,
62  [&predicate](long i) { return predicate[i]; }
63  );
64 #endif
65  return sep;
66 }
67 
74 template< typename ForwardIterator >
75 int iteratorDistance(ForwardIterator const first,
76  ForwardIterator const last)
77 {
78  return std::distance( first, last );
79 }
80 
92 {
93  public:
94  fillBufferFlag( WarpXParIter const& pti, amrex::iMultiFab const* bmasks,
95  amrex::Gpu::DeviceVector<int>& inexflag,
96  amrex::Geometry const& geom ) {
97 
98  // Extract simple structure that can be used directly on the GPU
99  m_particles = &(pti.GetArrayOfStructs()[0]);
100  m_buffer_mask = (*bmasks)[pti].array();
101  m_inexflag_ptr = inexflag.dataPtr();
102  m_domain = geom.Domain();
103  for (int idim=0; idim<AMREX_SPACEDIM; idim++) {
104  m_prob_lo[idim] = geom.ProbLo(idim);
105  m_inv_cell_size[idim] = geom.InvCellSize(idim);
106  }
107  }
108 
109 
110  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
111  void operator()( const long i ) const {
112  // Select a particle
113  auto const& p = m_particles[i];
114  // Find the index of the cell where this particle is located
115  amrex::IntVect const iv = amrex::getParticleCell( p,
117  // Find the value of the buffer flag in this cell and
118  // store it at the corresponding particle position in the array `inexflag`
120  }
121 
122  private:
123  amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> m_prob_lo;
124  amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> m_inv_cell_size;
125  amrex::Box m_domain;
128  amrex::Array4<int const> m_buffer_mask;
129 };
130 
146 {
147  public:
149  WarpXParIter const& pti,
150  amrex::iMultiFab const* bmasks,
151  amrex::Gpu::DeviceVector<int>& inexflag,
152  amrex::Geometry const& geom,
153  amrex::Gpu::DeviceVector<long> const& particle_indices,
154  long const start_index ) :
155  m_start_index(start_index) {
156 
157  // Extract simple structure that can be used directly on the GPU
158  m_particles = &(pti.GetArrayOfStructs()[0]);
159  m_buffer_mask = (*bmasks)[pti].array();
160  m_inexflag_ptr = inexflag.dataPtr();
161  m_indices_ptr = particle_indices.dataPtr();
162  m_domain = geom.Domain();
163  for (int idim=0; idim<AMREX_SPACEDIM; idim++) {
164  m_prob_lo[idim] = geom.ProbLo(idim);
165  m_inv_cell_size[idim] = geom.InvCellSize(idim);
166  }
167  }
168 
169 
170  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
171  void operator()( const long i ) const {
172  // Select a particle
173  auto const& p = m_particles[m_indices_ptr[i+m_start_index]];
174  // Find the index of the cell where this particle is located
175  amrex::IntVect const iv = amrex::getParticleCell( p,
177  // Find the value of the buffer flag in this cell and
178  // store it at the corresponding particle position in the array `inexflag`
179  m_inexflag_ptr[m_indices_ptr[i+m_start_index]] = m_buffer_mask(iv);
180  }
181 
182  private:
183  amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> m_prob_lo;
184  amrex::GpuArray<amrex::Real,AMREX_SPACEDIM> m_inv_cell_size;
185  amrex::Box m_domain;
188  amrex::Array4<int const> m_buffer_mask;
189  long const m_start_index;
190  long const* m_indices_ptr;
191 };
192 
200 template <typename T>
202 {
203  public:
205  amrex::Gpu::DeviceVector<T> const& src,
206  amrex::Gpu::DeviceVector<T>& dst,
207  amrex::Gpu::DeviceVector<long> const& indices ) {
208  // Extract simple structure that can be used directly on the GPU
209  m_src_ptr = src.dataPtr();
210  m_dst_ptr = dst.dataPtr();
211  m_indices_ptr = indices.dataPtr();
212  }
213 
214  AMREX_GPU_DEVICE AMREX_FORCE_INLINE
215  void operator()( const long ip ) const {
216  m_dst_ptr[ip] = m_src_ptr[ m_indices_ptr[ip] ];
217  }
218 
219  private:
220  T const* m_src_ptr;
222  long const* m_indices_ptr;
223 };
224 
225 #endif // WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i) const
Definition: SortingUtils.H:111
WarpXParticleContainer::ParticleType const * m_particles
Definition: SortingUtils.H:127
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:148
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i) const
Definition: SortingUtils.H:171
long const * m_indices_ptr
Definition: SortingUtils.H:190
WarpXParticleContainer::ParticleType const * m_particles
Definition: SortingUtils.H:187
data
Definition: run_alltests_1node.py:320
long const m_start_index
Definition: SortingUtils.H:189
Functor that fills the elements of the particle array inexflag with the value of the spatial array bm...
Definition: SortingUtils.H:145
T const * m_src_ptr
Definition: SortingUtils.H:220
int * m_inexflag_ptr
Definition: SortingUtils.H:126
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_inv_cell_size
Definition: SortingUtils.H:184
Functor that copies the elements of src into dst, while reordering them according to indices ...
Definition: SortingUtils.H:201
long const * m_indices_ptr
Definition: SortingUtils.H:222
copyAndReorder(amrex::Gpu::DeviceVector< T > const &src, amrex::Gpu::DeviceVector< T > &dst, amrex::Gpu::DeviceVector< long > const &indices)
Definition: SortingUtils.H:204
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_prob_lo
Definition: SortingUtils.H:183
amrex::Box m_domain
Definition: SortingUtils.H:185
i
Definition: check_interp_points_and_weights.py:171
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:45
amrex::Box m_domain
Definition: SortingUtils.H:125
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_prob_lo
Definition: SortingUtils.H:123
int iteratorDistance(ForwardIterator const first, ForwardIterator const last)
Return the number of elements between first and last
Definition: SortingUtils.H:75
int * m_inexflag_ptr
Definition: SortingUtils.H:186
amrex::Array4< int const > m_buffer_mask
Definition: SortingUtils.H:128
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_inv_cell_size
Definition: SortingUtils.H:124
T * m_dst_ptr
Definition: SortingUtils.H:221
Definition: WarpXParticleContainer.H:76
Functor that fills the elements of the particle array inexflag with the value of the spatial array bm...
Definition: SortingUtils.H:91
void fillWithConsecutiveIntegers(amrex::Gpu::DeviceVector< long > &v)
Fill the elements of the input vector with consecutive integer, starting from 0.
Definition: SortingUtils.H:22
WarpXParticleContainer::ParticleType ParticleType
Definition: CollisionType.cpp:46
amrex::Array4< int const > m_buffer_mask
Definition: SortingUtils.H:188
fillBufferFlag(WarpXParIter const &pti, amrex::iMultiFab const *bmasks, amrex::Gpu::DeviceVector< int > &inexflag, amrex::Geometry const &geom)
Definition: SortingUtils.H:94
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(const long ip) const
Definition: SortingUtils.H:215