ARTEMIS
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 {
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 
45 template< typename ForwardIterator >
46 ForwardIterator stablePartition(ForwardIterator const index_begin,
47  ForwardIterator const index_end,
48  amrex::Gpu::DeviceVector<int> const& predicate)
49 {
50 #ifdef AMREX_USE_GPU
51  // On GPU: Use amrex
52  int const* AMREX_RESTRICT predicate_ptr = predicate.dataPtr();
53  int N = static_cast<int>(std::distance(index_begin, index_end));
54  auto num_true = amrex::StablePartition(&(*index_begin), N,
55  [predicate_ptr] AMREX_GPU_DEVICE (long i) { return predicate_ptr[i]; });
56 
57  ForwardIterator sep = index_begin;
58  std::advance(sep, num_true);
59 #else
60  // On CPU: Use std library
61  ForwardIterator const sep = std::stable_partition(
62  index_begin, index_end,
63  [&predicate](long i) { return predicate[i]; }
64  );
65 #endif
66  return sep;
67 }
68 
76 template< typename ForwardIterator >
77 int iteratorDistance(ForwardIterator const first,
78  ForwardIterator const last)
79 {
80  return std::distance( first, last );
81 }
82 
94 {
95  public:
96  fillBufferFlag( WarpXParIter const& pti, amrex::iMultiFab const* bmasks,
98  amrex::Geometry const& geom ) {
99 
100  // Extract simple structure that can be used directly on the GPU
101  m_particles = &(pti.GetArrayOfStructs()[0]);
102  m_buffer_mask = (*bmasks)[pti].array();
103  m_inexflag_ptr = inexflag.dataPtr();
104  m_domain = geom.Domain();
105  for (int idim=0; idim<AMREX_SPACEDIM; idim++) {
106  m_prob_lo[idim] = geom.ProbLo(idim);
107  m_inv_cell_size[idim] = geom.InvCellSize(idim);
108  }
109  }
110 
111 
113  void operator()( const long i ) const {
114  // Select a particle
115  auto const& p = m_particles[i];
116  // Find the index of the cell where this particle is located
119  // Find the value of the buffer flag in this cell and
120  // store it at the corresponding particle position in the array `inexflag`
122  }
123 
124  private:
131 };
132 
148 {
149  public:
151  WarpXParIter const& pti,
152  amrex::iMultiFab const* bmasks,
154  amrex::Geometry const& geom,
155  amrex::Gpu::DeviceVector<long> const& particle_indices,
156  long const start_index ) :
157  m_start_index(start_index) {
158 
159  // Extract simple structure that can be used directly on the GPU
160  m_particles = &(pti.GetArrayOfStructs()[0]);
161  m_buffer_mask = (*bmasks)[pti].array();
162  m_inexflag_ptr = inexflag.dataPtr();
163  m_indices_ptr = particle_indices.dataPtr();
164  m_domain = geom.Domain();
165  for (int idim=0; idim<AMREX_SPACEDIM; idim++) {
166  m_prob_lo[idim] = geom.ProbLo(idim);
167  m_inv_cell_size[idim] = geom.InvCellSize(idim);
168  }
169  }
170 
171 
173  void operator()( const long i ) const {
174  // Select a particle
175  auto const& p = m_particles[m_indices_ptr[i+m_start_index]];
176  // Find the index of the cell where this particle is located
179  // Find the value of the buffer flag in this cell and
180  // store it at the corresponding particle position in the array `inexflag`
182  }
183 
184  private:
191  long const m_start_index;
192  long const* m_indices_ptr;
193 };
194 
202 template <typename T>
204 {
205  public:
207  amrex::Gpu::DeviceVector<T> const& src,
209  amrex::Gpu::DeviceVector<long> const& indices ) {
210  // Extract simple structure that can be used directly on the GPU
211  m_src_ptr = src.dataPtr();
212  m_dst_ptr = dst.dataPtr();
213  m_indices_ptr = indices.dataPtr();
214  }
215 
217  void operator()( const long ip ) const {
218  m_dst_ptr[ip] = m_src_ptr[ m_indices_ptr[ip] ];
219  }
220 
221  private:
222  T const* m_src_ptr;
224  long const* m_indices_ptr;
225 };
226 
227 #endif // WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
#define AMREX_FORCE_INLINE
#define AMREX_RESTRICT
#define AMREX_FOR_1D(...)
#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.H:22
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:46
int iteratorDistance(ForwardIterator const first, ForwardIterator const last)
Return the number of elements between first and last
Definition: SortingUtils.H:77
Definition: WarpXParticleContainer.H:52
const Real * InvCellSize() const noexcept
const Box & Domain() const noexcept
const Real * ProbLo() const noexcept
size_type size() const noexcept
iterator begin() noexcept
iterator end() noexcept
T * dataPtr() noexcept
T * data() noexcept
AoSRef GetArrayOfStructs() const
Functor that copies the elements of src into dst, while reordering them according to indices
Definition: SortingUtils.H:204
T const * m_src_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:206
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(const long ip) const
Definition: SortingUtils.H:217
long const * m_indices_ptr
Definition: SortingUtils.H:224
T * m_dst_ptr
Definition: SortingUtils.H:223
Functor that fills the elements of the particle array inexflag with the value of the spatial array bm...
Definition: SortingUtils.H:94
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_inv_cell_size
Definition: SortingUtils.H:126
fillBufferFlag(WarpXParIter const &pti, amrex::iMultiFab const *bmasks, amrex::Gpu::DeviceVector< int > &inexflag, amrex::Geometry const &geom)
Definition: SortingUtils.H:96
amrex::Array4< int const > m_buffer_mask
Definition: SortingUtils.H:130
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_prob_lo
Definition: SortingUtils.H:125
WarpXParticleContainer::ParticleType const * m_particles
Definition: SortingUtils.H:129
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i) const
Definition: SortingUtils.H:113
int * m_inexflag_ptr
Definition: SortingUtils.H:128
amrex::Box m_domain
Definition: SortingUtils.H:127
Functor that fills the elements of the particle array inexflag with the value of the spatial array bm...
Definition: SortingUtils.H:148
long const * m_indices_ptr
Definition: SortingUtils.H:192
amrex::Box m_domain
Definition: SortingUtils.H:187
long const m_start_index
Definition: SortingUtils.H:191
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_inv_cell_size
Definition: SortingUtils.H:186
amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > m_prob_lo
Definition: SortingUtils.H:185
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:150
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator()(const long i) const
Definition: SortingUtils.H:173
amrex::Array4< int const > m_buffer_mask
Definition: SortingUtils.H:190
WarpXParticleContainer::ParticleType const * m_particles
Definition: SortingUtils.H:189
int * m_inexflag_ptr
Definition: SortingUtils.H:188
list data
open text file and read data #####
Definition: Excitation_Flag_Generator.py:24
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) noexcept
int StablePartition(T *data, int beg, int end, F &&f)
i
Definition: check_interp_points_and_weights.py:174