8 #ifndef WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
9 #define WARPX_PARTICLES_SORTING_SORTINGUTILS_H_
34 template<
typename ForwardIterator >
36 ForwardIterator
const index_end,
42 int N =
static_cast<int>(std::distance(index_begin, index_end));
46 ForwardIterator sep = index_begin;
47 std::advance(sep, num_true);
50 ForwardIterator
const sep = std::stable_partition(
51 index_begin, index_end,
52 [&predicate](
long i) {
return predicate[
i]; }
65 template<
typename ForwardIterator >
67 ForwardIterator
const last)
69 return std::distance( first, last );
94 for (
int idim=0; idim<AMREX_SPACEDIM; idim++) {
145 long const start_index ) :
154 for (
int idim=0; idim<AMREX_SPACEDIM; idim++) {
191 template <
typename T>
#define AMREX_FORCE_INLINE
#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_ParticleType ParticleType
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