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 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_TPETRACRSMATRIX_DEF_HPP
47 #define XPETRA_TPETRACRSMATRIX_DEF_HPP
48 
49 #include <Xpetra_MultiVectorFactory.hpp>
51 #include "Tpetra_Details_residual.hpp"
52 
53 namespace Xpetra {
54 
55 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP<Teuchos::ParameterList> &params)
57  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), maxNumEntriesPerRow, params))) {}
58 
59 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
60 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)
61  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), NumEntriesPerRowToAlloc(), params))) {}
62 
63 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64 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)
65  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, params))) {}
66 
67 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68 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)
69  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc(), params))) {}
70 
71 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> &graph, const Teuchos::RCP<Teuchos::ParameterList> &params)
73  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), params))) {}
74 
75 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
76 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)
77  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), values, params))) {}
78 
79 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
81  const Import<LocalOrdinal, GlobalOrdinal, Node> &importer,
82  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
83  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
84  const Teuchos::RCP<Teuchos::ParameterList> &params) {
85  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
86  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
87  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
88 
89  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
90  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
91  mtx_ = Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(importer), myDomainMap, myRangeMap, params);
92  bool restrictComm = false;
93  if (!params.is_null()) restrictComm = params->get("Restrict Communicator", restrictComm);
94  if (restrictComm && mtx_->getRowMap().is_null()) mtx_ = Teuchos::null;
95 }
96 
97 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
99  const Export<LocalOrdinal, GlobalOrdinal, Node> &exporter,
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  mtx_ = Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(exporter), myDomainMap, myRangeMap, params);
110 }
111 
112 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
113 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
114  const Import<LocalOrdinal, GlobalOrdinal, Node> &RowImporter,
115  const Teuchos::RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> DomainImporter,
116  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
117  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
118  const Teuchos::RCP<Teuchos::ParameterList> &params) {
119  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
120  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
121  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
122 
123  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
124  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
125 
126  mtx_ = Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(RowImporter), toTpetra(*DomainImporter), myDomainMap, myRangeMap, params);
127  bool restrictComm = false;
128  if (!params.is_null()) restrictComm = params->get("Restrict Communicator", restrictComm);
129  if (restrictComm && mtx_->getRowMap().is_null()) mtx_ = Teuchos::null;
130 }
131 
132 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
133 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const CrsMatrix> &sourceMatrix,
134  const Export<LocalOrdinal, GlobalOrdinal, Node> &RowExporter,
135  const Teuchos::RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> DomainExporter,
136  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
137  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
138  const Teuchos::RCP<Teuchos::ParameterList> &params) {
139  typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
140  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
141  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
142 
143  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
144  RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
145 
146  mtx_ = Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(RowExporter), toTpetra(*DomainExporter), myDomainMap, myRangeMap, params);
147 }
148 
150 
151 #ifdef HAVE_XPETRA_TPETRA
152 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
153 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap,
154  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap,
155  const local_matrix_type &lclMatrix,
156  const Teuchos::RCP<Teuchos::ParameterList> &params)
157  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) {}
158 
159 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
160 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(
161  const local_matrix_type &lclMatrix,
162  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap,
163  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap,
164  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
165  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
166  const Teuchos::RCP<Teuchos::ParameterList> &params)
167  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) {}
168 
169 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(
171  const local_matrix_type &lclMatrix,
172  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rowMap,
173  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &colMap,
174  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &domainMap,
175  const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
176  const Teuchos::RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> &importer,
177  const Teuchos::RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> &exporter,
178  const Teuchos::RCP<Teuchos::ParameterList> &params)
179  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), toTpetra(importer), toTpetra(exporter), params))) {}
180 #endif
181 
182 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
184 
185 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
186 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
187  XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues");
188  mtx_->insertGlobalValues(globalRow, cols, vals);
189 }
190 
191 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
192 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
193  XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues");
194  mtx_->insertLocalValues(localRow, cols, vals);
195 }
196 
197 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
198 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
199  XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues");
200  mtx_->replaceGlobalValues(globalRow, cols, vals);
201 }
202 
203 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
205  const ArrayView<const LocalOrdinal> &cols,
206  const ArrayView<const Scalar> &vals) {
207  XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
208  typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
209  const LocalOrdinal numValid =
210  mtx_->replaceLocalValues(localRow, cols, vals);
211  TEUCHOS_TEST_FOR_EXCEPTION(
212  static_cast<size_type>(numValid) != cols.size(), std::runtime_error,
213  "replaceLocalValues returned " << numValid << " != cols.size() = " << cols.size() << ".");
214 }
215 
216 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
218  XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar");
219  mtx_->setAllToScalar(alpha);
220 }
221 
222 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224  XPETRA_MONITOR("TpetraCrsMatrix::scale");
225  mtx_->scale(alpha);
226 }
227 
228 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
229 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::allocateAllValues(size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) {
230  XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues");
231  rowptr.resize(getLocalNumRows() + 1);
232  colind.resize(numNonZeros);
233  values.resize(numNonZeros);
234 }
235 
236 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
237 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setAllValues(const ArrayRCP<size_t> &rowptr, const ArrayRCP<LocalOrdinal> &colind, const ArrayRCP<Scalar> &values) {
238  XPETRA_MONITOR("TpetraCrsMatrix::setAllValues");
239  mtx_->setAllValues(rowptr, colind, values);
240 }
241 
242 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
243 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values) const {
244  XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
245  // TODO: Change getAllValues interface to return Kokkos views,
246  // TODO: or add getLocalRowPtrsHost, getLocalIndicesHost, etc., interfaces
247  // TODO: and use them in MueLu
248  rowptr = Kokkos::Compat::persistingView(mtx_->getLocalRowPtrsHost());
249  colind = Kokkos::Compat::persistingView(mtx_->getLocalIndicesHost());
250  values = Teuchos::arcp_reinterpret_cast<const Scalar>(
251  Kokkos::Compat::persistingView(mtx_->getLocalValuesHost(
253 }
254 
255 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
257  XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
258  // TODO: Change getAllValues interface to return Kokkos view,
259  // TODO: or add getLocalValuesDevice() interfaces
260  values = Teuchos::arcp_reinterpret_cast<Scalar>(
261  Kokkos::Compat::persistingView(mtx_->getLocalValuesDevice(
263 }
264 
265 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
267 
268 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
270  XPETRA_MONITOR("TpetraCrsMatrix::resumeFill");
271  mtx_->resumeFill(params);
272 }
273 
274 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
275 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) {
276  XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
277  mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params);
278 }
279 
280 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
282  XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
283  mtx_->fillComplete(params);
284 }
285 
286 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
288  XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
289  XPETRA_DYNAMIC_CAST(const TpetraImportClass, *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
290  RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> myImport = tImporter.getTpetra_Import();
291  mtx_->replaceDomainMapAndImporter(toTpetra(newDomainMap), myImport);
292 }
293 
294 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
296  const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
297  const RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> &importer,
298  const RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> &exporter,
299  const RCP<ParameterList> &params) {
300  XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
301  RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> myImport;
302  RCP<const Tpetra::Export<LocalOrdinal, GlobalOrdinal, Node>> myExport;
303 
304  if (importer != Teuchos::null) {
305  XPETRA_DYNAMIC_CAST(const TpetraImportClass, *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
306  myImport = tImporter.getTpetra_Import();
307  }
308  if (exporter != Teuchos::null) {
309  XPETRA_DYNAMIC_CAST(const TpetraExportClass, *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
310  myExport = tExporter.getTpetra_Export();
311  }
312 
313  mtx_->expertStaticFillComplete(toTpetra(domainMap), toTpetra(rangeMap), myImport, myExport, params);
314 }
315 
316 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
317 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getRowMap() const {
318  XPETRA_MONITOR("TpetraCrsMatrix::getRowMap");
319  return toXpetra(mtx_->getRowMap());
320 }
321 
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
323 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getColMap() const {
324  XPETRA_MONITOR("TpetraCrsMatrix::getColMap");
325  return toXpetra(mtx_->getColMap());
326 }
327 
328 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329 RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getCrsGraph() const {
330  XPETRA_MONITOR("TpetraCrsMatrix::getCrsGraph");
331  return toXpetra(mtx_->getCrsGraph());
332 }
333 
334 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows");
337  return mtx_->getGlobalNumRows();
338 }
339 
340 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
342  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols");
343  return mtx_->getGlobalNumCols();
344 }
345 
346 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348  XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumRows");
349  return mtx_->getLocalNumRows();
350 }
351 
352 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354  XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumCols");
355  return mtx_->getLocalNumCols();
356 }
357 
358 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
360  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries");
361  return mtx_->getGlobalNumEntries();
362 }
363 
364 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
366  XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumEntries");
367  return mtx_->getLocalNumEntries();
368 }
369 
370 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372  XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow");
373  return mtx_->getNumEntriesInLocalRow(localRow);
374 }
375 
376 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
378  XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow");
379  return mtx_->getNumEntriesInGlobalRow(globalRow);
380 }
381 
382 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries");
385  return mtx_->getGlobalMaxNumRowEntries();
386 }
387 
388 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390  XPETRA_MONITOR("TpetraCrsMatrix::getLocalMaxNumRowEntries");
391  return mtx_->getLocalMaxNumRowEntries();
392 }
393 
394 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
396  XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed");
397  return mtx_->isLocallyIndexed();
398 }
399 
400 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
402  XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed");
403  return mtx_->isGloballyIndexed();
404 }
405 
406 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
408  XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete");
409  return mtx_->isFillComplete();
410 }
411 
412 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
414  XPETRA_MONITOR("TpetraCrsMatrix::isFillActive");
415  return mtx_->isFillActive();
416 }
417 
418 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
419 typename ScalarTraits<Scalar>::magnitudeType TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getFrobeniusNorm() const {
420  XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm");
421  return mtx_->getFrobeniusNorm();
422 }
423 
424 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
426  XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews");
427  return mtx_->supportsRowViews();
428 }
429 
430 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
431 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView<LocalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
432  XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy");
433  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_local_inds_host_view_type indices("indices", Indices.size());
434  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values("values", Values.size());
435 
436  mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
437  for (size_t i = 0; i < NumEntries; ++i) {
438  Indices[i] = indices(i);
439  Values[i] = values(i);
440  }
441 }
442 
443 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
444 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView<GlobalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
445  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy");
446  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_global_inds_host_view_type indices("indices", Indices.size());
447  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values("values", Values.size());
448 
449  mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
450  for (size_t i = 0; i < NumEntries; ++i) {
451  Indices[i] = indices(i);
452  Values[i] = values(i);
453  }
454 }
455 
456 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
457 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &Indices, ArrayView<const Scalar> &Values) const {
458  XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView");
459  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type indices;
460  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type values;
461 
462  mtx_->getGlobalRowView(GlobalRow, indices, values);
463  Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
464  Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
465 }
466 
467 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
468 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &Indices, ArrayView<const Scalar> &Values) const {
469  XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView");
470  typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type indices;
471  typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type values;
472 
473  mtx_->getLocalRowView(LocalRow, indices, values);
474  Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
475  Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
476 }
477 
478 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
479 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta) const {
480  XPETRA_MONITOR("TpetraCrsMatrix::apply");
481  mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
482 }
483 
484 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
485 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 {
486  XPETRA_MONITOR("TpetraCrsMatrix::apply(region)");
487  RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> regionInterfaceMap = regionInterfaceImporter->getTargetMap();
488  mtx_->localApply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
489  if (sumInterfaceValues) {
490  // preform communication to propagate local interface
491  // values to all the processor that share interfaces.
492  RCP<MultiVector> matvecInterfaceTmp = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(regionInterfaceMap, 1);
493  matvecInterfaceTmp->doImport(Y, *regionInterfaceImporter, INSERT);
494 
495  // sum all contributions to interface values
496  // on all ranks
497  ArrayRCP<Scalar> YData = Y.getDataNonConst(0);
498  ArrayRCP<Scalar> interfaceData = matvecInterfaceTmp->getDataNonConst(0);
499  for (LocalOrdinal interfaceIdx = 0; interfaceIdx < static_cast<LocalOrdinal>(interfaceData.size()); ++interfaceIdx) {
500  YData[regionInterfaceLIDs[interfaceIdx]] += interfaceData[interfaceIdx];
501  }
502  }
503 }
504 
505 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
506 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getDomainMap() const {
507  XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap");
508  return toXpetra(mtx_->getDomainMap());
509 }
510 
511 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
512 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getRangeMap() const {
513  XPETRA_MONITOR("TpetraCrsMatrix::getRangeMap");
514  return toXpetra(mtx_->getRangeMap());
515 }
516 
517 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
519  XPETRA_MONITOR("TpetraCrsMatrix::description");
520  return mtx_->description();
521 }
522 
523 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const {
525  XPETRA_MONITOR("TpetraCrsMatrix::describe");
526  mtx_->describe(out, verbLevel);
527 }
528 
529 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
531  XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
532  Teuchos::LabeledObject::setObjectLabel(objectLabel);
533  mtx_->setObjectLabel(objectLabel);
534 }
535 
536 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
538  : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(*(matrix.mtx_), Teuchos::Copy))) {}
539 
540 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
542  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
543  XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
544  mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
545 }
546 
547 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
549  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
550  mtx_->getLocalDiagOffsets(offsets);
551 }
552 
553 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
554 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagCopy(Vector &diag, const Teuchos::ArrayView<const size_t> &offsets) const {
555  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
556  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
557 }
558 
559 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
560 void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagCopy(Vector &diag, const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const {
561  XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
562  mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
563 }
564 
565 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
567  XPETRA_MONITOR("TpetraCrsMatrix::replaceDiag");
568  Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(toTpetra(diag)));
569 }
570 
571 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
573  XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
574  mtx_->leftScale(*(toTpetra(x)));
575 }
576 
577 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
579  XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
580  mtx_->rightScale(*(toTpetra(x)));
581 }
582 
583 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
584 Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getMap() const {
585  XPETRA_MONITOR("TpetraCrsMatrix::getMap");
586  return rcp(new TpetraMap<LocalOrdinal, GlobalOrdinal, Node>(mtx_->getMap()));
587 }
588 
589 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
592  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
593 
594  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
595  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSource.getTpetra_CrsMatrix();
596  // mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
597  mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
598 }
599 
601 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
604  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
605 
606  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
607  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tDest.getTpetra_CrsMatrix();
608  mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
609 }
610 
611 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
614  XPETRA_MONITOR("TpetraCrsMatrix::doImport");
615 
616  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
617  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSource.getTpetra_CrsMatrix();
618  mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
619 }
620 
621 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
624  XPETRA_MONITOR("TpetraCrsMatrix::doExport");
625 
626  XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
627  RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tDest.getTpetra_CrsMatrix();
628  mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
629 }
630 
631 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
633  XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
634  mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
635 }
636 
637 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
639  return !mtx_.is_null();
640 }
641 
642 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
643 TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &mtx)
644  : mtx_(mtx) {}
645 
646 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
647 RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getTpetra_CrsMatrix() const { return mtx_; }
648 
650 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
651 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getTpetra_CrsMatrixNonConst() const { return mtx_; } // TODO: remove
652 
654 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
656  const MultiVector &B,
657  MultiVector &R) const {
658  Tpetra::Details::residual(*mtx_, toTpetra(X), toTpetra(B), toTpetra(R));
659 }
660 
663 // End of TpetrCrsMatrix class definition //
666 
667 } // namespace Xpetra
668 
669 #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)