17 #ifndef KOKKOS_VECTORIMPORT_HPP
18 #define KOKKOS_VECTORIMPORT_HPP
26 #include <Kokkos_Core.hpp>
33 template<
class CommMessageType ,
class CommIdentType ,
class VectorType >
39 #if ! defined( KOKKOS_ENABLE_MPI )
47 template<
class CommMessageType ,
class CommIdentType ,
class VectorType >
55 const CommMessageType & ,
56 const CommMessageType & ,
57 const CommIdentType & ,
58 const unsigned arg_count_owned ,
59 const unsigned arg_count_receive )
83 template<
class CommMessageType ,
class CommIdentType ,
class VectorType >
89 ( VectorType::rank == 1 ) ||
90 std::is_same< typename VectorType::array_layout , Kokkos::LayoutRight >::value,
91 "Kokkos::Example::VectorImport Assert Fail: rank != 1 or array_layout != LayoutRight" );
93 typedef typename VectorType::HostMirror HostVectorType ;
94 typedef typename CommMessageType::HostMirror HostCommMessageType;
96 enum { ReceiveInPlace =
97 std::is_same<
typename VectorType::memory_space ,
98 typename HostVectorType::memory_space >::value };
100 const CommMessageType recv_msg ;
101 const CommMessageType send_msg ;
102 const CommIdentType send_nodeid ;
103 HostCommMessageType host_recv_msg ;
104 HostCommMessageType host_send_msg ;
105 VectorType send_buffer ;
106 HostVectorType host_send_buffer ;
107 HostVectorType host_recv_buffer ;
118 const CommIdentType index ;
119 const VectorType source ;
120 const VectorType buffer ;
122 KOKKOS_INLINE_FUNCTION
124 { buffer( i ) = source( index(i) ); }
126 Pack(
const CommIdentType & arg_index ,
127 const VectorType & arg_source ,
128 const VectorType & arg_buffer )
130 , source( arg_source )
131 , buffer( arg_buffer )
133 Kokkos::parallel_for( index.extent(0) , *this );
139 const CommMessageType & arg_recv_msg ,
140 const CommMessageType & arg_send_msg ,
141 const CommIdentType & arg_send_nodeid ,
142 const unsigned arg_count_owned ,
143 const unsigned arg_count_receive )
144 : recv_msg( arg_recv_msg )
145 , send_msg( arg_send_msg )
146 , send_nodeid( arg_send_nodeid )
160 if ( ! ReceiveInPlace ) {
161 host_recv_buffer = HostVectorType(
"recv_buffer",
count_receive);
164 unsigned send_count = 0 ;
165 for (
unsigned i = 0 ; i < send_msg.extent(0) ; ++i ) { send_count += host_send_msg(i,1); }
166 send_buffer = VectorType(
"send_buffer",send_count);
176 const Teuchos::MpiComm<int> & teuchos_mpi_comm =
dynamic_cast< const Teuchos::MpiComm<int> &
>( *comm );
178 MPI_Comm mpi_comm = * teuchos_mpi_comm.getRawMpiComm();
180 const int mpi_tag = 42 ;
181 const unsigned chunk = v.extent(1);
185 const VectorType recv_vector = Kokkos::subview( v , recv_range );
187 std::vector< MPI_Request > recv_request( recv_msg.extent(0) , MPI_REQUEST_NULL );
190 if (ReceiveInPlace) {
191 scalar_type * ptr = recv_vector.data();
193 for (
size_t i = 0 ; i < recv_msg.extent(0) ; ++i ) {
194 const int proc = host_recv_msg(i,0);
195 const int count = host_recv_msg(i,1) * chunk ;
197 MPI_Irecv( ptr , count *
sizeof(scalar_type) , MPI_BYTE ,
198 proc , mpi_tag , mpi_comm , & recv_request[i] );
204 host_scalar_type * ptr = host_recv_buffer.data();
206 for (
size_t i = 0 ; i < recv_msg.extent(0) ; ++i ) {
207 const int proc = host_recv_msg(i,0);
208 const int count = host_recv_msg(i,1) * chunk ;
210 MPI_Irecv( ptr , count *
sizeof(host_scalar_type) , MPI_BYTE ,
211 proc , mpi_tag , mpi_comm , & recv_request[i] );
218 MPI_Barrier( mpi_comm );
221 const Pack pack( send_nodeid , v , send_buffer );
225 host_scalar_type * ptr = host_send_buffer.data();
227 for (
size_t i = 0 ; i < send_msg.extent(0) ; ++i ) {
228 const int proc = host_send_msg(i,0);
229 const int count = host_send_msg(i,1) * chunk ;
239 count *
sizeof(host_scalar_type) , MPI_BYTE ,
240 proc , mpi_tag , mpi_comm );
248 for (
size_t i = 0 ; i < recv_msg.extent(0) ; ++i ) {
249 MPI_Status recv_status ;
253 MPI_Waitany( recv_msg.extent(0) , & recv_request[0] , & recv_which , & recv_status );
255 const int recv_proc = recv_status.MPI_SOURCE ;
257 MPI_Get_count( & recv_status , MPI_BYTE , & recv_size );
261 const int expected_proc = host_recv_msg(recv_which,0);
262 const int expected_size = host_recv_msg(recv_which,1) * chunk *
sizeof(
scalar_type);
264 if ( ( expected_proc != recv_proc ) ||
265 ( expected_size != recv_size ) ) {
269 MPI_Comm_rank( mpi_comm , & local_rank );
271 std::ostringstream msg ;
272 msg <<
"VectorImport error:"
273 <<
" P" << local_rank
274 <<
" received from P" << recv_proc
275 <<
" size " << recv_size
276 <<
" expected " << expected_size
277 <<
" from P" << expected_proc ;
278 throw std::runtime_error( msg.str() );
Kokkos::DefaultExecutionSpace execution_space
void operator()(const VectorType &) const
const Teuchos::RCP< const Teuchos::Comm< int > > comm
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
VectorImport(const Teuchos::RCP< const Teuchos::Comm< int > > arg_comm, const CommMessageType &, const CommMessageType &, const CommIdentType &, const unsigned arg_count_owned, const unsigned arg_count_receive)
Stokhos::CrsMatrix< ValueType, Device, Layout >::HostMirror create_mirror_view(const Stokhos::CrsMatrix< ValueType, Device, Layout > &A)
const unsigned count_receive
const unsigned count_owned