16#ifndef CABANA_MIGRATE_MPI_HPP
17#define CABANA_MIGRATE_MPI_HPP
22#include <Kokkos_Core.hpp>
23#include <Kokkos_Profiling_ScopedRegion.hpp>
40template <
class ExecutionSpace,
class Distributor_t,
class AoSoA_t>
42 Mpi, ExecutionSpace,
const Distributor_t& distributor,
const AoSoA_t& src,
44 typename std::enable_if<( ( is_distributor<Distributor_t>::value ) &&
45 is_aosoa<AoSoA_t>::value ),
48 Kokkos::Profiling::ScopedRegion region(
"Cabana::migrate" );
50 static_assert( is_accessible_from<
typename Distributor_t::memory_space,
56 MPI_Comm_rank( distributor.comm(), &my_rank );
59 int num_n = distributor.numNeighbor();
65 std::size_t num_stay =
66 ( num_n > 0 && distributor.neighborRank( 0 ) == my_rank )
67 ? distributor.numExport( 0 )
71 std::size_t num_send = distributor.totalNumExport() - num_stay;
72 Kokkos::View<
typename AoSoA_t::tuple_type*,
73 typename Distributor_t::memory_space>
74 send_buffer( Kokkos::ViewAllocateWithoutInitializing(
75 "distributor_send_buffer" ),
79 Kokkos::View<
typename AoSoA_t::tuple_type*,
80 typename Distributor_t::memory_space>
81 recv_buffer( Kokkos::ViewAllocateWithoutInitializing(
82 "distributor_recv_buffer" ),
83 distributor.totalNumImport() );
86 auto steering = distributor.getExportSteering();
92 auto build_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
94 auto tpl = src.getTuple( steering( i ) );
96 recv_buffer( i ) = tpl;
98 send_buffer( i - num_stay ) = tpl;
100 Kokkos::RangePolicy<ExecutionSpace> build_send_buffer_policy(
101 0, distributor.totalNumExport() );
102 Kokkos::parallel_for(
"Cabana::Impl::distributeData::build_send_buffer",
103 build_send_buffer_policy, build_send_buffer_func );
107 const int mpi_tag = 1234;
110 std::vector<MPI_Request> requests;
111 requests.reserve( num_n );
112 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
113 for (
int n = 0; n < num_n; ++n )
115 recv_range.second = recv_range.first + distributor.numImport( n );
117 if ( ( distributor.numImport( n ) > 0 ) &&
118 ( distributor.neighborRank( n ) != my_rank ) )
120 auto recv_subview = Kokkos::subview( recv_buffer, recv_range );
122 requests.push_back( MPI_Request() );
124 MPI_Irecv( recv_subview.data(),
125 recv_subview.size() *
126 sizeof(
typename AoSoA_t::tuple_type ),
127 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
128 distributor.comm(), &( requests.back() ) );
131 recv_range.first = recv_range.second;
135 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
136 for (
int n = 0; n < num_n; ++n )
138 if ( ( distributor.numExport( n ) > 0 ) &&
139 ( distributor.neighborRank( n ) != my_rank ) )
141 send_range.second = send_range.first + distributor.numExport( n );
143 auto send_subview = Kokkos::subview( send_buffer, send_range );
145 MPI_Send( send_subview.data(),
146 send_subview.size() *
147 sizeof(
typename AoSoA_t::tuple_type ),
148 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
149 distributor.comm() );
151 send_range.first = send_range.second;
156 std::vector<MPI_Status> status( requests.size() );
158 MPI_Waitall( requests.size(), requests.data(), status.data() );
159 if ( MPI_SUCCESS != ec )
160 throw std::logic_error(
161 "Cabana::Distributor: Failed MPI Communication" );
164 auto extract_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
166 dst.setTuple( i, recv_buffer( i ) );
168 Kokkos::RangePolicy<ExecutionSpace> extract_recv_buffer_policy(
169 0, distributor.totalNumImport() );
170 Kokkos::parallel_for(
"Cabana::Impl::distributeData::extract_recv_buffer",
171 extract_recv_buffer_policy,
172 extract_recv_buffer_func );
176 MPI_Barrier( distributor.comm() );
200template <
class ExecutionSpace,
class Distributor_t,
class Slice_t>
202 Mpi, ExecutionSpace,
const Distributor_t& distributor,
const Slice_t& src,
204 typename std::enable_if<( ( is_distributor<Distributor_t>::value ) &&
205 is_slice<Slice_t>::value ),
210 for (
size_t d = 2; d < src.viewRank(); ++d )
211 num_comp *= src.extent( d );
214 auto src_data = src.data();
215 auto dst_data = dst.data();
219 MPI_Comm_rank( distributor.comm(), &my_rank );
222 int num_n = distributor.numNeighbor();
228 std::size_t num_stay =
229 ( num_n > 0 && distributor.neighborRank( 0 ) == my_rank )
230 ? distributor.numExport( 0 )
235 std::size_t num_send = distributor.totalNumExport() - num_stay;
236 Kokkos::View<
typename Slice_t::value_type**, Kokkos::LayoutRight,
237 typename Distributor_t::memory_space>
238 send_buffer( Kokkos::ViewAllocateWithoutInitializing(
239 "distributor_send_buffer" ),
240 num_send, num_comp );
244 Kokkos::View<
typename Slice_t::value_type**, Kokkos::LayoutRight,
245 typename Distributor_t::memory_space>
246 recv_buffer( Kokkos::ViewAllocateWithoutInitializing(
247 "distributor_recv_buffer" ),
248 distributor.totalNumImport(), num_comp );
251 auto steering = distributor.getExportSteering();
256 auto build_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
258 auto s_src = Slice_t::index_type::s( steering( i ) );
259 auto a_src = Slice_t::index_type::a( steering( i ) );
260 std::size_t src_offset = s_src * src.stride( 0 ) + a_src;
262 for ( std::size_t n = 0; n < num_comp; ++n )
263 recv_buffer( i, n ) =
264 src_data[src_offset + n * Slice_t::vector_length];
266 for ( std::size_t n = 0; n < num_comp; ++n )
267 send_buffer( i - num_stay, n ) =
268 src_data[src_offset + n * Slice_t::vector_length];
270 Kokkos::RangePolicy<ExecutionSpace> build_send_buffer_policy(
271 0, distributor.totalNumExport() );
272 Kokkos::parallel_for(
"Cabana::migrate::build_send_buffer",
273 build_send_buffer_policy, build_send_buffer_func );
277 const int mpi_tag = 1234;
280 std::vector<MPI_Request> requests;
281 requests.reserve( num_n );
282 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
283 for (
int n = 0; n < num_n; ++n )
285 recv_range.second = recv_range.first + distributor.numImport( n );
287 if ( ( distributor.numImport( n ) > 0 ) &&
288 ( distributor.neighborRank( n ) != my_rank ) )
291 Kokkos::subview( recv_buffer, recv_range, Kokkos::ALL );
293 requests.push_back( MPI_Request() );
295 MPI_Irecv( recv_subview.data(),
296 recv_subview.size() *
297 sizeof(
typename Slice_t::value_type ),
298 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
299 distributor.comm(), &( requests.back() ) );
302 recv_range.first = recv_range.second;
306 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
307 for (
int n = 0; n < num_n; ++n )
309 if ( ( distributor.numExport( n ) > 0 ) &&
310 ( distributor.neighborRank( n ) != my_rank ) )
312 send_range.second = send_range.first + distributor.numExport( n );
315 Kokkos::subview( send_buffer, send_range, Kokkos::ALL );
317 MPI_Send( send_subview.data(),
318 send_subview.size() *
319 sizeof(
typename Slice_t::value_type ),
320 MPI_BYTE, distributor.neighborRank( n ), mpi_tag,
321 distributor.comm() );
323 send_range.first = send_range.second;
328 std::vector<MPI_Status> status( requests.size() );
330 MPI_Waitall( requests.size(), requests.data(), status.data() );
331 if ( MPI_SUCCESS != ec )
332 throw std::logic_error(
"Cabana::migrate: Failed MPI Communication" );
335 auto extract_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
337 auto s = Slice_t::index_type::s( i );
338 auto a = Slice_t::index_type::a( i );
339 std::size_t dst_offset = s * dst.stride( 0 ) + a;
340 for ( std::size_t n = 0; n < num_comp; ++n )
341 dst_data[dst_offset + n * Slice_t::vector_length] =
344 Kokkos::RangePolicy<ExecutionSpace> extract_recv_buffer_policy(
345 0, distributor.totalNumImport() );
346 Kokkos::parallel_for(
"Cabana::migrate::extract_recv_buffer",
347 extract_recv_buffer_policy,
348 extract_recv_buffer_func );
352 MPI_Barrier( distributor.comm() );
Array-of-Struct-of-Arrays particle data structure.
Slice a single particle property from an AoSoA.
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36