16#ifndef CABANA_DISTRIBUTOR_HPP
17#define CABANA_DISTRIBUTOR_HPP
23#include <Kokkos_Core.hpp>
24#include <Kokkos_Profiling_ScopedRegion.hpp>
61template <
class MemorySpace>
99 template <
class ViewType>
101 const std::vector<int>& neighbor_ranks )
105 element_export_ranks, neighbor_ranks );
136 template <
class ViewType>
148struct is_distributor_impl :
public std::false_type
152template <
typename MemorySpace>
153struct is_distributor_impl<
Distributor<MemorySpace>> :
public std::true_type
161 :
public is_distributor_impl<typename std::remove_cv<T>::type>::type
172template <
class ExecutionSpace,
class Distributor_t,
class AoSoA_t>
174 ExecutionSpace,
const Distributor_t& distributor,
const AoSoA_t& src,
180 Kokkos::Profiling::ScopedRegion region(
"Cabana::migrate" );
188 MPI_Comm_rank( distributor.comm(), &my_rank );
191 int num_n = distributor.numNeighbor();
197 std::size_t num_stay =
198 ( num_n > 0 && distributor.neighborRank( 0 ) == my_rank )
199 ? distributor.numExport( 0 )
203 std::size_t num_send = distributor.totalNumExport() - num_stay;
204 Kokkos::View<
typename AoSoA_t::tuple_type*,
205 typename Distributor_t::memory_space>
206 send_buffer( Kokkos::ViewAllocateWithoutInitializing(
207 "distributor_send_buffer" ),
211 Kokkos::View<
typename AoSoA_t::tuple_type*,
212 typename Distributor_t::memory_space>
213 recv_buffer( Kokkos::ViewAllocateWithoutInitializing(
214 "distributor_recv_buffer" ),
215 distributor.totalNumImport() );
218 auto steering = distributor.getExportSteering();
224 auto build_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
226 auto tpl = src.getTuple( steering( i ) );
228 recv_buffer( i ) = tpl;
230 send_buffer( i - num_stay ) = tpl;
232 Kokkos::RangePolicy<ExecutionSpace> build_send_buffer_policy(
233 0, distributor.totalNumExport() );
234 Kokkos::parallel_for(
"Cabana::Impl::distributeData::build_send_buffer",
235 build_send_buffer_policy, build_send_buffer_func );
239 const int mpi_tag = 1234;
242 std::vector<MPI_Request> requests;
243 requests.reserve( num_n );
244 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
245 for (
int n = 0; n < num_n; ++n )
247 recv_range.second = recv_range.first + distributor.numImport( n );
249 if ( ( distributor.numImport( n ) > 0 ) &&
250 ( distributor.neighborRank( n ) != my_rank ) )
252 auto recv_subview = Kokkos::subview( recv_buffer, recv_range );
254 requests.push_back( MPI_Request() );
256 MPI_Irecv( recv_subview.data(),
257 recv_subview.size() *
258 sizeof(
typename AoSoA_t::tuple_type ),
259 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
260 distributor.comm(), &( requests.back() ) );
263 recv_range.first = recv_range.second;
267 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
268 for (
int n = 0; n < num_n; ++n )
270 if ( ( distributor.numExport( n ) > 0 ) &&
271 ( distributor.neighborRank( n ) != my_rank ) )
273 send_range.second = send_range.first + distributor.numExport( n );
275 auto send_subview = Kokkos::subview( send_buffer, send_range );
277 MPI_Send( send_subview.data(),
278 send_subview.size() *
279 sizeof(
typename AoSoA_t::tuple_type ),
280 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
281 distributor.comm() );
283 send_range.first = send_range.second;
288 std::vector<MPI_Status> status( requests.size() );
290 MPI_Waitall( requests.size(), requests.data(), status.data() );
291 if ( MPI_SUCCESS != ec )
292 throw std::logic_error(
293 "Cabana::Distributor: Failed MPI Communication" );
296 auto extract_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
298 dst.setTuple( i, recv_buffer( i ) );
300 Kokkos::RangePolicy<ExecutionSpace> extract_recv_buffer_policy(
301 0, distributor.totalNumImport() );
302 Kokkos::parallel_for(
"Cabana::Impl::distributeData::extract_recv_buffer",
303 extract_recv_buffer_policy,
304 extract_recv_buffer_func );
308 MPI_Barrier( distributor.comm() );
335template <
class ExecutionSpace,
class Distributor_t,
class AoSoA_t>
336void migrate( ExecutionSpace exec_space,
const Distributor_t& distributor,
337 const AoSoA_t& src, AoSoA_t& dst,
343 if ( src.size() != distributor.exportSize() )
344 throw std::runtime_error(
"Cabana::migrate: Source is "
345 "the wrong size for migration! (Label: " +
347 if ( dst.size() != distributor.totalNumImport() )
348 throw std::runtime_error(
"Cabana::migrate: Destination "
349 "is the wrong size for migration! (Label: " +
353 Impl::distributeData( exec_space, distributor, src, dst );
373template <
class Distributor_t,
class AoSoA_t>
374void migrate(
const Distributor_t& distributor,
const AoSoA_t& src,
380 migrate(
typename Distributor_t::execution_space{}, distributor, src, dst );
407template <
class ExecutionSpace,
class Distributor_t,
class AoSoA_t>
408void migrate( ExecutionSpace exec_space,
const Distributor_t& distributor,
415 if ( aosoa.size() != distributor.exportSize() )
416 throw std::runtime_error(
417 "Cabana::migrate (in-place): "
418 "AoSoA is the wrong size for migration! (Label: " +
419 aosoa.label() +
")" );
424 ( distributor.totalNumImport() > distributor.exportSize() );
429 aosoa.resize( distributor.totalNumImport() );
432 Impl::distributeData( exec_space, distributor, aosoa, aosoa );
436 if ( !dst_is_bigger )
437 aosoa.resize( distributor.totalNumImport() );
461template <
class Distributor_t,
class AoSoA_t>
462void migrate(
const Distributor_t& distributor, AoSoA_t& aosoa,
467 migrate(
typename Distributor_t::execution_space{}, distributor, aosoa );
491template <
class ExecutionSpace,
class Distributor_t,
class Slice_t>
492void migrate( ExecutionSpace,
const Distributor_t& distributor,
493 const Slice_t& src, Slice_t& dst,
499 if ( src.size() != distributor.exportSize() )
500 throw std::runtime_error(
"Cabana::migrate: Source Slice is the wrong "
501 "size for migration! (Label: " +
503 if ( dst.size() != distributor.totalNumImport() )
504 throw std::runtime_error(
"Cabana::migrate: Destination Slice is the "
505 "wrong size for migration! (Label: " +
510 for (
size_t d = 2; d < src.viewRank(); ++d )
511 num_comp *= src.extent( d );
514 auto src_data = src.data();
515 auto dst_data = dst.data();
519 MPI_Comm_rank( distributor.comm(), &my_rank );
522 int num_n = distributor.numNeighbor();
528 std::size_t num_stay =
529 ( num_n > 0 && distributor.neighborRank( 0 ) == my_rank )
530 ? distributor.numExport( 0 )
535 std::size_t num_send = distributor.totalNumExport() - num_stay;
536 Kokkos::View<
typename Slice_t::value_type**, Kokkos::LayoutRight,
537 typename Distributor_t::memory_space>
538 send_buffer( Kokkos::ViewAllocateWithoutInitializing(
539 "distributor_send_buffer" ),
540 num_send, num_comp );
544 Kokkos::View<
typename Slice_t::value_type**, Kokkos::LayoutRight,
545 typename Distributor_t::memory_space>
546 recv_buffer( Kokkos::ViewAllocateWithoutInitializing(
547 "distributor_recv_buffer" ),
548 distributor.totalNumImport(), num_comp );
551 auto steering = distributor.getExportSteering();
556 auto build_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
558 auto s_src = Slice_t::index_type::s( steering( i ) );
559 auto a_src = Slice_t::index_type::a( steering( i ) );
560 std::size_t src_offset = s_src * src.stride( 0 ) + a_src;
562 for ( std::size_t n = 0; n < num_comp; ++n )
563 recv_buffer( i, n ) =
564 src_data[src_offset + n * Slice_t::vector_length];
566 for ( std::size_t n = 0; n < num_comp; ++n )
567 send_buffer( i - num_stay, n ) =
568 src_data[src_offset + n * Slice_t::vector_length];
570 Kokkos::RangePolicy<ExecutionSpace> build_send_buffer_policy(
571 0, distributor.totalNumExport() );
572 Kokkos::parallel_for(
"Cabana::migrate::build_send_buffer",
573 build_send_buffer_policy, build_send_buffer_func );
577 const int mpi_tag = 1234;
580 std::vector<MPI_Request> requests;
581 requests.reserve( num_n );
582 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
583 for (
int n = 0; n < num_n; ++n )
585 recv_range.second = recv_range.first + distributor.numImport( n );
587 if ( ( distributor.numImport( n ) > 0 ) &&
588 ( distributor.neighborRank( n ) != my_rank ) )
591 Kokkos::subview( recv_buffer, recv_range, Kokkos::ALL );
593 requests.push_back( MPI_Request() );
595 MPI_Irecv( recv_subview.data(),
596 recv_subview.size() *
597 sizeof(
typename Slice_t::value_type ),
598 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
599 distributor.comm(), &( requests.back() ) );
602 recv_range.first = recv_range.second;
606 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
607 for (
int n = 0; n < num_n; ++n )
609 if ( ( distributor.numExport( n ) > 0 ) &&
610 ( distributor.neighborRank( n ) != my_rank ) )
612 send_range.second = send_range.first + distributor.numExport( n );
615 Kokkos::subview( send_buffer, send_range, Kokkos::ALL );
617 MPI_Send( send_subview.data(),
618 send_subview.size() *
619 sizeof(
typename Slice_t::value_type ),
620 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
621 distributor.comm() );
623 send_range.first = send_range.second;
628 std::vector<MPI_Status> status( requests.size() );
630 MPI_Waitall( requests.size(), requests.data(), status.data() );
631 if ( MPI_SUCCESS != ec )
632 throw std::logic_error(
"Cabana::migrate: Failed MPI Communication" );
635 auto extract_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
637 auto s = Slice_t::index_type::s( i );
638 auto a = Slice_t::index_type::a( i );
639 std::size_t dst_offset = s * dst.stride( 0 ) + a;
640 for ( std::size_t n = 0; n < num_comp; ++n )
641 dst_data[dst_offset + n * Slice_t::vector_length] =
644 Kokkos::RangePolicy<ExecutionSpace> extract_recv_buffer_policy(
645 0, distributor.totalNumImport() );
646 Kokkos::parallel_for(
"Cabana::migrate::extract_recv_buffer",
647 extract_recv_buffer_policy,
648 extract_recv_buffer_func );
652 MPI_Barrier( distributor.comm() );
674template <
class Distributor_t,
class Slice_t>
675void migrate(
const Distributor_t& distributor,
const Slice_t& src,
681 migrate(
typename Distributor_t::execution_space{}, distributor, src, dst );
Array-of-Struct-of-Arrays particle data structure.
Multi-node communication patterns.
Slice a single particle property from an AoSoA.
Kokkos::View< size_type *, memory_space > createFromExportsOnly(ExecutionSpace exec_space, const ViewType &element_export_ranks)
Export rank creator. Use this when you don't know who you will receiving from - only who you are send...
Definition Cabana_CommunicationPlan.hpp:758
void createExportSteering(const PackViewType &neighbor_ids, const RankViewType &element_export_ranks)
Create the export steering vector.
Definition Cabana_CommunicationPlan.hpp:946
MPI_Comm comm() const
Get the MPI communicator.
Definition Cabana_CommunicationPlan.hpp:468
Kokkos::View< size_type *, memory_space > createFromExportsAndTopology(ExecutionSpace exec_space, const ViewType &element_export_ranks, const std::vector< int > &neighbor_ranks)
Neighbor and export rank creator. Use this when you already know which ranks neighbor each other (i....
Definition Cabana_CommunicationPlan.hpp:596
CommunicationPlan(MPI_Comm comm)
Constructor.
Definition Cabana_CommunicationPlan.hpp:446
A communication plan for migrating data from one uniquely-owned decomposition to another uniquely own...
Definition Cabana_Distributor.hpp:63
Distributor(MPI_Comm comm, const ViewType &element_export_ranks)
Export rank constructor. Use this when you don't know who you will be receiving from - only who you a...
Definition Cabana_Distributor.hpp:137
Distributor(MPI_Comm comm, const ViewType &element_export_ranks, const std::vector< int > &neighbor_ranks)
Topology and export rank constructor. Use this when you already know which ranks neighbor each other ...
Definition Cabana_Distributor.hpp:100
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
void migrate(ExecutionSpace exec_space, const Distributor_t &distributor, const AoSoA_t &src, AoSoA_t &dst, typename std::enable_if<(is_distributor< Distributor_t >::value &&is_aosoa< AoSoA_t >::value), int >::type *=0)
Synchronously migrate data between two different decompositions using the distributor forward communi...
Definition Cabana_Distributor.hpp:336
Definition Cabana_Types.hpp:88
AoSoA static type checker.
Definition Cabana_AoSoA.hpp:61
Distributor static type checker.
Definition Cabana_Distributor.hpp:162
Slice static type checker.
Definition Cabana_Slice.hpp:861