14 #ifndef _ZOLTAN2_TPETRACRSMATRIXADAPTER_HPP_
15 #define _ZOLTAN2_TPETRACRSMATRIXADAPTER_HPP_
23 #include <Tpetra_CrsMatrix.hpp>
49 template <
typename User,
typename UserCoord = User>
53 #ifndef DOXYGEN_SHOULD_SKIP_THIS
60 using tmatrix_t = Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t>;
61 using device_t =
typename node_t::device_type;
62 using host_t =
typename Kokkos::HostSpace::memory_space;
64 using userCoord_t = UserCoord;
76 int nWeightsPerRow=0);
82 template <
typename Adapter>
86 template <
typename Adapter>
95 template <
typename User,
typename UserCoord>
97 const RCP<const User> &inmatrix,
int nWeightsPerRow):
98 RowMatrix(nWeightsPerRow, inmatrix) {
100 auto colIdsHost = inmatrix->getLocalIndicesHost();
102 auto colIdsGlobalHost =
103 typename Base::IdsHostView(
"colIdsGlobalHost", colIdsHost.extent(0));
104 auto colMap = inmatrix->getColMap();
107 Kokkos::parallel_for(
"colIdsGlobalHost",
108 Kokkos::RangePolicy<Kokkos::HostSpace::execution_space>(
109 0, colIdsGlobalHost.extent(0)),
111 colIdsGlobalHost(i) =
112 colMap->getGlobalElement(colIdsHost(i));
115 auto colIdsDevice = Kokkos::create_mirror_view_and_copy(
119 this->
offsDevice_ = inmatrix->getLocalRowPtrsDevice();
124 "rowWeightsDevice_", inmatrix->getLocalNumRows(),
137 template <
typename User,
typename UserCoord>
138 template <
typename Adapter>
140 const User &in, User *&out,
145 ArrayRCP<gno_t> importList;
149 (solution,
this, importList);
154 RCP<User> outPtr = this->doMigration(in, numNewRows,importList.getRawPtr());
155 out =
const_cast<User *
>(outPtr.get());
160 template <
typename User,
typename UserCoord>
161 template <
typename Adapter>
163 const User &in, RCP<User> &out,
168 ArrayRCP<gno_t> importList;
172 (solution,
this, importList);
177 out = this->doMigration(in, numNewRows, importList.getRawPtr());
Helper functions for Partitioning Problems.
Base::ConstIdsDeviceView colIdsDevice_
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
typename InputTraits< User >::scalar_t scalar_t
MatrixAdapter defines the adapter interface for matrices.
typename Kokkos::HostSpace::memory_space host_t
TpetraCrsMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
Provides access for Zoltan2 to Tpetra::CrsMatrix data.
typename node_t::device_type device_t
typename InputTraits< User >::part_t part_t
size_t getImportList(const PartitioningSolution< SolutionAdapter > &solution, const DataAdapter *const data, ArrayRCP< typename DataAdapter::gno_t > &imports)
From a PartitioningSolution, get a list of IDs to be imported. Assumes part numbers in PartitioningSo...
Provides access for Zoltan2 to Tpetra::RowMatrix data.
typename InputTraits< User >::node_t node_t
A PartitioningSolution is a solution to a partitioning problem.
typename InputTraits< User >::gno_t gno_t
Base::ConstOffsetsDeviceView offsDevice_
Defines the TpetraRowMatrixAdapter class.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
RCP< const User > matrix_
typename InputTraits< User >::offset_t offset_t
Defines the MatrixAdapter interface.
Kokkos::View< bool *, host_t > numNzWeight_
RCP< const User > getUserMatrix() const
Access to user's matrix.
Base::WeightsDeviceView rowWeightsDevice_
typename InputTraits< User >::lno_t lno_t
This file defines the StridedData class.