Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_CrsMatrixWrap_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 // WARNING: This code is experimental. Backwards compatibility should not be expected.
11 
12 #ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
13 #define XPETRA_CRSMATRIXWRAP_DEF_HPP
14 
15 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
16 
17 #include <Teuchos_SerialDenseMatrix.hpp>
18 #include <Teuchos_Hashtable.hpp>
19 
20 #include "Xpetra_ConfigDefs.hpp"
21 #include "Xpetra_Exceptions.hpp"
22 
23 #include "Xpetra_MultiVector.hpp"
24 #include "Xpetra_CrsGraph.hpp"
25 #include "Xpetra_CrsMatrix.hpp"
27 
28 #include "Xpetra_Matrix.hpp"
29 
31 
32 namespace Xpetra {
33 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35  : finalDefaultView_(false) {
38 }
39 
40 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42  size_t maxNumEntriesPerRow)
43  : finalDefaultView_(false) {
44  matrixData_ = CrsMatrixFactory::Build(rowMap, maxNumEntriesPerRow);
46 }
47 
48 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
50  const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc)
51  : finalDefaultView_(false) {
52  matrixData_ = CrsMatrixFactory::Build(rowMap, NumEntriesPerRowToAlloc);
54 }
55 
56 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57 CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow)
58  : finalDefaultView_(false) {
59  // Set matrix data
60  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, maxNumEntriesPerRow);
61 
62  // Default view
64 }
65 
66 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67 CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc)
68  : finalDefaultView_(false) {
69  // Set matrix data
70  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, NumEntriesPerRowToAlloc);
71 
72  // Default view
74 }
75 
76 #ifdef HAVE_XPETRA_TPETRA
77 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78 CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params)
79  : finalDefaultView_(false) {
80  // Set matrix data
81  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
82 
83  // Default view
85 }
86 
87 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
88 CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
89  const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap,
90  const Teuchos::RCP<Teuchos::ParameterList> &params)
91  : finalDefaultView_(false) {
92  // Set matrix data
93  matrixData_ = CrsMatrixFactory::Build(lclMatrix, rowMap, colMap, domainMap, rangeMap, params);
94 
95  // Default view
97 }
98 #else
99 #ifdef __GNUC__
100 #warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
101 #endif
102 #endif
103 
104 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
106  : finalDefaultView_(matrix->isFillComplete()) {
107  // Set matrix data
108  matrixData_ = matrix;
109 
110  // Default view
112 }
113 
114 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115 CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList)
116  : finalDefaultView_(false) {
117  // Set matrix data
118  matrixData_ = CrsMatrixFactory::Build(graph, paramList);
119 
120  // Default view
122 }
123 
124 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
127 
128  const RCP<ParameterList> &paramList)
129  : finalDefaultView_(false) {
130  // Set matrix data
131  matrixData_ = CrsMatrixFactory::Build(graph, values, paramList);
132 
133  // Default view
135 }
136 
137 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139 
140 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
141 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
142  matrixData_->insertGlobalValues(globalRow, cols, vals);
143 }
144 
145 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
146 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
147  matrixData_->insertLocalValues(localRow, cols, vals);
148 }
149 
150 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152  const ArrayView<const GlobalOrdinal> &cols,
153  const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
154 
155 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
157  const ArrayView<const LocalOrdinal> &cols,
158  const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
159 
160 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
161 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
162 
163 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
165  matrixData_->scale(alpha);
166 }
167 
168 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170  matrixData_->resumeFill(params);
171 }
172 
173 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params) {
175  matrixData_->fillComplete(domainMap, rangeMap, params);
176 
177  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
178  updateDefaultView();
179 }
180 
181 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
183  matrixData_->fillComplete(params);
184 
185  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
186  updateDefaultView();
187 }
188 
189 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191  return matrixData_->getGlobalNumRows();
192 }
193 
194 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196  return matrixData_->getGlobalNumCols();
197 }
198 
199 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201  return matrixData_->getLocalNumRows();
202 }
203 
204 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
206  return matrixData_->getGlobalNumEntries();
207 }
208 
209 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
211  return matrixData_->getLocalNumEntries();
212 }
213 
214 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
216  return matrixData_->getNumEntriesInLocalRow(localRow);
217 }
218 
219 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221  return matrixData_->getNumEntriesInGlobalRow(globalRow);
222 }
223 
224 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226  return matrixData_->getGlobalMaxNumRowEntries();
227 }
228 
229 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231  return matrixData_->getLocalMaxNumRowEntries();
232 }
233 
234 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236  return matrixData_->isLocallyIndexed();
237 }
238 
239 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241  return matrixData_->isGloballyIndexed();
242 }
243 
244 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246  return matrixData_->isFillComplete();
247 }
248 
249 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
251  const ArrayView<LocalOrdinal> &Indices,
252  const ArrayView<Scalar> &Values,
253  size_t &NumEntries) const {
254  matrixData_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
255 }
256 
257 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const {
259  matrixData_->getGlobalRowView(GlobalRow, indices, values);
260 }
261 
262 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
263 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const {
264  matrixData_->getLocalRowView(LocalRow, indices, values);
265 }
266 
267 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
269  matrixData_->getLocalDiagCopy(diag);
270 }
271 
272 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
274  matrixData_->getLocalDiagOffsets(offsets);
275 }
276 
277 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
279  matrixData_->getLocalDiagCopy(diag, offsets);
280 }
281 
282 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
283 typename ScalarTraits<Scalar>::magnitudeType CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getFrobeniusNorm() const {
284  return matrixData_->getFrobeniusNorm();
285 }
286 
287 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289  matrixData_->leftScale(x);
290 }
291 
292 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
294  matrixData_->rightScale(x);
295 }
296 
297 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
299  return matrixData_->haveGlobalConstants();
300 }
301 
302 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305  Teuchos::ETransp mode,
306  Scalar alpha,
307  Scalar beta) const {
308  matrixData_->apply(X, Y, mode, alpha, beta);
309 }
310 
311 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
314  Teuchos::ETransp mode,
315  Scalar alpha,
316  Scalar beta,
317  bool sumInterfaceValues,
318  const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter,
319  const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const {
320  matrixData_->apply(X, Y, mode, alpha, beta, sumInterfaceValues, regionInterfaceImporter, regionInterfaceLIDs);
321 }
322 
323 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
324 const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getDomainMap() const {
325  return matrixData_->getDomainMap();
326 }
327 
328 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329 const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getRangeMap() const {
330  return matrixData_->getRangeMap();
331 }
332 
333 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
334 const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> &CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getColMap() const { return getColMap(Matrix::GetCurrentViewLabel()); }
335 
336 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
337 const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> &CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getColMap(viewLabel_t viewLabel) const {
338  TEUCHOS_TEST_FOR_EXCEPTION(Matrix::operatorViewTable_.containsKey(viewLabel) == false, Xpetra::Exceptions::RuntimeError, "Xpetra::Matrix.GetColMap(): view '" + viewLabel + "' does not exist.");
339  updateDefaultView(); // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated.
340  return Matrix::operatorViewTable_.get(viewLabel)->GetColMap();
341 }
342 
343 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
345  matrixData_->removeEmptyProcessesInPlace(newMap);
346  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetRowMap(matrixData_->getRowMap());
347  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetColMap(matrixData_->getColMap());
348 }
349 
350 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
351 const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getMap() const {
352  return matrixData_->getMap();
353 }
354 
355 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
358  const CrsMatrixWrap &sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
359  matrixData_->doImport(*sourceWrp.getCrsMatrix(), importer, CM);
360 }
361 
362 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
365  const CrsMatrixWrap &destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
366  matrixData_->doExport(*destWrp.getCrsMatrix(), importer, CM);
367 }
368 
369 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372  const CrsMatrixWrap &sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
373  matrixData_->doImport(*sourceWrp.getCrsMatrix(), exporter, CM);
374 }
375 
376 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
379  const CrsMatrixWrap &destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
380  matrixData_->doExport(*destWrp.getCrsMatrix(), exporter, CM);
381 }
382 
383 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
385  return "Xpetra::CrsMatrixWrap";
386 }
387 
388 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
389 void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const {
390  // Teuchos::EVerbosityLevel vl = verbLevel;
391  // if (vl == VERB_DEFAULT) vl = VERB_LOW;
392  // RCP<const Comm<int> > comm = this->getComm();
393  // const int myImageID = comm->getRank(),
394  // numImages = comm->getSize();
395 
396  // if (myImageID == 0) out << this->description() << std::endl;
397 
398  matrixData_->describe(out, verbLevel);
399 
400  // Teuchos::OSTab tab(out);
401 }
402 
403 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
405  Teuchos::LabeledObject::setObjectLabel(objectLabel);
406  matrixData_->setObjectLabel(objectLabel);
407 }
408 
409 #ifdef HAVE_XPETRA_TPETRA
410 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
413  return matrixData_->getLocalMatrixHost();
414 }
415 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
418  return matrixData_->getLocalMatrixDevice();
419 }
420 #else
421 #ifdef __GNUC__
422 #warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
423 #endif
424 #endif
425 
426 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
428 
429 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
430 RCP<const Xpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getCrsGraph() const { return matrixData_->getCrsGraph(); }
431 
432 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
433 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getCrsMatrix() const { return matrixData_; }
434 
435 // Default view is created after fillComplete()
436 // Because ColMap might not be available before fillComplete().
437 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
439  // Create default view
440  this->defaultViewLabel_ = "point";
441  this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
442 
443  // Set current view
444  this->currentViewLabel_ = this->GetDefaultViewLabel();
445 }
446 
447 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
449  if ((finalDefaultView_ == false) && matrixData_->isFillComplete()) {
450  // Update default view with the colMap
451  Matrix::operatorViewTable_.get(Matrix::GetDefaultViewLabel())->SetColMap(matrixData_->getColMap());
452  finalDefaultView_ = true;
453  }
454 }
455 
456 // Expert only
457 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
459  // Clear the old view table
460  Teuchos::Hashtable<viewLabel_t, RCP<MatrixView>> dummy_table;
461  Matrix::operatorViewTable_ = dummy_table;
462 
463  finalDefaultView_ = M->isFillComplete();
464  // Set matrix data
465  matrixData_ = M;
466 
467  // Default view
468  CreateDefaultView();
469 }
470 
471 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
473  return matrixData_->GetStorageBlockSize();
474 }
475 
476 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
481  matrixData_->residual(X, B, R);
482 }
483 
484 } // namespace Xpetra
485 
486 #endif // ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
void setObjectLabel(const std::string &objectLabel)
RCP< CrsMatrix > getCrsMatrix() const
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
const viewLabel_t & GetCurrentViewLabel() const
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
virtual local_matrix_type getLocalMatrixDevice() const
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
static RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap)
Constructor for empty matrix (intended use is an import/export target - can&#39;t insert entries directly...
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
Exception throws to report errors in the internal logical of the program.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
const viewLabel_t & GetDefaultViewLabel() const
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
std::string description() const
Return a simple one-line description of this object.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
CrsMatrix::local_matrix_type local_matrix_type
void getLocalDiagCopy(Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
CrsMatrixWrap(const RCP< const Map > &rowMap)
Constructor for a dynamic profile matrix (Epetra only)
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
size_t global_size_t
Global size_t object.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< Teuchos::ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
virtual ~CrsMatrixWrap()
Destructor.
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
void resumeFill(const RCP< ParameterList > &params=null)
void replaceCrsMatrix(RCP< CrsMatrix > &M)
Expert only.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Concrete implementation of Xpetra::Matrix.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
CombineMode
Xpetra::Combine Mode enumerable type.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
virtual void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
virtual local_matrix_type::HostMirror getLocalMatrixHost() const
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
Xpetra-specific matrix class.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...