Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_CrsMatrixWrap_decl.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_DECL_HPP
13 #define XPETRA_CRSMATRIXWRAP_DECL_HPP
14 
15 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
16 
17 #include "Xpetra_ConfigDefs.hpp"
18 #include "Xpetra_Exceptions.hpp"
19 
21 #include "Xpetra_CrsGraph.hpp"
22 #include "Xpetra_CrsMatrix.hpp"
24 
25 #include "Xpetra_Matrix.hpp"
26 
27 #include <Teuchos_SerialDenseMatrix.hpp>
28 #include <Teuchos_Hashtable.hpp>
29 
34 namespace Xpetra {
35 
36 typedef std::string viewLabel_t;
37 
42 template <class Scalar,
43  class LocalOrdinal,
44  class GlobalOrdinal,
45  class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
46 class CrsMatrixWrap : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
51 #ifdef HAVE_XPETRA_TPETRA
53 #endif
56 #ifdef HAVE_XPETRA_TPETRA
58 #endif
59 
60  public:
62 
63 
65  CrsMatrixWrap(const RCP<const Map> &rowMap);
66 
68  CrsMatrixWrap(const RCP<const Map> &rowMap,
69  size_t maxNumEntriesPerRow);
70 
72  CrsMatrixWrap(const RCP<const Map> &rowMap,
73  const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
74 
76  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow);
77 
79  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
80 
81 #ifdef HAVE_XPETRA_TPETRA
82  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params = null);
84 
86  CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
87  const RCP<const Map> &domainMap = Teuchos::null, const RCP<const Map> &rangeMap = Teuchos::null,
88  const Teuchos::RCP<Teuchos::ParameterList> &params = null);
89 #else
90 #ifdef __GNUC__
91 #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."
92 #endif
93 #endif
94 
95  CrsMatrixWrap(RCP<CrsMatrix> matrix);
96 
97  CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList = Teuchos::null);
98 
99  CrsMatrixWrap(const RCP<const CrsGraph> &graph,
101  const RCP<ParameterList> &paramList = Teuchos::null);
102 
104  virtual ~CrsMatrixWrap();
105 
107 
109 
110 
112 
123  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
124 
126 
133  void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
134 
136 
141  void replaceGlobalValues(GlobalOrdinal globalRow,
142  const ArrayView<const GlobalOrdinal> &cols,
143  const ArrayView<const Scalar> &vals);
144 
146 
149  void replaceLocalValues(LocalOrdinal localRow,
150  const ArrayView<const LocalOrdinal> &cols,
151  const ArrayView<const Scalar> &vals);
152 
154  virtual void setAllToScalar(const Scalar &alpha);
155 
157  void scale(const Scalar &alpha);
158 
160 
162 
163 
172  void resumeFill(const RCP<ParameterList> &params = null);
173 
185  void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
186 
200  // TODO : Get ride of "Tpetra"::OptimizeOption
201  void fillComplete(const RCP<ParameterList> &params = null);
202 
204 
206 
209 
211 
214 
216  size_t getLocalNumRows() const;
217 
220 
222  size_t getLocalNumEntries() const;
223 
225 
226  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
227 
229 
230  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
231 
233 
235  size_t getGlobalMaxNumRowEntries() const;
236 
238 
240  size_t getLocalMaxNumRowEntries() const;
241 
243  bool isLocallyIndexed() const;
244 
246  bool isGloballyIndexed() const;
247 
249  bool isFillComplete() const;
250 
252 
264  void getLocalRowCopy(LocalOrdinal LocalRow,
265  const ArrayView<LocalOrdinal> &Indices,
266  const ArrayView<Scalar> &Values,
267  size_t &NumEntries) const;
268 
270 
279  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
280 
282 
291  void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
292 
294 
297 
299  void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const;
300 
302  void getLocalDiagCopy(Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &diag, const Teuchos::ArrayView<const size_t> &offsets) const;
303 
305  typename ScalarTraits<Scalar>::magnitudeType getFrobeniusNorm() const;
306 
309 
312 
314  bool haveGlobalConstants() const;
315 
317 
319 
320 
322 
331  // TODO virtual=0 // TODO: Add default parameters ?
332  // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
333  // matrixData_->multiply(X, Y, trans, alpha, beta);
334  // }
335 
337 
339 
340 
342 
347  Teuchos::ETransp mode = Teuchos::NO_TRANS,
348  Scalar alpha = ScalarTraits<Scalar>::one(),
349  Scalar beta = ScalarTraits<Scalar>::zero()) const;
350 
354  Teuchos::ETransp mode,
355  Scalar alpha,
356  Scalar beta,
357  bool sumInterfaceValues,
358  const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> > &regionInterfaceImporter,
359  const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const;
360 
363  const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getDomainMap() const;
364 
367  const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getRangeMap() const;
368 
371  const RCP<const Map> &getColMap() const;
372 
374  const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
375 
376  void removeEmptyProcessesInPlace(const Teuchos::RCP<const Map> &newMap);
377 
379 
381  //{@
382 
384  const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getMap() const;
385 
387  void doImport(const Matrix &source,
389 
391  void doExport(const Matrix &dest,
393 
395  void doImport(const Matrix &source,
397 
399  void doExport(const Matrix &dest,
401 
402  // @}
403 
405 
406 
408  std::string description() const;
409 
411  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const;
412 
414 
415  void setObjectLabel(const std::string &objectLabel);
417 
418 #ifdef HAVE_XPETRA_TPETRA
419  virtual local_matrix_type getLocalMatrixDevice() const;
420  virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const;
421 #else
422 #ifdef __GNUC__
423 #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."
424 #endif
425 #endif
426 
427  // JG: Added:
428 
429  bool hasCrsGraph() const;
430 
432  RCP<const CrsGraph> getCrsGraph() const;
433 
434  RCP<CrsMatrix> getCrsMatrix() const;
435 
437  LocalOrdinal GetStorageBlockSize() const;
438 
443 
445  void replaceCrsMatrix(RCP<CrsMatrix> &M);
446 
448  private:
449  // Default view is created after fillComplete()
450  // Because ColMap might not be available before fillComplete().
451  void CreateDefaultView();
452 
453  private:
454  // The colMap can be <tt>null</tt> until fillComplete() is called. The default view of the Matrix have to be updated when fillComplete() is called.
455  // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
456  void updateDefaultView() const;
457  // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
458  // See also CrsMatrixWrap::updateDefaultView()
459  mutable bool finalDefaultView_;
460 
461  // The underlying matrix object
462  RCP<CrsMatrix> matrixData_;
463 
464 }; // class CrsMatrixWrap
465 
466 } // namespace Xpetra
467 
468 #define XPETRA_CRSMATRIXWRAP_SHORT
469 #endif // XPETRA_CRSMATRIXWRAP_DECL_HPP
470 
471 // NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
472 // TODO: getUnderlyingMatrix() method
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.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Xpetra::TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrix
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.
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...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
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.
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.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
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.
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
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.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
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.
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
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.
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
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...