15 #ifndef __IFPACK2_CRSARRAYS_DECL_HPP__
16 #define __IFPACK2_CRSARRAYS_DECL_HPP__
18 #include <Tpetra_RowMatrix.hpp>
19 #include <Tpetra_CrsMatrix.hpp>
20 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
21 #include <Tpetra_BlockCrsMatrix_Helpers_decl.hpp>
22 #include <KokkosSparse_CrsMatrix.hpp>
23 #include <Ifpack2_LocalFilter.hpp>
24 #include <Ifpack2_ReorderFilter.hpp>
31 template <
typename Scalar,
typename ImplScalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
32 struct CrsArrayReader {
33 typedef typename Node::device_type device_type;
34 typedef typename device_type::execution_space execution_space;
35 typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TRowMatrix;
36 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TCrsMatrix;
37 typedef Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> TBcrsMatrix;
40 typedef KokkosSparse::CrsMatrix<ImplScalar, LocalOrdinal, execution_space> KCrsMatrix;
41 typedef Kokkos::View<LocalOrdinal*, execution_space> OrdinalArray;
42 typedef Kokkos::View<ImplScalar*, execution_space> ScalarArray;
43 typedef typename OrdinalArray::HostMirror OrdinalArrayHost;
45 typedef Kokkos::Serial functor_space;
46 typedef Kokkos::RangePolicy<functor_space, int> RangePol;
52 static void getValues(
const TRowMatrix* A, ScalarArray& vals, OrdinalArrayHost& rowptrs) {
53 auto Acrs =
dynamic_cast<const TCrsMatrix*
>(A);
54 auto Abcrs =
dynamic_cast<const TBcrsMatrix*
>(A);
56 getValuesCrs(Acrs, vals);
60 getValuesBcrs(Abcrs, vals);
63 using range_type = Kokkos::pair<int, int>;
64 using local_inds_host_view_type =
typename TRowMatrix::nonconst_local_inds_host_view_type;
65 using values_host_view_type =
typename TRowMatrix::nonconst_values_host_view_type;
66 using scalar_type =
typename values_host_view_type::value_type;
68 LocalOrdinal nrows = A->getLocalNumRows();
69 size_t nnz = A->getLocalNumEntries();
70 size_t maxNnz = A->getLocalMaxNumRowEntries();
72 vals = ScalarArray(
"Values", nnz);
73 auto valsHost = Kokkos::create_mirror(vals);
74 local_inds_host_view_type lclColInds(
"lclColinds", maxNnz);
77 for (LocalOrdinal i = 0; i < nrows; i++) {
78 size_t NumEntries = A->getNumEntriesInLocalRow(i);
79 auto constLclValues = Kokkos::subview(valsHost, range_type(nnz, nnz + NumEntries));
80 values_host_view_type lclValues(const_cast<scalar_type*>(constLclValues.data()), NumEntries);
82 A->getLocalRowCopy(i, lclColInds, lclValues, NumEntries);
85 Kokkos::deep_copy(vals, valsHost);
93 static void getStructure(
const TRowMatrix* A, OrdinalArrayHost& rowptrsHost, OrdinalArray& rowptrs, OrdinalArray& colinds) {
94 auto Acrs =
dynamic_cast<const TCrsMatrix*
>(A);
95 auto Abcrs =
dynamic_cast<const TBcrsMatrix*
>(A);
97 getStructureCrs(Acrs, rowptrsHost, rowptrs, colinds);
101 getStructureBcrs(Abcrs, rowptrsHost, rowptrs, colinds);
107 LocalOrdinal nrows = A->getLocalNumRows();
108 rowptrsHost = OrdinalArrayHost(
"RowPtrs (host)", nrows + 1);
110 using range_type = Kokkos::pair<int, int>;
111 using values_host_view_type =
typename TRowMatrix::nonconst_values_host_view_type;
112 using local_inds_host_view_type =
typename TRowMatrix::nonconst_local_inds_host_view_type;
113 using local_ind_type =
typename local_inds_host_view_type::value_type;
114 size_t nnz = A->getLocalNumEntries();
115 size_t maxNnz = A->getLocalMaxNumRowEntries();
117 colinds = OrdinalArray(
"ColInds", nnz);
118 auto colindsHost = Kokkos::create_mirror(colinds);
119 values_host_view_type lclValues(
"lclValues", maxNnz);
122 rowptrsHost[0] = nnz;
123 for (LocalOrdinal i = 0; i < nrows; i++) {
124 size_t NumEntries = A->getNumEntriesInLocalRow(i);
125 auto constLclValues = Kokkos::subview(colindsHost, range_type(nnz, nnz + NumEntries));
126 local_inds_host_view_type lclColInds(const_cast<local_ind_type*>(constLclValues.data()), NumEntries);
127 A->getLocalRowCopy(i, lclColInds, lclValues, NumEntries);
130 rowptrsHost[i + 1] = nnz;
133 rowptrs = OrdinalArray(
"RowPtrs", nrows + 1);
134 Kokkos::deep_copy(rowptrs, rowptrsHost);
135 Kokkos::deep_copy(colinds, colindsHost);
140 static void getValuesCrs(
const TCrsMatrix* A, ScalarArray& values_) {
141 auto localA = A->getLocalMatrixDevice();
142 auto values = localA.values;
143 auto nnz = values.extent(0);
144 values_ = ScalarArray(
"Values", nnz);
145 Kokkos::deep_copy(values_, values);
149 static void getStructureCrs(
const TCrsMatrix* A, OrdinalArrayHost& rowptrsHost_, OrdinalArray& rowptrs_, OrdinalArray& colinds_) {
151 auto localA = A->getLocalMatrixDevice();
152 auto rowptrs = localA.graph.row_map;
153 auto colinds = localA.graph.entries;
154 auto numRows = A->getLocalNumRows();
155 auto nnz = colinds.extent(0);
157 rowptrs_ = OrdinalArray(
"RowPtrs", numRows + 1);
158 colinds_ = OrdinalArray(
"ColInds", nnz);
159 Kokkos::deep_copy(rowptrs_, rowptrs);
160 Kokkos::deep_copy(colinds_, colinds);
162 rowptrsHost_ = Kokkos::create_mirror(rowptrs_);
163 Kokkos::deep_copy(rowptrsHost_, rowptrs_);
167 static void getValuesBcrs(
const TBcrsMatrix* A, ScalarArray& values_) {
168 auto localA = A->getLocalMatrixDevice();
169 auto values = localA.values;
170 auto nnz = values.extent(0);
171 values_ = ScalarArray(
"Values", nnz);
172 Kokkos::deep_copy(values_, values);
176 static void getStructureBcrs(
const TBcrsMatrix* A, OrdinalArrayHost& rowptrsHost_, OrdinalArray& rowptrs_, OrdinalArray& colinds_) {
178 auto localA = A->getLocalMatrixDevice();
179 auto rowptrs = localA.graph.row_map;
180 auto colinds = localA.graph.entries;
181 auto numRows = A->getLocalNumRows();
182 auto nnz = colinds.extent(0);
184 rowptrs_ = OrdinalArray(
"RowPtrs", numRows + 1);
185 colinds_ = OrdinalArray(
"ColInds", nnz);
186 Kokkos::deep_copy(rowptrs_, rowptrs);
187 Kokkos::deep_copy(colinds_, colinds);
189 rowptrsHost_ = Kokkos::create_mirror(rowptrs_);
190 Kokkos::deep_copy(rowptrsHost_, rowptrs_);
Wraps a Tpetra::RowMatrix in a filter that reorders local rows and columns.
Definition: Ifpack2_ReorderFilter_decl.hpp:36
Access only local rows and columns of a sparse matrix.
Definition: Ifpack2_LocalFilter_decl.hpp:127