16#ifndef CABANA_COMMUNICATIONPLANBASE_HPP
17#define CABANA_COMMUNICATIONPLANBASE_HPP
19#include <Cabana_Core_Config.hpp>
23#include <Kokkos_Core.hpp>
24#include <Kokkos_ScatterView.hpp>
44struct CountSendsAndCreateSteeringDuplicated
47struct CountSendsAndCreateSteeringAtomic
53template <
class ExecutionSpace>
54struct CountSendsAndCreateSteeringAlgorithm;
57#ifdef KOKKOS_ENABLE_CUDA
59struct CountSendsAndCreateSteeringAlgorithm<Kokkos::Cuda>
61 using type = CountSendsAndCreateSteeringAtomic;
64#ifdef KOKKOS_ENABLE_HIP
66struct CountSendsAndCreateSteeringAlgorithm<Kokkos::HIP>
68 using type = CountSendsAndCreateSteeringAtomic;
71#ifdef KOKKOS_ENABLE_SYCL
73struct CountSendsAndCreateSteeringAlgorithm<Kokkos::Experimental::SYCL>
75 using type = CountSendsAndCreateSteeringAtomic;
78#ifdef KOKKOS_ENABLE_OPENMPTARGET
80struct CountSendsAndCreateSteeringAlgorithm<Kokkos::Experimental::OpenMPTarget>
82 using type = CountSendsAndCreateSteeringAtomic;
87template <
class ExecutionSpace>
88struct CountSendsAndCreateSteeringAlgorithm
90 using type = CountSendsAndCreateSteeringDuplicated;
95template <
class ExecutionSpace,
class ExportRankView>
96auto countSendsAndCreateSteering( ExecutionSpace,
97 const ExportRankView element_export_ranks,
99 CountSendsAndCreateSteeringAtomic )
100 -> std::pair<Kokkos::View<int*, typename ExportRankView::memory_space>,
101 Kokkos::View<
typename ExportRankView::size_type*,
102 typename ExportRankView::memory_space>>
104 using memory_space =
typename ExportRankView::memory_space;
105 using size_type =
typename ExportRankView::size_type;
108 Kokkos::View<int*, memory_space> neighbor_counts(
"neighbor_counts",
110 Kokkos::View<size_type*, memory_space> neighbor_ids(
111 Kokkos::ViewAllocateWithoutInitializing(
"neighbor_ids" ),
112 element_export_ranks.size() );
118 if ( comm_size <= 64 )
120 constexpr int team_size = 256;
121 Kokkos::TeamPolicy<ExecutionSpace> team_policy(
122 ( element_export_ranks.size() + team_size - 1 ) / team_size,
124 team_policy = team_policy.set_scratch_size(
126 Kokkos::PerTeam(
sizeof(
int ) * ( team_size + 2 * comm_size ) ) );
127 Kokkos::parallel_for(
128 "Cabana::CommunicationPlan::countSendsAndCreateSteeringShared",
131 const typename Kokkos::TeamPolicy<ExecutionSpace>::member_type&
139 int* scratch = (
int*)team.team_shmem().get_shmem(
140 ( team.team_size() + 2 * comm_size ) *
sizeof(
int ), 0 );
144 int* local_neighbor_ids = scratch;
146 int* histo = local_neighbor_ids + team.team_size();
148 int* global_offset = histo + comm_size;
151 team.team_rank() + team.league_rank() * team.team_size();
153 const int local_id = team.team_rank();
155 const int num_elements = element_export_ranks.size();
157 const int my_element_export_rank =
158 ( tid < num_elements ? element_export_ranks( tid ) : -1 );
162 const bool in_bounds =
163 tid < num_elements && my_element_export_rank >= 0;
165 Kokkos::parallel_for(
166 Kokkos::TeamThreadRange( team, comm_size ),
167 [&](
const int i ) { histo[i] = 0; } );
176 local_neighbor_ids[local_id] = Kokkos::atomic_fetch_add(
177 &histo[my_element_export_rank], 1 );
184 Kokkos::parallel_for(
185 Kokkos::TeamThreadRange( team, comm_size ),
189 global_offset[i] = Kokkos::atomic_fetch_add(
190 &neighbor_counts( i ), histo[i] );
199 neighbor_ids( tid ) =
200 global_offset[my_element_export_rank] +
201 local_neighbor_ids[local_id];
210 Kokkos::parallel_for(
211 "Cabana::CommunicationPlan::countSendsAndCreateSteering",
212 Kokkos::RangePolicy<ExecutionSpace>( 0,
213 element_export_ranks.size() ),
214 KOKKOS_LAMBDA(
const size_type i ) {
215 if ( element_export_ranks( i ) >= 0 )
216 neighbor_ids( i ) = Kokkos::atomic_fetch_add(
217 &neighbor_counts( element_export_ranks( i ) ), 1 );
223 return std::make_pair( neighbor_counts, neighbor_ids );
228template <
class ExecutionSpace,
class ExportRankView>
229auto countSendsAndCreateSteering( ExecutionSpace,
230 const ExportRankView element_export_ranks,
232 CountSendsAndCreateSteeringDuplicated )
233 -> std::pair<Kokkos::View<int*, typename ExportRankView::memory_space>,
234 Kokkos::View<
typename ExportRankView::size_type*,
235 typename ExportRankView::memory_space>>
237 using memory_space =
typename ExportRankView::memory_space;
238 using size_type =
typename ExportRankView::size_type;
241 Kokkos::Experimental::UniqueToken<
242 ExecutionSpace, Kokkos::Experimental::UniqueTokenScope::Global>
246 Kokkos::View<int*, memory_space> neighbor_counts(
247 Kokkos::ViewAllocateWithoutInitializing(
"neighbor_counts" ),
249 Kokkos::View<size_type*, memory_space> neighbor_ids(
250 Kokkos::ViewAllocateWithoutInitializing(
"neighbor_ids" ),
251 element_export_ranks.size() );
252 Kokkos::View<int**, memory_space> neighbor_counts_dup(
253 "neighbor_counts", unique_token.size(), comm_size );
254 Kokkos::View<size_type**, memory_space> neighbor_ids_dup(
255 "neighbor_ids", unique_token.size(), element_export_ranks.size() );
258 Kokkos::parallel_for(
259 "Cabana::CommunicationPlan::intialCount",
260 Kokkos::RangePolicy<ExecutionSpace>( 0, element_export_ranks.size() ),
261 KOKKOS_LAMBDA(
const size_type i ) {
262 if ( element_export_ranks( i ) >= 0 )
265 auto thread_id = unique_token.acquire();
275 neighbor_ids_dup( thread_id, i ) = ++neighbor_counts_dup(
276 thread_id, element_export_ranks( i ) );
279 unique_token.release( thread_id );
286 Kokkos::TeamPolicy<ExecutionSpace, Kokkos::Schedule<Kokkos::Dynamic>>;
287 using index_type =
typename team_policy::index_type;
291 Kokkos::parallel_for(
292 "Cabana::CommunicationPlan::finalCount",
293 team_policy( neighbor_counts.extent( 0 ), Kokkos::AUTO ),
294 KOKKOS_LAMBDA(
const typename team_policy::member_type& team ) {
296 auto i = team.league_rank();
299 int thread_counts = 0;
300 Kokkos::parallel_reduce(
301 Kokkos::TeamThreadRange( team,
302 neighbor_counts_dup.extent( 0 ) ),
303 [&](
const index_type thread_id,
int& result )
304 { result += neighbor_counts_dup( thread_id, i ); },
306 neighbor_counts( i ) = thread_counts;
312 Kokkos::parallel_for(
313 "Cabana::CommunicationPlan::createSteering",
314 team_policy( element_export_ranks.size(), Kokkos::AUTO ),
315 KOKKOS_LAMBDA(
const typename team_policy::member_type& team ) {
317 auto i = team.league_rank();
320 if ( element_export_ranks( i ) >= 0 )
325 index_type dup_thread = 0;
326 Kokkos::parallel_reduce(
327 Kokkos::TeamThreadRange( team,
328 neighbor_ids_dup.extent( 0 ) ),
329 [&](
const index_type thread_id, index_type& result )
331 if ( neighbor_ids_dup( thread_id, i ) > 0 )
341 size_type thread_offset = 0;
342 Kokkos::parallel_reduce(
343 Kokkos::TeamThreadRange( team, dup_thread ),
344 [&](
const index_type thread_id, size_type& result ) {
345 result += neighbor_counts_dup(
346 thread_id, element_export_ranks( i ) );
354 thread_offset + neighbor_ids_dup( dup_thread, i ) - 1;
360 return std::make_pair( neighbor_counts, neighbor_ids );
375 std::vector<int> topology )
377 auto remove_end = std::remove( topology.begin(), topology.end(), -1 );
378 std::sort( topology.begin(), remove_end );
379 auto unique_end = std::unique( topology.begin(), remove_end );
380 topology.resize( std::distance( topology.begin(), unique_end ) );
384 MPI_Comm_rank( comm, &my_rank );
385 for (
auto& n : topology )
389 std::swap( n, topology[0] );
422template <
class MemorySpace>
428 static_assert( Kokkos::is_memory_space<MemorySpace>() );
434 using size_type =
typename memory_space::memory_space::size_type;
449 auto p = std::make_unique<MPI_Comm>();
450 MPI_Comm_dup(
comm, p.get() );
568 template <
class PackViewType,
class RankViewType>
570 const RankViewType& element_export_ranks )
573 createSteering(
true, neighbor_ids, element_export_ranks,
574 element_export_ranks );
596 template <
class PackViewType,
class RankViewType,
class IdViewType>
598 const RankViewType& element_export_ranks,
599 const IdViewType& element_export_ids )
601 createSteering(
false, neighbor_ids, element_export_ranks,
602 element_export_ids );
607 template <
class ExecutionSpace,
class PackViewType,
class RankViewType,
609 void createSteering( ExecutionSpace,
const bool use_iota,
610 const PackViewType& neighbor_ids,
611 const RankViewType& element_export_ranks,
612 const IdViewType& element_export_ids )
617 ( element_export_ids.size() != element_export_ranks.size() ) )
618 throw std::runtime_error(
619 "Cabana::CommunicationPlan::createSteering: Export ids and "
620 "ranks different sizes!" );
629 std::vector<std::size_t> offsets( num_n, 0.0 );
630 for (
int n = 1; n < num_n; ++n )
634 Kokkos::View<std::size_t*, Kokkos::HostSpace> rank_offsets_host(
635 Kokkos::ViewAllocateWithoutInitializing(
"rank_map" ), comm_size );
636 for (
int n = 0; n < num_n; ++n )
637 rank_offsets_host(
_neighbors[n] ) = offsets[n];
638 auto rank_offsets = Kokkos::create_mirror_view_and_copy(
645 Kokkos::ViewAllocateWithoutInitializing(
"export_steering" ),
648 Kokkos::parallel_for(
649 "Cabana::CommunicationPlan::createSteering",
651 KOKKOS_LAMBDA(
const int i ) {
652 if ( element_export_ranks( i ) >= 0 )
653 steer_vec( rank_offsets( element_export_ranks( i ) ) +
654 neighbor_ids( i ) ) =
655 ( use_iota ) ? i : element_export_ids( i );
660 template <
class PackViewType,
class RankViewType,
class IdViewType>
661 void createSteering(
const bool use_iota,
const PackViewType& neighbor_ids,
662 const RankViewType& element_export_ranks,
663 const IdViewType& element_export_ids )
667 element_export_ranks, element_export_ids );
696template <
class AoSoAType>
706 using data_type =
typename particle_data_type::tuple_type;
708 using buffer_type =
typename Kokkos::View<data_type*, memory_space>;
718 Kokkos::ViewAllocateWithoutInitializing(
"send_buffer" ), 0 );
720 Kokkos::ViewAllocateWithoutInitializing(
"recv_buffer" ), 0 );
747template <
class SliceType>
757 using data_type =
typename particle_data_type::value_type;
760 typename Kokkos::View<data_type**, Kokkos::LayoutRight, memory_space>;
772 Kokkos::ViewAllocateWithoutInitializing(
"send_buffer" ), 0, 0 );
774 Kokkos::ViewAllocateWithoutInitializing(
"recv_buffer" ), 0, 0 );
792 for ( std::size_t d = 2; d <
_particles.viewRank(); ++d )
810template <
class CommPlanType,
class CommDataType>
840 const double overallocation = 1.0 )
886 if ( use_overallocation )
891 _comm_data.reallocateSend( shrunk_send_size );
892 _comm_data.reallocateReceive( shrunk_recv_size );
899 template <
class ExecutionSpace>
903 void reserveImpl(
const CommPlanType& comm_plan,
905 const std::size_t total_send,
906 const std::size_t total_recv,
907 const double overallocation )
909 if ( overallocation < 1.0 )
910 throw std::runtime_error(
"Cabana::CommunicationPlan: "
911 "Cannot allocate buffers with less space "
912 "than data to communicate!" );
915 reserveImpl( comm_plan, particles, total_send, total_recv );
917 void reserveImpl(
const CommPlanType& comm_plan,
919 const std::size_t total_send,
920 const std::size_t total_recv )
926 auto new_send_size =
static_cast<std::size_t
>(
928 if ( new_send_size > send_capacity )
932 auto new_recv_size =
static_cast<std::size_t
>(
934 if ( new_recv_size > recv_capacity )
935 _comm_data.reallocateReceive( new_recv_size );
959template <
class MemorySpace,
class CommSpaceType = Mpi>
963template <
class CommPlanType,
class CommDataType,
class CommSpaceType = Mpi>
969#ifdef Cabana_ENABLE_MPI
Multi-node communication patterns. Uses vanilla MPI as the communication backend.
auto sendSize()
Current send buffer size.
Definition Cabana_CommunicationPlanBase.hpp:866
std::size_t _send_size
Send sizes.
Definition Cabana_CommunicationPlanBase.hpp:953
Kokkos::RangePolicy< execution_space > policy_type
Kokkos execution policy.
Definition Cabana_CommunicationPlanBase.hpp:819
auto receiveSize()
Current receive buffer size.
Definition Cabana_CommunicationPlanBase.hpp:872
typename comm_data_type::memory_space memory_space
Kokkos memory space.
Definition Cabana_CommunicationPlanBase.hpp:825
auto receiveCapacity()
Current allocated receive buffer space.
Definition Cabana_CommunicationPlanBase.hpp:876
comm_data_type _comm_data
Communication plan.
Definition Cabana_CommunicationPlanBase.hpp:949
void setData(const particle_data_type &particles)
Update particles to communicate.
Definition Cabana_CommunicationPlanBase.hpp:856
plan_type _comm_plan
Communication plan.
Definition Cabana_CommunicationPlanBase.hpp:947
particle_data_type getData() const
Get the particles to communicate.
Definition Cabana_CommunicationPlanBase.hpp:854
typename comm_data_type::data_type data_type
Communication data type.
Definition Cabana_CommunicationPlanBase.hpp:827
CommPlanType plan_type
Communication plan type (Halo, Distributor)
Definition Cabana_CommunicationPlanBase.hpp:815
buffer_type getSendBuffer() const
Get the communication send buffer.
Definition Cabana_CommunicationPlanBase.hpp:849
CommDataType comm_data_type
Communication data type.
Definition Cabana_CommunicationPlanBase.hpp:821
typename plan_type::execution_space execution_space
Kokkos execution space.
Definition Cabana_CommunicationPlanBase.hpp:817
typename comm_data_type::particle_data_type particle_data_type
Particle data type.
Definition Cabana_CommunicationPlanBase.hpp:823
void apply(ExecutionSpace)
Perform the communication (migrate, gather, scatter).
double _overallocation
Overallocation factor.
Definition Cabana_CommunicationPlanBase.hpp:951
CommunicationDataBase(const CommPlanType &comm_plan, const particle_data_type &particles, const double overallocation=1.0)
Definition Cabana_CommunicationPlanBase.hpp:838
void shrinkToFit(const bool use_overallocation=false)
Reduce communication buffers to current send/receive sizes.
Definition Cabana_CommunicationPlanBase.hpp:882
typename comm_data_type::buffer_type buffer_type
Communication buffer type.
Definition Cabana_CommunicationPlanBase.hpp:829
auto sendCapacity()
Current allocated send buffer space.
Definition Cabana_CommunicationPlanBase.hpp:874
std::size_t _recv_size
Receive sizes.
Definition Cabana_CommunicationPlanBase.hpp:955
virtual void apply()=0
Perform the communication (migrate, gather, scatter).
buffer_type getReceiveBuffer() const
Get the communication receive buffer.
Definition Cabana_CommunicationPlanBase.hpp:851
auto getSliceComponents()
Get the total number of components in the slice.
Definition Cabana_CommunicationPlanBase.hpp:944
Definition Cabana_CommunicationPlanBase.hpp:964
std::size_t exportSize() const
Get the number of export elements.
Definition Cabana_CommunicationPlanBase.hpp:534
std::size_t _total_num_export
Number of elements exported.
Definition Cabana_CommunicationPlanBase.hpp:677
std::vector< std::size_t > _num_export
Number of elements exported to each neighbor.
Definition Cabana_CommunicationPlanBase.hpp:681
Kokkos::View< std::size_t *, memory_space > _export_steering
Export steering vector.
Definition Cabana_CommunicationPlanBase.hpp:689
CommunicationPlanBase(MPI_Comm comm)
Constructor.
Definition Cabana_CommunicationPlanBase.hpp:442
int neighborRank(const int neighbor) const
Given a local neighbor id get its rank in the MPI communicator.
Definition Cabana_CommunicationPlanBase.hpp:479
std::size_t numImport(const int neighbor) const
Get the number of elements this rank will import from a given neighbor.
Definition Cabana_CommunicationPlanBase.hpp:512
std::size_t totalNumExport() const
Get the total number of exports this rank will do.
Definition Cabana_CommunicationPlanBase.hpp:502
std::size_t totalNumImport() const
Get the total number of imports this rank will do.
Definition Cabana_CommunicationPlanBase.hpp:522
std::shared_ptr< MPI_Comm > _comm_ptr
Shared pointer to Mpi communicator.
Definition Cabana_CommunicationPlanBase.hpp:673
void createExportSteering(const PackViewType &neighbor_ids, const RankViewType &element_export_ranks, const IdViewType &element_export_ids)
Create the export steering vector.
Definition Cabana_CommunicationPlanBase.hpp:597
std::size_t _total_num_import
Number of elements imported.
Definition Cabana_CommunicationPlanBase.hpp:679
void createExportSteering(const PackViewType &neighbor_ids, const RankViewType &element_export_ranks)
Create the export steering vector.
Definition Cabana_CommunicationPlanBase.hpp:569
typename memory_space::memory_space::size_type size_type
Size type.
Definition Cabana_CommunicationPlanBase.hpp:434
int numNeighbor() const
Get the number of neighbor ranks that this rank will communicate with.
Definition Cabana_CommunicationPlanBase.hpp:472
std::size_t _num_export_element
Definition Cabana_CommunicationPlanBase.hpp:687
std::vector< std::size_t > _num_import
Number of elements imported from each neighbor.
Definition Cabana_CommunicationPlanBase.hpp:683
std::vector< int > _neighbors
List of Mpi neighbors.
Definition Cabana_CommunicationPlanBase.hpp:675
std::size_t numExport(const int neighbor) const
Get the number of elements this rank will export to a given neighbor.
Definition Cabana_CommunicationPlanBase.hpp:492
Kokkos::View< std::size_t *, memory_space > getExportSteering() const
Get the steering vector for the exports.
Definition Cabana_CommunicationPlanBase.hpp:545
typename memory_space::execution_space execution_space
Default execution space.
Definition Cabana_CommunicationPlanBase.hpp:431
MemorySpace memory_space
Kokkos memory space.
Definition Cabana_CommunicationPlanBase.hpp:427
MPI_Comm comm() const
Get the MPI communicator.
Definition Cabana_CommunicationPlanBase.hpp:465
Definition Cabana_CommunicationPlanBase.hpp:960
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
std::vector< int > getUniqueTopology(MPI_Comm comm, std::vector< int > topology)
Return unique neighbor ranks, with the current rank first.
Definition Cabana_CommunicationPlanBase.hpp:374
buffer_type _send_buffer
Send buffer.
Definition Cabana_CommunicationPlanBase.hpp:735
void reallocateSend(const std::size_t num_send)
Resize the send buffer.
Definition Cabana_CommunicationPlanBase.hpp:724
typename particle_data_type::tuple_type data_type
Communication data type.
Definition Cabana_CommunicationPlanBase.hpp:706
AoSoAType particle_data_type
Particle data type.
Definition Cabana_CommunicationPlanBase.hpp:702
buffer_type _recv_buffer
Receive buffer.
Definition Cabana_CommunicationPlanBase.hpp:737
typename particle_data_type::memory_space memory_space
Kokkos memory space.
Definition Cabana_CommunicationPlanBase.hpp:704
particle_data_type _particles
Particle AoSoA.
Definition Cabana_CommunicationPlanBase.hpp:739
CommunicationDataAoSoA(particle_data_type particles)
Definition Cabana_CommunicationPlanBase.hpp:714
typename Kokkos::View< data_type *, memory_space > buffer_type
Communication buffer type.
Definition Cabana_CommunicationPlanBase.hpp:708
std::size_t _num_comp
Slice components.
Definition Cabana_CommunicationPlanBase.hpp:741
void reallocateReceive(const std::size_t num_recv)
Resize the receive buffer.
Definition Cabana_CommunicationPlanBase.hpp:729
void setSliceComponents()
Get the total number of components in the slice.
Definition Cabana_CommunicationPlanBase.hpp:789
void reallocateSend(const std::size_t num_send)
Resize the send buffer.
Definition Cabana_CommunicationPlanBase.hpp:778
buffer_type _send_buffer
Send buffer.
Definition Cabana_CommunicationPlanBase.hpp:797
buffer_type _recv_buffer
Receive buffer.
Definition Cabana_CommunicationPlanBase.hpp:799
void reallocateReceive(const std::size_t num_recv)
Resize the receive buffer.
Definition Cabana_CommunicationPlanBase.hpp:783
typename Kokkos::View< data_type **, Kokkos::LayoutRight, memory_space > buffer_type
Communication buffer type.
Definition Cabana_CommunicationPlanBase.hpp:759
typename particle_data_type::value_type data_type
Communication data type.
Definition Cabana_CommunicationPlanBase.hpp:757
typename particle_data_type::memory_space memory_space
Kokkos memory space.
Definition Cabana_CommunicationPlanBase.hpp:755
SliceType particle_data_type
Particle data type.
Definition Cabana_CommunicationPlanBase.hpp:753
std::size_t _num_comp
Slice components.
Definition Cabana_CommunicationPlanBase.hpp:803
particle_data_type _particles
Particle slice.
Definition Cabana_CommunicationPlanBase.hpp:801
CommunicationDataSlice(particle_data_type particles)
Definition Cabana_CommunicationPlanBase.hpp:766
Definition Cabana_Types.hpp:88
AoSoA static type checker.
Definition Cabana_AoSoA.hpp:61
Slice static type checker.
Definition Cabana_Slice.hpp:868