Xpetra Version of the Day
Loading...
Searching...
No Matches
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
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
32namespace Xpetra {
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
39
40template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42 size_t maxNumEntriesPerRow)
43 : finalDefaultView_(false) {
44 matrixData_ = CrsMatrixFactory::Build(rowMap, maxNumEntriesPerRow);
46}
47
48template <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
56template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
58 : finalDefaultView_(false) {
59 // Set matrix data
60 matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, maxNumEntriesPerRow);
61
62 // Default view
64}
65
66template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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
77template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79 : finalDefaultView_(false) {
80 // Set matrix data
81 matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
82
83 // Default view
85}
86
87template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
89 const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap,
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
104template <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
114template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116 : finalDefaultView_(false) {
117 // Set matrix data
118 matrixData_ = CrsMatrixFactory::Build(graph, paramList);
119
120 // Default view
122}
123
124template <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
137template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139
140template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
142 matrixData_->insertGlobalValues(globalRow, cols, vals);
143}
144
145template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
147 matrixData_->insertLocalValues(localRow, cols, vals);
148}
149
150template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
153 const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
154
155template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
158 const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
159
160template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
161void CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>::setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
162
163template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
165 matrixData_->scale(alpha);
166}
167
168template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170 matrixData_->resumeFill(params);
171}
172
173template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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
181template <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
189template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191 return matrixData_->getGlobalNumRows();
192}
193
194template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196 return matrixData_->getGlobalNumCols();
197}
198
199template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
201 return matrixData_->getLocalNumRows();
202}
203
204template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
206 return matrixData_->getGlobalNumEntries();
207}
208
209template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
211 return matrixData_->getLocalNumEntries();
212}
213
214template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
216 return matrixData_->getNumEntriesInLocalRow(localRow);
217}
218
219template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221 return matrixData_->getNumEntriesInGlobalRow(globalRow);
222}
223
224template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226 return matrixData_->getGlobalMaxNumRowEntries();
227}
228
229template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231 return matrixData_->getLocalMaxNumRowEntries();
232}
233
234template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236 return matrixData_->isLocallyIndexed();
237}
238
239template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241 return matrixData_->isGloballyIndexed();
242}
243
244template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246 return matrixData_->isFillComplete();
247}
248
249template <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
257template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
259 matrixData_->getGlobalRowView(GlobalRow, indices, values);
260}
261
262template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
264 matrixData_->getLocalRowView(LocalRow, indices, values);
265}
266
267template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271
272template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
274 matrixData_->getLocalDiagOffsets(offsets);
275}
276
277template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281
282template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
286
287template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
291
292template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
296
297template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
299 return matrixData_->haveGlobalConstants();
300}
301
302template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
310
311template <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
323template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
327
328template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
332
333template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
335
336template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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
343template <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
350template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354
355template <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
362template <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
369template <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
376template <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
383template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
385 return "Xpetra::CrsMatrixWrap";
386}
387
388template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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
403template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
406 matrixData_->setObjectLabel(objectLabel);
407}
408
409#ifdef HAVE_XPETRA_TPETRA
410template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
411#if KOKKOS_VERSION >= 40799
413#else
415#endif
419template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
422 return matrixData_->getLocalMatrixDevice();
423}
424#else
425#ifdef __GNUC__
426#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."
427#endif
428#endif
429
430template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
432
433template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
435
436template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
438
439// Default view is created after fillComplete()
440// Because ColMap might not be available before fillComplete().
441template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
443 // Create default view
444 this->defaultViewLabel_ = "point";
445 this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
446
447 // Set current view
448 this->currentViewLabel_ = this->GetDefaultViewLabel();
449}
450
451template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
453 if ((finalDefaultView_ == false) && matrixData_->isFillComplete()) {
454 // Update default view with the colMap
455 Matrix::operatorViewTable_.get(Matrix::GetDefaultViewLabel())->SetColMap(matrixData_->getColMap());
456 finalDefaultView_ = true;
457 }
458}
459
460// Expert only
461template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
463 // Clear the old view table
465 Matrix::operatorViewTable_ = dummy_table;
466
467 finalDefaultView_ = M->isFillComplete();
468 // Set matrix data
469 matrixData_ = M;
470
471 // Default view
472 CreateDefaultView();
473}
474
475template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
477 return matrixData_->GetStorageBlockSize();
478}
479
480template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
487
488} // namespace Xpetra
489
490#endif // ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
virtual void setObjectLabel(const std::string &objectLabel)
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't insert entries directly...
Concrete implementation of Xpetra::Matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
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.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
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.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
std::string description() const
Return a simple one-line description of this object.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void resumeFill(const RCP< ParameterList > &params=null)
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
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.
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
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...
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.
void setObjectLabel(const std::string &objectLabel)
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
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...
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
virtual local_matrix_type getLocalMatrixDevice() const
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void replaceCrsMatrix(RCP< CrsMatrix > &M)
Expert only.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
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.
virtual ~CrsMatrixWrap()
Destructor.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
virtual local_matrix_type::HostMirror getLocalMatrixHost() const
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.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
CrsMatrixWrap(const RCP< const Map > &rowMap)
Constructor for a dynamic profile matrix (Epetra only)
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
RCP< CrsMatrix > getCrsMatrix() const
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.
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.
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...
Exception throws to report errors in the internal logical of the program.
Xpetra-specific matrix class.
const viewLabel_t & GetCurrentViewLabel() const
const viewLabel_t & GetDefaultViewLabel() const
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.