Xpetra Version of the Day
Loading...
Searching...
No Matches
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"
27
29#include <Teuchos_Hashtable.hpp>
30
35namespace Xpetra {
36
37typedef std::string viewLabel_t;
38
43template <class Scalar,
44 class LocalOrdinal,
45 class GlobalOrdinal,
46 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
47class CrsMatrixWrap : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
52#ifdef HAVE_XPETRA_TPETRA
54#endif
57#ifdef HAVE_XPETRA_TPETRA
59#endif
60
61 public:
63
64
66 CrsMatrixWrap(const RCP<const Map> &rowMap);
67
69 CrsMatrixWrap(const RCP<const Map> &rowMap,
70 size_t maxNumEntriesPerRow);
71
73 CrsMatrixWrap(const RCP<const Map> &rowMap,
74 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
75
77 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow);
78
80 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
81
82#ifdef HAVE_XPETRA_TPETRA
84 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params = null);
85
87 CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
88 const RCP<const Map> &domainMap = Teuchos::null, const RCP<const Map> &rangeMap = Teuchos::null,
89 const Teuchos::RCP<Teuchos::ParameterList> &params = null);
90#else
91#ifdef __GNUC__
92#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."
93#endif
94#endif
95
97
98 CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList = Teuchos::null);
99
102 const RCP<ParameterList> &paramList = Teuchos::null);
103
105 virtual ~CrsMatrixWrap();
106
108
110
111
113
124 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
125
127
134 void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
135
137
142 void replaceGlobalValues(GlobalOrdinal globalRow,
144 const ArrayView<const Scalar> &vals);
145
147
150 void replaceLocalValues(LocalOrdinal localRow,
152 const ArrayView<const Scalar> &vals);
153
155 virtual void setAllToScalar(const Scalar &alpha);
156
158 void scale(const Scalar &alpha);
159
161
163
164
173 void resumeFill(const RCP<ParameterList> &params = null);
174
186 void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
187
201 // TODO : Get ride of "Tpetra"::OptimizeOption
202 void fillComplete(const RCP<ParameterList> &params = null);
203
205
207
210
212
215
217 size_t getLocalNumRows() const;
218
221
223 size_t getLocalNumEntries() const;
224
226
227 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
228
230
231 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
232
234
236 size_t getGlobalMaxNumRowEntries() const;
237
239
241 size_t getLocalMaxNumRowEntries() const;
242
244 bool isLocallyIndexed() const;
245
247 bool isGloballyIndexed() const;
248
250 bool isFillComplete() const;
251
253
265 void getLocalRowCopy(LocalOrdinal LocalRow,
266 const ArrayView<LocalOrdinal> &Indices,
267 const ArrayView<Scalar> &Values,
268 size_t &NumEntries) const;
269
271
280 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
281
283
292 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
293
295
298
301
304
307
310
313
315 bool haveGlobalConstants() const;
316
318
320
321
323
332 // TODO virtual=0 // TODO: Add default parameters ?
333 // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
334 // matrixData_->multiply(X, Y, trans, alpha, beta);
335 // }
336
338
340
341
343
349 Scalar alpha = ScalarTraits<Scalar>::one(),
350 Scalar beta = ScalarTraits<Scalar>::zero()) const;
351
355 Teuchos::ETransp mode,
356 Scalar alpha,
357 Scalar beta,
358 bool sumInterfaceValues,
359 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter,
360 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const;
361
365
369
372 const RCP<const Map> &getColMap() const;
373
375 const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
376
378
380
382 //{@
383
386
388 void doImport(const Matrix &source,
390
392 void doExport(const Matrix &dest,
394
396 void doImport(const Matrix &source,
398
400 void doExport(const Matrix &dest,
402
403 // @}
404
406
407
409 std::string description() const;
410
413
415
416 void setObjectLabel(const std::string &objectLabel);
418
419#ifdef HAVE_XPETRA_TPETRA
421 virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const;
422#else
423#ifdef __GNUC__
424#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."
425#endif
426#endif
427
428 // JG: Added:
429
430 bool hasCrsGraph() const;
431
434
436
438 LocalOrdinal GetStorageBlockSize() const;
439
444
447
449 private:
450 // Default view is created after fillComplete()
451 // Because ColMap might not be available before fillComplete().
452 void CreateDefaultView();
453
454 private:
455 // 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.
456 // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
457 void updateDefaultView() const;
458 // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
459 // See also CrsMatrixWrap::updateDefaultView()
460 mutable bool finalDefaultView_;
461
462 // The underlying matrix object
464
465}; // class CrsMatrixWrap
466
467template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
470 return Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
471}
472
473template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
476 return Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
477}
478
479template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
484
485template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
490
491template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
496
497template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
502
503template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
504Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
508
509template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
510Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
514
515template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
518 return Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrix();
519}
520
521template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524 return Teuchos::rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrixNonConst();
525}
526
527template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
528Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
532
533template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
534const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
538
539template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
544
545template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
550
551template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
552Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
556
557template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
558const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
562
563} // namespace Xpetra
564
565#define XPETRA_CRSMATRIXWRAP_SHORT
566#endif // XPETRA_CRSMATRIXWRAP_DECL_HPP
567
568// NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
569// TODO: getUnderlyingMatrix() method
static const EVerbosityLevel verbLevel_default
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...
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
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.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
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.
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
Xpetra::TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrix
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
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.
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
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.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
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
Xpetra-specific matrix class.
RCP< const Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrix() const
Get the underlying Tpetra matrix.
RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrixNonConst() const
Get the underlying Tpetra matrix.
basic_FancyOStream< char > FancyOStream
Teuchos::RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > toCrsMatrix(const Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &matrix)
size_t global_size_t
Global size_t object.
Teuchos::RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > toTpetraBlock(const Teuchos::RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.