108 const RankViewType& element_export_ranks,
109 const std::vector<int>& neighbor_ranks )
122 MPI_Comm_size( this->
comm(), &comm_size );
126 MPI_Comm_rank( this->
comm(), &my_rank );
130 const int mpi_tag = 1221;
138 auto counts_and_ids = Impl::countSendsAndCreateSteering(
139 exec_space, element_export_ranks, comm_size,
140 typename Impl::CountSendsAndCreateSteeringAlgorithm<
141 ExecutionSpace>::type() );
144 auto neighbor_counts_host = Kokkos::create_mirror_view_and_copy(
145 Kokkos::HostSpace(), counts_and_ids.first );
148 for (
int n = 0; n < num_n; ++n )
152 std::vector<MPI_Request> requests;
153 requests.reserve( num_n );
154 for (
int n = 0; n < num_n; ++n )
157 requests.push_back( MPI_Request() );
158 MPI_Irecv( &this->
_num_import[n], 1, MPI_UNSIGNED_LONG,
160 &( requests.back() ) );
166 for (
int n = 0; n < num_n; ++n )
168 MPI_Send( &this->
_num_export[n], 1, MPI_UNSIGNED_LONG,
172 std::vector<MPI_Status> status( requests.size() );
174 MPI_Waitall( requests.size(), requests.data(), status.data() );
175 if ( MPI_SUCCESS != ec )
176 throw std::logic_error(
177 "Cabana::CommunicationPlan::createFromExportsAndTopology: "
178 "Failed MPI Communication" );
182 std::accumulate( this->
_num_export.begin(), this->_num_export.end(),
185 std::accumulate( this->
_num_import.begin(), this->_num_import.end(),
192 return counts_and_ids.second;
273 const RankViewType& element_export_ranks )
282 MPI_Comm_size( this->
comm(), &comm_size );
286 MPI_Comm_rank( this->
comm(), &my_rank );
290 const int mpi_tag = 1221;
294 auto counts_and_ids = Impl::countSendsAndCreateSteering(
295 exec_space, element_export_ranks, comm_size,
296 typename Impl::CountSendsAndCreateSteeringAlgorithm<
297 ExecutionSpace>::type() );
300 auto neighbor_counts_host = Kokkos::create_mirror_view_and_copy(
301 Kokkos::HostSpace(), counts_and_ids.first );
308 for (
int r = 0; r < comm_size; ++r )
309 if ( neighbor_counts_host( r ) > 0 )
312 this->
_num_export.push_back( neighbor_counts_host( r ) );
314 neighbor_counts_host( r ) = 1;
319 int num_export_rank = this->
_neighbors.size();
324 bool self_send =
false;
325 for (
int n = 0; n < num_export_rank; ++n )
336 int num_import_rank = -1;
337 std::vector<int> recv_counts( comm_size, 1 );
338 MPI_Reduce_scatter( neighbor_counts_host.data(), &num_import_rank,
339 recv_counts.data(), MPI_INT, MPI_SUM,
346 std::vector<std::size_t> import_sizes( num_import_rank );
347 std::vector<MPI_Request> requests( num_import_rank );
348 for (
int n = 0; n < num_import_rank; ++n )
349 MPI_Irecv( &import_sizes[n], 1, MPI_UNSIGNED_LONG, MPI_ANY_SOURCE,
350 mpi_tag, this->
comm(), &requests[n] );
353 int self_offset = ( self_send ) ? 1 : 0;
354 for (
int n = self_offset; n < num_export_rank; ++n )
355 MPI_Send( &this->
_num_export[n], 1, MPI_UNSIGNED_LONG,
359 std::vector<MPI_Status> status( requests.size() );
361 MPI_Waitall( requests.size(), requests.data(), status.data() );
362 if ( MPI_SUCCESS != ec )
363 throw std::logic_error(
364 "Cabana::CommunicationPlan::createFromExportsOnly: Failed MPI "
369 std::accumulate( import_sizes.begin(), import_sizes.end(),
370 ( self_send ) ? this->_num_import[0] : 0 );
374 for (
int i = 0; i < num_import_rank; ++i )
377 const auto source = status[i].MPI_SOURCE;
381 auto found_neighbor = std::find( this->
_neighbors.begin(),
382 this->_neighbors.end(), source );
386 if ( found_neighbor == std::end( this->
_neighbors ) )
399 std::distance( this->
_neighbors.begin(), found_neighbor );
406 MPI_Barrier( this->
comm() );
409 return counts_and_ids.second;
490 const RankViewType& element_import_ranks,
491 const IdViewType& element_import_ids,
492 const std::vector<int>& neighbor_ranks )
493 -> std::tuple<Kokkos::View<
typename RankViewType::size_type*,
494 typename RankViewType::memory_space>,
495 Kokkos::View<int*, typename RankViewType::memory_space>,
496 Kokkos::View<int*, typename IdViewType::memory_space>>
500 if ( element_import_ids.size() != element_import_ranks.size() )
501 throw std::runtime_error(
"Export ids and ranks different sizes!" );
509 MPI_Comm_size( this->
comm(), &comm_size );
513 MPI_Comm_rank( this->
comm(), &my_rank );
517 const int mpi_tag = 1221;
525 auto counts_and_ids = Impl::countSendsAndCreateSteering(
526 exec_space, element_import_ranks, comm_size,
527 typename Impl::CountSendsAndCreateSteeringAlgorithm<
528 ExecutionSpace>::type() );
531 auto neighbor_counts_host = Kokkos::create_mirror_view_and_copy(
532 Kokkos::HostSpace(), counts_and_ids.first );
535 for ( std::size_t n = 0; n < num_n; ++n )
541 std::vector<MPI_Request> requests;
542 requests.reserve( num_n * 2 );
543 for ( std::size_t n = 0; n < num_n; ++n )
546 requests.push_back( MPI_Request() );
547 MPI_Irecv( &this->
_num_export[n], 1, MPI_UNSIGNED_LONG,
549 &( requests.back() ) );
557 for ( std::size_t n = 0; n < num_n; ++n )
560 requests.push_back( MPI_Request() );
561 MPI_Isend( &this->
_num_import[n], 1, MPI_UNSIGNED_LONG,
563 &( requests.back() ) );
567 std::vector<MPI_Status> status( requests.size() );
569 MPI_Waitall( requests.size(), requests.data(), status.data() );
570 if ( MPI_SUCCESS != ec )
571 throw std::logic_error(
"Failed MPI Communication" );
575 this->_num_export.end(), 0 );
577 this->_num_import.end(), 0 );
582 Kokkos::View<int*, memory_space> export_indices(
587 std::vector<MPI_Request> mpi_requests( num_messages );
588 std::vector<MPI_Status> mpi_statuses( num_messages );
593 for ( std::size_t i = 0; i < num_n; i++ )
595 for ( std::size_t j = 0; j < this->
_num_export[i]; j++ )
597 MPI_Irecv( export_indices.data() + idx, 1, MPI_INT,
598 this->_neighbors[i], mpi_tag + 1, this->comm(),
599 &mpi_requests[idx] );
605 for ( std::size_t i = 0; i < element_import_ranks.extent( 0 ); i++ )
607 MPI_Isend( element_import_ids.data() + i, 1, MPI_INT,
608 *( element_import_ranks.data() + i ), mpi_tag + 1,
609 this->comm(), &mpi_requests[idx++] );
613 const int ec1 = MPI_Waitall( num_messages, mpi_requests.data(),
614 mpi_statuses.data() );
615 if ( MPI_SUCCESS != ec1 )
616 throw std::logic_error(
"Failed MPI Communication" );
621 Kokkos::View<int*, Kokkos::HostSpace> element_export_ranks_h(
625 element_export_ranks_h[i] = mpi_statuses[i].MPI_SOURCE;
627 auto element_export_ranks = Kokkos::create_mirror_view_and_copy(
630 auto counts_and_ids2 = Impl::countSendsAndCreateSteering(
631 exec_space, element_export_ranks, comm_size,
632 typename Impl::CountSendsAndCreateSteeringAlgorithm<
633 ExecutionSpace>::type() );
639 return std::tuple{ counts_and_ids2.second, element_export_ranks,
725 const RankViewType& element_import_ranks,
726 const IdViewType& element_import_ids )
727 -> std::tuple<Kokkos::View<
typename RankViewType::size_type*,
728 typename RankViewType::memory_space>,
729 Kokkos::View<int*, typename RankViewType::memory_space>,
730 Kokkos::View<int*, typename IdViewType::memory_space>>
734 if ( element_import_ids.size() != element_import_ranks.size() )
735 throw std::runtime_error(
"Export ids and ranks different sizes!" );
739 MPI_Comm_size( this->
comm(), &comm_size );
743 MPI_Comm_rank( this->
comm(), &rank );
747 const int mpi_tag = 1221;
750 Kokkos::View<int*, memory_space> importing_ranks(
"importing_ranks",
752 Kokkos::deep_copy( importing_ranks, 0 );
753 Kokkos::parallel_for(
754 "Cabana::storeImportRanks",
755 Kokkos::RangePolicy<ExecutionSpace>(
756 0, element_import_ranks.extent( 0 ) ),
757 KOKKOS_LAMBDA(
const int i ) {
758 int import_rank = element_import_ranks( i );
759 Kokkos::atomic_store( &importing_ranks( import_rank ), 1 );
762 auto importing_ranks_h = Kokkos::create_mirror_view_and_copy(
763 Kokkos::HostSpace(), importing_ranks );
766 Kokkos::View<int*, Kokkos::HostSpace> num_ranks_communicate(
767 "num_ranks_communicate", comm_size );
768 MPI_Allreduce( importing_ranks_h.data(), num_ranks_communicate.data(),
769 comm_size, MPI_INT, MPI_SUM, this->comm() );
773 int num_recvs = num_ranks_communicate( rank );
774 Kokkos::View<int*, Kokkos::HostSpace> send_counts(
"send_counts",
776 Kokkos::View<int*, Kokkos::HostSpace> send_to(
"send_to", num_recvs );
778 std::vector<MPI_Request> mpi_requests( num_recvs );
779 std::vector<MPI_Status> mpi_statuses( num_recvs );
782 for (
int i = 0; i < num_recvs; i++ )
784 MPI_Irecv( &send_counts( i ), 1, MPI_INT, MPI_ANY_SOURCE, mpi_tag,
785 this->
comm(), &mpi_requests[i] );
790 auto counts_and_ids = Impl::countSendsAndCreateSteering(
791 exec_space, element_import_ranks, comm_size,
792 typename Impl::CountSendsAndCreateSteeringAlgorithm<
793 ExecutionSpace>::type() );
796 auto neighbor_counts_host = Kokkos::create_mirror_view_and_copy(
797 Kokkos::HostSpace(), counts_and_ids.first );
804 for ( std::size_t i = 0; i < neighbor_counts_host.extent( 0 ); i++ )
806 if ( neighbor_counts_host( i ) != 0 )
809 MPI_Send( &neighbor_counts_host( i ), 1, MPI_INT, i, mpi_tag,
814 this->
_num_import.push_back( neighbor_counts_host( i ) );
822 MPI_Waitall( num_recvs, mpi_requests.data(), mpi_statuses.data() );
823 if ( MPI_SUCCESS != ec0 )
824 throw std::logic_error(
"Failed MPI Communication" );
829 for (
int i = 0; i < num_recvs; i++ )
831 send_to( i ) = mpi_statuses[i].MPI_SOURCE;
837 for (
int r = 0; r < num_recvs; ++r )
839 int export_to = send_to( r );
840 if ( export_to > -1 )
844 auto found_neighbor =
845 std::find( this->
_neighbors.begin(), this->_neighbors.end(),
850 if ( found_neighbor == std::end( this->
_neighbors ) )
862 auto n = std::distance( this->
_neighbors.begin(),
871 throw std::runtime_error(
872 "CommunicationPlan::createFromImportsOnly: "
873 "mpi_statuses[i].MPI_SOURCE returned a value >= -1" );
878 for ( std::size_t n = 0; n < this->
_neighbors.size(); ++n )
894 Kokkos::View<int*, memory_space> export_indices(
897 mpi_requests.clear();
898 mpi_statuses.clear();
901 mpi_requests.resize( num_messages );
902 mpi_statuses.resize( num_messages );
907 for (
int i = 0; i < num_recvs; i++ )
909 for (
int j = 0; j < send_counts( i ); j++ )
911 MPI_Irecv( export_indices.data() + idx, 1, MPI_INT,
912 send_to( i ), mpi_tag + 1, this->
comm(),
913 &mpi_requests[idx] );
919 for ( std::size_t i = 0; i < element_import_ranks.extent( 0 ); i++ )
921 MPI_Isend( element_import_ids.data() + i, 1, MPI_INT,
922 *( element_import_ranks.data() + i ), mpi_tag + 1,
923 this->comm(), &mpi_requests[idx++] );
927 const int ec1 = MPI_Waitall( num_messages, mpi_requests.data(),
928 mpi_statuses.data() );
929 if ( MPI_SUCCESS != ec1 )
930 throw std::logic_error(
"Failed MPI Communication" );
935 Kokkos::View<int*, Kokkos::HostSpace> element_export_ranks_h(
939 element_export_ranks_h[i] = mpi_statuses[i].MPI_SOURCE;
941 auto element_export_ranks = Kokkos::create_mirror_view_and_copy(
944 auto counts_and_ids2 = Impl::countSendsAndCreateSteering(
945 exec_space, element_export_ranks, comm_size,
946 typename Impl::CountSendsAndCreateSteeringAlgorithm<
947 ExecutionSpace>::type() );
951 MPI_Barrier( this->
comm() );
953 return std::tuple{ counts_and_ids2.second, element_export_ranks,