12 #ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
13 #define XPETRA_CRSMATRIXWRAP_DEF_HPP
15 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
17 #include <Teuchos_SerialDenseMatrix.hpp>
18 #include <Teuchos_Hashtable.hpp>
23 #include "Xpetra_MultiVector.hpp"
28 #include "Xpetra_Matrix.hpp"
33 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
35 : finalDefaultView_(false) {
40 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
42 size_t maxNumEntriesPerRow)
43 : finalDefaultView_(false) {
48 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
50 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc)
51 : finalDefaultView_(false) {
56 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
58 : finalDefaultView_(false) {
66 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 : finalDefaultView_(false) {
76 #ifdef HAVE_XPETRA_TPETRA
77 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
79 : finalDefaultView_(false) {
87 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
89 const RCP<const Map> &domainMap,
const RCP<const Map> &rangeMap,
90 const Teuchos::RCP<Teuchos::ParameterList> ¶ms)
91 : finalDefaultView_(false) {
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."
104 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
106 : finalDefaultView_(matrix->isFillComplete()) {
114 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
116 : finalDefaultView_(false) {
124 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 const RCP<ParameterList> ¶mList)
129 : finalDefaultView_(false) {
137 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
142 matrixData_->insertGlobalValues(globalRow, cols, vals);
145 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
147 matrixData_->insertLocalValues(localRow, cols, vals);
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); }
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); }
160 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
163 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
165 matrixData_->scale(alpha);
168 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
170 matrixData_->resumeFill(params);
173 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
175 matrixData_->fillComplete(domainMap, rangeMap, params);
181 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
183 matrixData_->fillComplete(params);
189 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
191 return matrixData_->getGlobalNumRows();
194 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
196 return matrixData_->getGlobalNumCols();
199 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
201 return matrixData_->getLocalNumRows();
204 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 return matrixData_->getGlobalNumEntries();
209 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
211 return matrixData_->getLocalNumEntries();
214 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
216 return matrixData_->getNumEntriesInLocalRow(localRow);
219 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
221 return matrixData_->getNumEntriesInGlobalRow(globalRow);
224 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
226 return matrixData_->getGlobalMaxNumRowEntries();
229 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
231 return matrixData_->getLocalMaxNumRowEntries();
234 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236 return matrixData_->isLocallyIndexed();
239 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 return matrixData_->isGloballyIndexed();
244 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
246 return matrixData_->isFillComplete();
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);
257 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
259 matrixData_->getGlobalRowView(GlobalRow, indices, values);
262 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
264 matrixData_->getLocalRowView(LocalRow, indices, values);
267 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
269 matrixData_->getLocalDiagCopy(diag);
272 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
274 matrixData_->getLocalDiagOffsets(offsets);
277 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
279 matrixData_->getLocalDiagCopy(diag, offsets);
282 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
284 return matrixData_->getFrobeniusNorm();
287 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
289 matrixData_->leftScale(x);
292 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
294 matrixData_->rightScale(x);
297 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
299 return matrixData_->haveGlobalConstants();
302 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
305 Teuchos::ETransp mode,
308 matrixData_->apply(X, Y, mode, alpha, beta);
311 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
314 Teuchos::ETransp mode,
317 bool sumInterfaceValues,
319 const Teuchos::ArrayRCP<LocalOrdinal> ®ionInterfaceLIDs)
const {
320 matrixData_->apply(X, Y, mode, alpha, beta, sumInterfaceValues, regionInterfaceImporter, regionInterfaceLIDs);
323 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
325 return matrixData_->getDomainMap();
328 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
330 return matrixData_->getRangeMap();
333 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
336 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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());
350 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
352 return matrixData_->getMap();
355 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
362 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
369 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
376 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
383 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
385 return "Xpetra::CrsMatrixWrap";
388 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
398 matrixData_->describe(out, verbLevel);
403 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
405 Teuchos::LabeledObject::setObjectLabel(objectLabel);
406 matrixData_->setObjectLabel(objectLabel);
409 #ifdef HAVE_XPETRA_TPETRA
410 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
413 return matrixData_->getLocalMatrixHost();
415 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
418 return matrixData_->getLocalMatrixDevice();
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."
426 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
429 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
432 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
437 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
440 this->defaultViewLabel_ =
"point";
441 this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
444 this->currentViewLabel_ = this->GetDefaultViewLabel();
447 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
449 if ((finalDefaultView_ ==
false) && matrixData_->isFillComplete()) {
452 finalDefaultView_ =
true;
457 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
460 Teuchos::Hashtable<viewLabel_t, RCP<MatrixView>> dummy_table;
463 finalDefaultView_ = M->isFillComplete();
471 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
473 return matrixData_->GetStorageBlockSize();
476 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
481 matrixData_->residual(X, B, R);
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'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_
void updateDefaultView() const
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 > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps.
RCP< CrsMatrix > matrixData_
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 > ¶ms=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...