Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_TpetraCrsMatrix_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Xpetra: A linear algebra interface package
4 //
5 // Copyright 2012 NTESS and the Xpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef XPETRA_TPETRACRSMATRIX_DEF_HPP
11 #define XPETRA_TPETRACRSMATRIX_DEF_HPP
12 
13 #include <Xpetra_MultiVectorFactory.hpp>
15 #include "Tpetra_Details_residual.hpp"
16 
17 namespace Xpetra {
18 
19 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
20 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP<Teuchos::ParameterList> &params)
21  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), maxNumEntriesPerRow, params))) {}
22 
23 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
24 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc, const Teuchos::RCP<Teuchos::ParameterList> &params)
25  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), NumEntriesPerRowToAlloc(), params))) {}
26 
27 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap, const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap, size_t maxNumEntriesPerRow, const Teuchos::RCP<Teuchos::ParameterList> &params)
29  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, params))) {}
30 
31 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
32 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap, const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc, const Teuchos::RCP<Teuchos::ParameterList> &params)
33  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc(), params))) {}
34 
35 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> &graph, const Teuchos::RCP<Teuchos::ParameterList> &params)
37  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), params))) {}
38 
39 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> &graph, typename Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_type::values_type &values, const Teuchos::RCP<Teuchos::ParameterList> &params)
41  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), values, params))) {}
42 
43 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
44 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
45  const Import<LocalOrdinal, GlobalOrdinal, Node> &importer,
46  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
47  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
48  const Teuchos::RCP<Teuchos::ParameterList> &params) {
49  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
50  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
51  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
52 
53  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
54  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
55  mtx_ = Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(importer), myDomainMap, myRangeMap, params);
56  bool restrictComm = false;
57  if (!params.is_null()) restrictComm = params->get("Restrict Communicator", restrictComm);
58  if (restrictComm && mtx_->getRowMap().is_null()) mtx_ = Teuchos::null;
59 }
60 
61 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
62 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
63  const Export<LocalOrdinal, GlobalOrdinal, Node> &exporter,
64  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
65  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
66  const Teuchos::RCP<Teuchos::ParameterList> &params) {
67  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
68  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
69  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
70 
71  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
72  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
73  mtx_ = Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(exporter), myDomainMap, myRangeMap, params);
74 }
75 
76 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
77 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
78  const Import<LocalOrdinal, GlobalOrdinal, Node> &RowImporter,
79  const Teuchos::RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> DomainImporter,
80  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
81  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
82  const Teuchos::RCP<Teuchos::ParameterList> &params) {
83  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
84  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
85  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
86 
87  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
88  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
89 
90  mtx_ = Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(RowImporter), toTpetra(*DomainImporter), myDomainMap, myRangeMap, params);
91  bool restrictComm = false;
92  if (!params.is_null()) restrictComm = params->get("Restrict Communicator", restrictComm);
93  if (restrictComm && mtx_->getRowMap().is_null()) mtx_ = Teuchos::null;
94 }
95 
96 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
97 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
98  const Export<LocalOrdinal, GlobalOrdinal, Node> &RowExporter,
99  const Teuchos::RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> DomainExporter,
100  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
101  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
102  const Teuchos::RCP<Teuchos::ParameterList> &params) {
103  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
104  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
105  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
106 
107  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
108  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
109 
110  mtx_ = Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(RowExporter), toTpetra(*DomainExporter), myDomainMap, myRangeMap, params);
111 }
112 
114 
115 #ifdef HAVE_XPETRA_TPETRA
116 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
117 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap,
118  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap,
119  const local_matrix_type &lclMatrix,
120  const Teuchos::RCP<Teuchos::ParameterList> &params)
121  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) {}
122 
123 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
124 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(
125  const local_matrix_type &lclMatrix,
126  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap,
127  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap,
128  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
129  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
130  const Teuchos::RCP<Teuchos::ParameterList> &params)
131  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) {}
132 
133 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(
135  const local_matrix_type &lclMatrix,
136  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap,
137  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap,
138  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
139  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
140  const Teuchos::RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> &importer,
141  const Teuchos::RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> &exporter,
142  const Teuchos::RCP<Teuchos::ParameterList> &params)
143  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), toTpetra(importer), toTpetra(exporter), params))) {}
144 #endif
145 
146 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148 
149 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
150 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
151  XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues");
152  mtx_->insertGlobalValues(globalRow, cols, vals);
153 }
154 
155 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
156 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
157  XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues");
158  mtx_->insertLocalValues(localRow, cols, vals);
159 }
160 
161 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
162 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
163  XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues");
164  mtx_->replaceGlobalValues(globalRow, cols, vals);
165 }
166 
167 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
169  const ArrayView<const LocalOrdinal> &cols,
170  const ArrayView<const Scalar> &vals) {
171  XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
172  typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
173  const LocalOrdinal numValid =
174  mtx_->replaceLocalValues(localRow, cols, vals);
175  TEUCHOS_TEST_FOR_EXCEPTION(
176  static_cast<size_type>(numValid) != cols.size(), std::runtime_error,
177  "replaceLocalValues returned " << numValid << " != cols.size() = " << cols.size() << ".");
178 }
179 
180 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182  XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar");
183  mtx_->setAllToScalar(alpha);
184 }
185 
186 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
188  XPETRA_MONITOR("TpetraCrsMatrix::scale");
189  mtx_->scale(alpha);
190 }
191 
192 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
193 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::allocateAllValues(size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) {
194  XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues");
195  rowptr.resize(getLocalNumRows() + 1);
196  colind.resize(numNonZeros);
197  values.resize(numNonZeros);
198 }
199 
200 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setAllValues(const ArrayRCP<size_t> &rowptr, const ArrayRCP<LocalOrdinal> &colind, const ArrayRCP<Scalar> &values) {
202  XPETRA_MONITOR("TpetraCrsMatrix::setAllValues");
203  mtx_->setAllValues(rowptr, colind, values);
204 }
205 
206 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
207 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values) const {
208  XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
209  // TODO: Change getAllValues interface to return Kokkos views,
210  // TODO: or add getLocalRowPtrsHost, getLocalIndicesHost, etc., interfaces
211  // TODO: and use them in MueLu
212  rowptr = Kokkos::Compat::persistingView(mtx_->getLocalRowPtrsHost());
213  colind = Kokkos::Compat::persistingView(mtx_->getLocalIndicesHost());
214  values = Teuchos::arcp_reinterpret_cast<const Scalar>(
215  Kokkos::Compat::persistingView(mtx_->getLocalValuesHost(
217 }
218 
219 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221  XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
222  // TODO: Change getAllValues interface to return Kokkos view,
223  // TODO: or add getLocalValuesDevice() interfaces
224  values = Teuchos::arcp_reinterpret_cast<Scalar>(
225  Kokkos::Compat::persistingView(mtx_->getLocalValuesDevice(
227 }
228 
229 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231 
232 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
234  XPETRA_MONITOR("TpetraCrsMatrix::resumeFill");
235  mtx_->resumeFill(params);
236 }
237 
238 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
239 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::fillComplete(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap, const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap, const RCP<ParameterList> &params) {
240  XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
241  mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params);
242 }
243 
244 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246  XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
247  mtx_->fillComplete(params);
248 }
249 
250 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
252  XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
253  XPETRA_DYNAMIC_CAST(const TpetraImportClass, *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
254  RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> myImport = tImporter.getTpetra_Import();
255  mtx_->replaceDomainMapAndImporter(toTpetra(newDomainMap), myImport);
256 }
257 
258 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260  const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
261  const RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> &importer,
262  const RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> &exporter,
263  const RCP<ParameterList> &params) {
264  XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
265  RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> myImport;
266  RCP<const Tpetra::Export<LocalOrdinal, GlobalOrdinal, Node>> myExport;
267 
268  if (importer != Teuchos::null) {
269  XPETRA_DYNAMIC_CAST(const TpetraImportClass, *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
270  myImport = tImporter.getTpetra_Import();
271  }
272  if (exporter != Teuchos::null) {
273  XPETRA_DYNAMIC_CAST(const TpetraExportClass, *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
274  myExport = tExporter.getTpetra_Export();
275  }
276 
277  mtx_->expertStaticFillComplete(toTpetra(domainMap), toTpetra(rangeMap), myImport, myExport, params);
278 }
279 
280 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getRowMap() const {
282  XPETRA_MONITOR("TpetraCrsMatrix::getRowMap");
283  return toXpetra(mtx_->getRowMap());
284 }
285 
286 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
287 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getColMap() const {
288  XPETRA_MONITOR("TpetraCrsMatrix::getColMap");
289  return toXpetra(mtx_->getColMap());
290 }
291 
292 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
293 RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getCrsGraph() const {
294  XPETRA_MONITOR("TpetraCrsMatrix::getCrsGraph");
295  return toXpetra(mtx_->getCrsGraph());
296 }
297 
298 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows");
301  return mtx_->getGlobalNumRows();
302 }
303 
304 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
306  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols");
307  return mtx_->getGlobalNumCols();
308 }
309 
310 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
312  XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumRows");
313  return mtx_->getLocalNumRows();
314 }
315 
316 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
318  XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumCols");
319  return mtx_->getLocalNumCols();
320 }
321 
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
324  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries");
325  return mtx_->getGlobalNumEntries();
326 }
327 
328 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
330  XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumEntries");
331  return mtx_->getLocalNumEntries();
332 }
333 
334 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336  XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow");
337  return mtx_->getNumEntriesInLocalRow(localRow);
338 }
339 
340 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
342  XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow");
343  return mtx_->getNumEntriesInGlobalRow(globalRow);
344 }
345 
346 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries");
349  return mtx_->getGlobalMaxNumRowEntries();
350 }
351 
352 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354  XPETRA_MONITOR("TpetraCrsMatrix::getLocalMaxNumRowEntries");
355  return mtx_->getLocalMaxNumRowEntries();
356 }
357 
358 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
360  XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed");
361  return mtx_->isLocallyIndexed();
362 }
363 
364 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
366  XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed");
367  return mtx_->isGloballyIndexed();
368 }
369 
370 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372  XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete");
373  return mtx_->isFillComplete();
374 }
375 
376 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
378  XPETRA_MONITOR("TpetraCrsMatrix::isFillActive");
379  return mtx_->isFillActive();
380 }
381 
382 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
383 typename ScalarTraits<Scalar>::magnitudeType TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getFrobeniusNorm() const {
384  XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm");
385  return mtx_->getFrobeniusNorm();
386 }
387 
388 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390  XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews");
391  return mtx_->supportsRowViews();
392 }
393 
394 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
395 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView<LocalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
396  XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy");
397  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_local_inds_host_view_type indices("indices", Indices.size());
398  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values("values", Values.size());
399 
400  mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
401  for (size_t i = 0; i < NumEntries; ++i) {
402  Indices[i] = indices(i);
403  Values[i] = values(i);
404  }
405 }
406 
407 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
408 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView<GlobalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
409  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy");
410  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_global_inds_host_view_type indices("indices", Indices.size());
411  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values("values", Values.size());
412 
413  mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
414  for (size_t i = 0; i < NumEntries; ++i) {
415  Indices[i] = indices(i);
416  Values[i] = values(i);
417  }
418 }
419 
420 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
421 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &Indices, ArrayView<const Scalar> &Values) const {
422  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView");
423  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type indices;
424  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type values;
425 
426  mtx_->getGlobalRowView(GlobalRow, indices, values);
427  Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
428  Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
429 }
430 
431 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
432 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &Indices, ArrayView<const Scalar> &Values) const {
433  XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView");
434  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type indices;
435  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type values;
436 
437  mtx_->getLocalRowView(LocalRow, indices, values);
438  Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
439  Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
440 }
441 
442 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
443 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta) const {
444  XPETRA_MONITOR("TpetraCrsMatrix::apply");
445  mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
446 }
447 
448 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
449 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter, const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const {
450  XPETRA_MONITOR("TpetraCrsMatrix::apply(region)");
451  RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> regionInterfaceMap = regionInterfaceImporter->getTargetMap();
452  mtx_->localApply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
453  if (sumInterfaceValues) {
454  // preform communication to propagate local interface
455  // values to all the processor that share interfaces.
456  RCP<MultiVector> matvecInterfaceTmp = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(regionInterfaceMap, 1);
457  matvecInterfaceTmp->doImport(Y, *regionInterfaceImporter, INSERT);
458 
459  // sum all contributions to interface values
460  // on all ranks
461  ArrayRCP<Scalar> YData = Y.getDataNonConst(0);
462  ArrayRCP<Scalar> interfaceData = matvecInterfaceTmp->getDataNonConst(0);
463  for (LocalOrdinal interfaceIdx = 0; interfaceIdx < static_cast<LocalOrdinal>(interfaceData.size()); ++interfaceIdx) {
464  YData[regionInterfaceLIDs[interfaceIdx]] += interfaceData[interfaceIdx];
465  }
466  }
467 }
468 
469 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
470 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getDomainMap() const {
471  XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap");
472  return toXpetra(mtx_->getDomainMap());
473 }
474 
475 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
476 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getRangeMap() const {
477  XPETRA_MONITOR("TpetraCrsMatrix::getRangeMap");
478  return toXpetra(mtx_->getRangeMap());
479 }
480 
481 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
483  XPETRA_MONITOR("TpetraCrsMatrix::description");
484  return mtx_->description();
485 }
486 
487 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
488 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const {
489  XPETRA_MONITOR("TpetraCrsMatrix::describe");
490  mtx_->describe(out, verbLevel);
491 }
492 
493 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
495  XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
496  Teuchos::LabeledObject::setObjectLabel(objectLabel);
497  mtx_->setObjectLabel(objectLabel);
498 }
499 
500 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
502  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(*(matrix.mtx_), Teuchos::Copy))) {}
503 
504 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
506  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
507  XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
508  mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
509 }
510 
511 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
513  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
514  mtx_->getLocalDiagOffsets(offsets);
515 }
516 
517 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
518 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagCopy(Vector &diag, const Teuchos::ArrayView<const size_t> &offsets) const {
519  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
520  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
521 }
522 
523 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagCopy(Vector &diag, const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const {
525  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
526  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
527 }
528 
529 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
531  XPETRA_MONITOR("TpetraCrsMatrix::replaceDiag");
532  Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(toTpetra(diag)));
533 }
534 
535 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
537  XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
538  mtx_->leftScale(*(toTpetra(x)));
539 }
540 
541 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
543  XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
544  mtx_->rightScale(*(toTpetra(x)));
545 }
546 
547 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
548 Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getMap() const {
549  XPETRA_MONITOR("TpetraCrsMatrix::getMap");
550  return rcp(new TpetraMap<LocalOrdinal, GlobalOrdinal, Node>(mtx_->getMap()));
551 }
552 
553 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
556  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
557 
558  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
559  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSource.getTpetra_CrsMatrix();
560  // mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
561  mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
562 }
563 
565 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
568  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
569 
570  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
571  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tDest.getTpetra_CrsMatrix();
572  mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
573 }
574 
575 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
578  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
579 
580  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
581  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSource.getTpetra_CrsMatrix();
582  mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
583 }
584 
585 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
588  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
589 
590  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
591  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tDest.getTpetra_CrsMatrix();
592  mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
593 }
594 
595 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
597  XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
598  mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
599 }
600 
601 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
603  return !mtx_.is_null();
604 }
605 
606 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
607 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &mtx)
608  : mtx_(mtx) {}
609 
610 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
611 RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getTpetra_CrsMatrix() const { return mtx_; }
612 
614 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
615 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getTpetra_CrsMatrixNonConst() const { return mtx_; } // TODO: remove
616 
618 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
620  const MultiVector &B,
621  MultiVector &R) const {
622  Tpetra::Details::residual(*mtx_, toTpetra(X), toTpetra(B), toTpetra(R));
623 }
624 
627 // End of TpetrCrsMatrix class definition //
630 
631 } // namespace Xpetra
632 
633 #endif // XPETRA_TPETRACRSMATRIX_DEF_HPP
void replaceDiag(const Vector &diag)
Replace the diagonal entries of the matrix.
void leftScale(const Vector &x)
Left scale operator with given vector values.
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &map, size_t NumVectors, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
constexpr struct ReadOnlyStruct ReadOnly
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)=0
View of the local values in a particular vector of this multivector.
size_t global_size_t
Global size_t object.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
bool hasMatrix() const
Does this have an underlying matrix.
CombineMode
Xpetra::Combine Mode enumerable type.
constexpr struct ReadWriteStruct ReadWrite
#define XPETRA_MONITOR(funcName)
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
void rightScale(const Vector &x)
Right scale operator with given vector values.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)