#ifndef WARPX_PARTICLES_SORTING_SORTINGUTILS_H_ #define WARPX_PARTICLES_SORTING_SORTINGUTILS_H_ #include #include #include #ifdef AMREX_USE_GPU #include #include #endif /* \brief Fill the elements of the input vector with consecutive integer, * starting from 0 * * \param[inout] v Vector of integers, to be filled by this routine */ void fillWithConsecutiveIntegers( amrex::Gpu::DeviceVector& v ) { #ifdef AMREX_USE_GPU // On GPU: Use thrust thrust::sequence( v.begin(), v.end() ); #else // On CPU: Use std library std::iota( v.begin(), v.end(), 0L ); #endif } /* \brief Find the indices that would reorder the elements of `predicate` * so that the elements with non-zero value precede the other elements * * \param[in, out] index_begin Point to the beginning of the vector which is * to be filled with these indices * \param[in, out] index_begin Point to the end of the vector which is * to be filled with these indices * \param[in] Vector that indicates the elements that need to be reordered first */ template< typename ForwardIterator > ForwardIterator stablePartition(ForwardIterator const index_begin, ForwardIterator const index_end, amrex::Gpu::DeviceVector const& predicate) { #ifdef AMREX_USE_GPU // On GPU: Use thrust int const* AMREX_RESTRICT predicate_ptr = predicate.dataPtr(); ForwardIterator const sep = thrust::stable_partition( thrust::cuda::par(amrex::Cuda::The_ThrustCachedAllocator()), index_begin, index_end, [predicate_ptr] AMREX_GPU_DEVICE (long i) { return predicate_ptr[i]; } ); #else // On CPU: Use std library ForwardIterator const sep = std::stable_partition( index_begin, index_end, [&predicate](long i) { return predicate[i]; } ); #endif return sep; } /* \brief Return the number of elements between `first` and `last` * * \param[in] fist Points to a position in a vector * \param[in] last Points to another position in a vector * \return The number of elements between `first` and `last` */ template< typename ForwardIterator > int iteratorDistance(ForwardIterator const first, ForwardIterator const last) { #ifdef AMREX_USE_GPU // On GPU: Use thrust return thrust::distance( first, last ); #else return std::distance( first, last ); #endif } /* \brief Functor that fills the elements of the particle array `inexflag` * with the value of the spatial array `bmasks`, at the corresponding particle position. * * \param[in] pti Contains information on the particle positions * \param[in] bmasks Spatial array, that contains a flag indicating * whether each cell is part of the gathering/deposition buffers * \param[out] inexflag Vector to be filled with the value of `bmasks` * \param[in] geom Geometry object, necessary to locate particles within the array `bmasks` * */ class fillBufferFlag { public: fillBufferFlag( WarpXParIter const& pti, amrex::iMultiFab const* bmasks, amrex::Gpu::DeviceVector& inexflag, amrex::Geometry const& geom ) { // Extract simple structure that can be used directly on the GPU m_particles = &(pti.GetArrayOfStructs()[0]); m_buffer_mask = (*bmasks)[pti].array(); m_inexflag_ptr = inexflag.dataPtr(); m_domain = geom.Domain(); for (int idim=0; idim m_prob_lo; amrex::GpuArray m_inv_cell_size; amrex::Box m_domain; int* m_inexflag_ptr; WarpXParticleContainer::ParticleType const* m_particles; amrex::Array4 m_buffer_mask; }; /* \brief Functor that fills the elements of the particle array `inexflag` * with the value of the spatial array `bmasks`, at the corresponding particle position. * * Contrary to `fillBufferFlag`, here this is done only for the particles that * the last elements of `particle_indices` point to (from the element at * index `start_index` in `particle_indices`, to the last element of `particle_indices`) * * \param[in] pti Contains information on the particle positions * \param[in] bmasks Spatial array, that contains a flag indicating * whether each cell is part of the gathering/deposition buffers * \param[out] inexflag Vector to be filled with the value of `bmasks` * \param[in] geom Geometry object, necessary to locate particles within the array `bmasks` * \param[in] start_index Index that which elements start to be modified */ class fillBufferFlagRemainingParticles { public: fillBufferFlagRemainingParticles( WarpXParIter const& pti, amrex::iMultiFab const* bmasks, amrex::Gpu::DeviceVector& inexflag, amrex::Geometry const& geom, amrex::Gpu::DeviceVector const& particle_indices, long const start_index ) : m_start_index(start_index) { // Extract simple structure that can be used directly on the GPU m_particles = &(pti.GetArrayOfStructs()[0]); m_buffer_mask = (*bmasks)[pti].array(); m_inexflag_ptr = inexflag.dataPtr(); m_indices_ptr = particle_indices.dataPtr(); m_domain = geom.Domain(); for (int idim=0; idim m_prob_lo; amrex::GpuArray m_inv_cell_size; amrex::Box m_domain; int* m_inexflag_ptr; WarpXParticleContainer::ParticleType const* m_particles; amrex::Array4 m_buffer_mask; long const m_start_index; long const* m_indices_ptr; }; /* \brief Functor that copies the elements of `src` into `dst`, * while reordering them according to `indices` * * \param[in] src Source vector * \param[out] dst Destination vector * \param[in] indices Array of indices that indicate how to reorder elements */ template class copyAndReorder { public: copyAndReorder( amrex::Gpu::ManagedDeviceVector const& src, amrex::Gpu::ManagedDeviceVector& dst, amrex::Gpu::DeviceVector const& indices ) { // Extract simple structure that can be used directly on the GPU m_src_ptr = src.dataPtr(); m_dst_ptr = dst.dataPtr(); m_indices_ptr = indices.dataPtr(); }; AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()( const long ip ) const { m_dst_ptr[ip] = m_src_ptr[ m_indices_ptr[ip] ]; }; private: T const* m_src_ptr; T* m_dst_ptr; long const* m_indices_ptr; }; #endif // WARPX_PARTICLES_SORTING_SORTINGUTILS_H_