16#ifndef CABANA_HALO_MPI_HPP
17#define CABANA_HALO_MPI_HPP
23#include <Kokkos_Core.hpp>
24#include <Kokkos_Profiling_ScopedRegion.hpp>
41template <
class HaloType,
class AoSoAType>
42template <
class ExecutionSpace,
class CommSpaceType>
43std::enable_if_t<std::is_same<CommSpaceType, Mpi>::value,
void>
45 typename std::enable_if<is_aosoa<AoSoAType>::value>::type>::
46 applyImpl( ExecutionSpace, CommSpaceType )
48 Kokkos::Profiling::ScopedRegion region(
"Cabana::gather" );
51 auto send_buffer = this->getSendBuffer();
52 auto recv_buffer = this->getReceiveBuffer();
53 auto aosoa = this->getData();
56 auto steering = _comm_plan.getExportSteering();
58 auto gather_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
60 send_buffer( i ) = aosoa.getTuple( steering( i ) );
62 Kokkos::RangePolicy<ExecutionSpace> send_policy( 0, _send_size );
63 Kokkos::parallel_for(
"Cabana::gather::gather_send_buffer", send_policy,
64 gather_send_buffer_func );
68 const int mpi_tag = 2345;
71 int num_n = _comm_plan.numNeighbor();
72 std::vector<MPI_Request> requests( num_n );
73 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
74 for (
int n = 0; n < num_n; ++n )
76 recv_range.second = recv_range.first + _comm_plan.numImport( n );
78 auto recv_subview = Kokkos::subview( recv_buffer, recv_range );
80 MPI_Irecv( recv_subview.data(),
81 recv_subview.size() *
sizeof( data_type ), MPI_BYTE,
82 _comm_plan.neighborRank( n ), mpi_tag, _comm_plan.comm(),
85 recv_range.first = recv_range.second;
89 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
90 for (
int n = 0; n < num_n; ++n )
92 send_range.second = send_range.first + _comm_plan.numExport( n );
94 auto send_subview = Kokkos::subview( send_buffer, send_range );
96 MPI_Send( send_subview.data(),
97 send_subview.size() *
sizeof( data_type ), MPI_BYTE,
98 _comm_plan.neighborRank( n ), mpi_tag, _comm_plan.comm() );
100 send_range.first = send_range.second;
104 std::vector<MPI_Status> status( num_n );
106 MPI_Waitall( requests.size(), requests.data(), status.data() );
107 if ( MPI_SUCCESS != ec )
108 throw std::logic_error(
109 "Cabana::Gather::apply: Failed MPI Communication" );
112 std::size_t num_local = _comm_plan.numLocal();
113 auto extract_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
115 std::size_t ghost_idx = i + num_local;
116 aosoa.setTuple( ghost_idx, recv_buffer( i ) );
118 Kokkos::RangePolicy<ExecutionSpace> recv_policy( 0, _recv_size );
119 Kokkos::parallel_for(
"Cabana::gather::apply::extract_recv_buffer",
120 recv_policy, extract_recv_buffer_func );
124 MPI_Barrier( _comm_plan.comm() );
132template <
class HaloType,
class SliceType>
133template <
class ExecutionSpace,
class CommSpaceType>
134std::enable_if_t<std::is_same<CommSpaceType, Mpi>::value,
void>
135Gather<HaloType, SliceType,
136 typename std::enable_if<is_slice<SliceType>::value>::type>::
137 applyImpl( ExecutionSpace, CommSpaceType )
139 Kokkos::Profiling::ScopedRegion region(
"Cabana::gather" );
142 auto send_buffer = this->getSendBuffer();
143 auto recv_buffer = this->getReceiveBuffer();
144 auto slice = this->getData();
147 std::size_t num_comp = this->getSliceComponents();
150 auto slice_data =
slice.data();
153 auto steering = _comm_plan.getExportSteering();
156 auto gather_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
158 auto s = SliceType::index_type::s( steering( i ) );
159 auto a = SliceType::index_type::a( steering( i ) );
160 std::size_t slice_offset = s *
slice.stride( 0 ) + a;
161 for ( std::size_t n = 0; n < num_comp; ++n )
162 send_buffer( i, n ) =
163 slice_data[slice_offset + n * SliceType::vector_length];
165 Kokkos::RangePolicy<ExecutionSpace> send_policy( 0, _send_size );
166 Kokkos::parallel_for(
"Cabana::gather::gather_send_buffer", send_policy,
167 gather_send_buffer_func );
171 const int mpi_tag = 2345;
174 int num_n = _comm_plan.numNeighbor();
175 std::vector<MPI_Request> requests( num_n );
176 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
177 for (
int n = 0; n < num_n; ++n )
179 recv_range.second = recv_range.first + _comm_plan.numImport( n );
182 Kokkos::subview( recv_buffer, recv_range, Kokkos::ALL );
184 MPI_Irecv( recv_subview.data(),
185 recv_subview.size() *
sizeof( data_type ), MPI_BYTE,
186 _comm_plan.neighborRank( n ), mpi_tag, _comm_plan.comm(),
189 recv_range.first = recv_range.second;
193 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
194 for (
int n = 0; n < num_n; ++n )
196 send_range.second = send_range.first + _comm_plan.numExport( n );
199 Kokkos::subview( send_buffer, send_range, Kokkos::ALL );
201 MPI_Send( send_subview.data(),
202 send_subview.size() *
sizeof( data_type ), MPI_BYTE,
203 _comm_plan.neighborRank( n ), mpi_tag, _comm_plan.comm() );
205 send_range.first = send_range.second;
209 std::vector<MPI_Status> status( num_n );
211 MPI_Waitall( requests.size(), requests.data(), status.data() );
212 if ( MPI_SUCCESS != ec )
213 throw std::logic_error(
214 "Cabana::gather::apply (SliceType): Failed MPI Communication" );
217 std::size_t num_local = _comm_plan.numLocal();
218 auto extract_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
220 std::size_t ghost_idx = i + num_local;
221 auto s = SliceType::index_type::s( ghost_idx );
222 auto a = SliceType::index_type::a( ghost_idx );
223 std::size_t slice_offset = s *
slice.stride( 0 ) + a;
224 for ( std::size_t n = 0; n < num_comp; ++n )
225 slice_data[slice_offset + SliceType::vector_length * n] =
228 Kokkos::RangePolicy<ExecutionSpace> recv_policy( 0, _recv_size );
229 Kokkos::parallel_for(
"Cabana::gather::extract_recv_buffer", recv_policy,
230 extract_recv_buffer_func );
234 MPI_Barrier( _comm_plan.comm() );
246template <
class HaloType,
class SliceType>
247template <
class ExecutionSpace,
class CommSpaceType>
248std::enable_if_t<std::is_same<CommSpaceType, Mpi>::value,
void>
251 Kokkos::Profiling::ScopedRegion region(
"Cabana::scatter" );
254 auto send_buffer = this->getSendBuffer();
255 auto recv_buffer = this->getReceiveBuffer();
256 auto slice = this->getData();
259 std::size_t num_comp = this->getSliceComponents();
263 Kokkos::View<data_type*, memory_space,
264 Kokkos::MemoryTraits<Kokkos::Unmanaged>>
268 std::size_t num_local = _comm_plan.numLocal();
269 auto extract_send_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
271 std::size_t ghost_idx = i + num_local;
272 auto s = SliceType::index_type::s( ghost_idx );
273 auto a = SliceType::index_type::a( ghost_idx );
274 std::size_t slice_offset = s *
slice.stride( 0 ) + a;
275 for ( std::size_t n = 0; n < num_comp; ++n )
276 send_buffer( i, n ) =
277 slice_data( slice_offset + SliceType::vector_length * n );
279 Kokkos::RangePolicy<ExecutionSpace> send_policy( 0, _send_size );
280 Kokkos::parallel_for(
"Cabana::scatter::extract_send_buffer", send_policy,
281 extract_send_buffer_func );
285 const int mpi_tag = 2345;
288 int num_n = _comm_plan.numNeighbor();
289 std::vector<MPI_Request> requests( num_n );
290 std::pair<std::size_t, std::size_t> recv_range = { 0, 0 };
291 for (
int n = 0; n < num_n; ++n )
293 recv_range.second = recv_range.first + _comm_plan.numExport( n );
296 Kokkos::subview( recv_buffer, recv_range, Kokkos::ALL );
298 MPI_Irecv( recv_subview.data(),
299 recv_subview.size() *
sizeof( data_type ), MPI_BYTE,
300 _comm_plan.neighborRank( n ), mpi_tag, _comm_plan.comm(),
303 recv_range.first = recv_range.second;
307 std::pair<std::size_t, std::size_t> send_range = { 0, 0 };
308 for (
int n = 0; n < num_n; ++n )
310 send_range.second = send_range.first + _comm_plan.numImport( n );
313 Kokkos::subview( send_buffer, send_range, Kokkos::ALL );
315 MPI_Send( send_subview.data(),
316 send_subview.size() *
sizeof( data_type ), MPI_BYTE,
317 _comm_plan.neighborRank( n ), mpi_tag, _comm_plan.comm() );
319 send_range.first = send_range.second;
323 std::vector<MPI_Status> status( num_n );
325 MPI_Waitall( requests.size(), requests.data(), status.data() );
326 if ( MPI_SUCCESS != ec )
327 throw std::logic_error(
"Cabana::scatter::apply (SliceType): "
328 "Failed MPI Communication" );
331 auto steering = _comm_plan.getExportSteering();
334 auto scatter_recv_buffer_func = KOKKOS_LAMBDA(
const std::size_t i )
336 auto s = SliceType::index_type::s( steering( i ) );
337 auto a = SliceType::index_type::a( steering( i ) );
338 std::size_t slice_offset = s *
slice.stride( 0 ) + a;
339 for ( std::size_t n = 0; n < num_comp; ++n )
341 &slice_data( slice_offset + SliceType::vector_length * n ),
342 recv_buffer( i, n ) );
344 Kokkos::RangePolicy<ExecutionSpace> recv_policy( 0, _recv_size );
345 Kokkos::parallel_for(
"Cabana::scatter::apply::scatter_recv_buffer",
346 recv_policy, scatter_recv_buffer_func );
350 MPI_Barrier( _comm_plan.comm() );
Array-of-Struct-of-Arrays particle data structure.
Multi-node communication patterns. Base class.
Slice a single particle property from an AoSoA.
Definition Cabana_Halo.hpp:393
std::enable_if_t< std::is_same< CommSpaceType, Mpi >::value, void > applyImpl(ExecutionSpace, CommSpaceType)
Vanilla Mpi implementation of the scatter operation.
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
AoSoA_t::template member_slice_type< M > slice(const AoSoA_t &aosoa, const std::string &slice_label="")
Create a slice from an AoSoA.
Definition Cabana_AoSoA.hpp:77