42 #ifndef TPETRA_DETAILS_UNPACKCRSMATRIXANDCOMBINE_DECL_HPP
43 #define TPETRA_DETAILS_UNPACKCRSMATRIXANDCOMBINE_DECL_HPP
45 #include "TpetraCore_config.h"
47 #include "Kokkos_DualView.hpp"
72 #ifndef DOXYGEN_SHOULD_SKIP_THIS
75 template<
class T>
class Array;
77 template<
class T>
class ArrayView;
79 #endif // DOXYGEN_SHOULD_SKIP_THIS
83 #ifndef DOXYGEN_SHOULD_SKIP_THIS
86 #endif // DOXYGEN_SHOULD_SKIP_THIS
133 template<
typename ST,
typename LO,
typename GO,
typename NT>
136 const Teuchos::ArrayView<const char>& imports,
137 const Teuchos::ArrayView<const size_t>& numPacketsPerLID,
138 const Teuchos::ArrayView<const LO>& importLIDs,
139 size_t constantNumPackets,
140 Distributor & distor,
144 template<
typename ST,
typename LO,
typename GO,
typename NT>
146 unpackCrsMatrixAndCombineNew (
const CrsMatrix<ST, LO, GO, NT>& sourceMatrix,
147 const Kokkos::DualView<
const char*,
149 const Kokkos::DualView<
const size_t*,
151 const Kokkos::DualView<
const LO*,
153 const size_t constantNumPackets,
154 Distributor & distor,
213 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
216 const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> & sourceMatrix,
217 const Teuchos::ArrayView<const LocalOrdinal> &importLIDs,
218 const Teuchos::ArrayView<const char> &imports,
219 const Teuchos::ArrayView<const size_t>& numPacketsPerLID,
220 size_t constantNumPackets,
224 const Teuchos::ArrayView<const LocalOrdinal>& permuteToLIDs,
225 const Teuchos::ArrayView<const LocalOrdinal>& permuteFromLIDs);
241 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
244 const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> & sourceMatrix,
245 const Teuchos::ArrayView<const LocalOrdinal>& importLIDs,
246 const Teuchos::ArrayView<const char>& imports,
247 const Teuchos::ArrayView<const size_t>& numPacketsPerLID,
248 const size_t constantNumPackets,
251 const size_t numSameIDs,
252 const Teuchos::ArrayView<const LocalOrdinal>& permuteToLIDs,
253 const Teuchos::ArrayView<const LocalOrdinal>& permuteFromLIDs,
254 size_t TargetNumRows,
255 size_t TargetNumNonzeros,
256 const int MyTargetPID,
257 const Teuchos::ArrayView<size_t>& CRS_rowptr,
258 const Teuchos::ArrayView<GlobalOrdinal>& CRS_colind,
259 const Teuchos::ArrayView<
typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::impl_scalar_type>& CRS_vals,
260 const Teuchos::ArrayView<const int>& SourcePids,
261 Teuchos::Array<int>& TargetPids);
266 #endif // TPETRA_DETAILS_UNPACKCRSMATRIXANDCOMBINE_DECL_HPP
size_t unpackAndCombineWithOwningPIDsCount(const CrsGraph< LO, GO, NT > &sourceGraph, const Teuchos::ArrayView< const LO > &importLIDs, const Teuchos::ArrayView< const typename CrsGraph< LO, GO, NT >::packet_type > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode combineMode, size_t numSameIDs, const Teuchos::ArrayView< const LO > &permuteToLIDs, const Teuchos::ArrayView< const LO > &permuteFromLIDs)
Special version of Tpetra::Details::unpackCrsGraphAndCombine that also unpacks owning process ranks...
Declaration of Tpetra::CombineMode enum, and a function for setting a Tpetra::CombineMode parameter i...
CombineMode
Rule for combining data in an Import or Export.
void unpackAndCombineIntoCrsArrays(const CrsGraph< LO, GO, NT > &sourceGraph, const Teuchos::ArrayView< const LO > &importLIDs, const Teuchos::ArrayView< const typename CrsGraph< LO, GO, NT >::packet_type > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, const size_t constantNumPackets, Distributor &distor, const CombineMode combineMode, const size_t numSameIDs, const Teuchos::ArrayView< const LO > &permuteToLIDs, const Teuchos::ArrayView< const LO > &permuteFromLIDs, size_t TargetNumRows, size_t TargetNumNonzeros, const int MyTargetPID, const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< GO > &CRS_colind, const Teuchos::ArrayView< const int > &SourcePids, Teuchos::Array< int > &TargetPids)
unpackAndCombineIntoCrsArrays
Forward declaration of Tpetra::CrsMatrix.
Declaration of the Tpetra::DistObject class.
void unpackCrsMatrixAndCombine(const CrsMatrix< ST, LO, GO, NT > &sourceMatrix, const Teuchos::ArrayView< const char > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, const Teuchos::ArrayView< const LO > &importLIDs, size_t constantNumPackets, Distributor &distor, CombineMode combineMode, const bool atomic)
Unpack the imported column indices and values, and combine into matrix.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type