14#ifndef _ZOLTAN2_TPETRACRSMATRIXADAPTER_HPP_
15#define _ZOLTAN2_TPETRACRSMATRIXADAPTER_HPP_
23#include <Tpetra_CrsMatrix.hpp>
49template <
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);
78 void init(
const RCP<const User> &inmatrix);
84 template <
typename Adapter>
88 template <
typename Adapter>
97 template <
typename User,
typename UserCoord>
99 auto colIdsDevice = inmatrix->getLocalIndicesDevice();
101 auto colIdsGlobalDevice =
103 auto colMap = inmatrix->getColMap();
104 auto lclColMap = colMap->getLocalMap();
107 Kokkos::parallel_for(
"colIdsGlobalDevice",
108 Kokkos::RangePolicy<typename User::node_type::execution_space>(
109 0, colIdsGlobalDevice.extent(0)),
110 KOKKOS_LAMBDA(
const int i) {
111 colIdsGlobalDevice(i) =
112 lclColMap.getGlobalElement(colIdsDevice(i));
115 this->colIdsDevice_ = colIdsGlobalDevice;
116 this->offsDevice_ = inmatrix->getLocalRowPtrsDevice();
119 template <
typename User,
typename UserCoord>
121 const RCP<const User> &inmatrix,
int nWeightsPerRow):
122 RowMatrix(nWeightsPerRow, inmatrix) {
129 "rowWeightsDevice_", inmatrix->getLocalNumRows(),
130 this->nWeightsPerRow_);
142template <
typename User,
typename UserCoord>
143 template <
typename Adapter>
145 const User &in, User *&out,
150 ArrayRCP<gno_t> importList;
154 (solution,
this, importList);
159 RCP<User> outPtr = this->doMigration(in, numNewRows,importList.getRawPtr());
160 out =
const_cast<User *
>(outPtr.get());
165template <
typename User,
typename UserCoord>
166 template <
typename Adapter>
168 const User &in, RCP<User> &out,
173 ArrayRCP<gno_t> importList;
177 (solution,
this, importList);
182 out = this->doMigration(in, numNewRows, importList.getRawPtr());
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the MatrixAdapter interface.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
Defines the TpetraRowMatrixAdapter class.
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::lno_t lno_t
typename InputTraits< User >::scalar_t scalar_t
typename InputTraits< User >::gno_t gno_t
typename Kokkos::HostSpace::memory_space host_t
typename InputTraits< User >::offset_t offset_t
typename InputTraits< User >::part_t part_t
typename node_t::device_type device_t
MatrixAdapter defines the adapter interface for matrices.
A PartitioningSolution is a solution to a partitioning problem.
Provides access for Zoltan2 to Tpetra::CrsMatrix data.
void init(const RCP< const User > &inmatrix)
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
RCP< const User > getUserMatrix() const
Access to user's matrix.
TpetraCrsMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
Provides access for Zoltan2 to Tpetra::RowMatrix data.
Base::WeightsDeviceView rowWeightsDevice_
Kokkos::View< bool *, host_t > numNzWeight_
RCP< const User > matrix_
Created by mbenlioglu on Aug 31, 2020.
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...