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#include "Xpetra_TpetraBlockCrsMatrix.hpp"
29#include "Xpetra_TpetraCrsMatrix.hpp"
30
32#include <Teuchos_Hashtable.hpp>
33
38namespace Xpetra {
39
40typedef std::string viewLabel_t;
41
46template <class Scalar,
47 class LocalOrdinal,
48 class GlobalOrdinal,
49 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
50class CrsMatrixWrap : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
55#ifdef HAVE_XPETRA_TPETRA
57#endif
60#ifdef HAVE_XPETRA_TPETRA
62#endif
63
64 public:
66
67
69 CrsMatrixWrap(const RCP<const Map> &rowMap);
70
72 CrsMatrixWrap(const RCP<const Map> &rowMap,
73 size_t maxNumEntriesPerRow);
74
76 CrsMatrixWrap(const RCP<const Map> &rowMap,
77 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
78
80 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow);
81
83 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
84
85#ifdef HAVE_XPETRA_TPETRA
87 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params = null);
88
90 CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
91 const RCP<const Map> &domainMap = Teuchos::null, const RCP<const Map> &rangeMap = Teuchos::null,
92 const Teuchos::RCP<Teuchos::ParameterList> &params = null);
93#else
94#ifdef __GNUC__
95#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."
96#endif
97#endif
98
100
101 CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList = Teuchos::null);
102
105 const RCP<ParameterList> &paramList = Teuchos::null);
106
108 virtual ~CrsMatrixWrap();
109
111
113
114
116
127 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
128
130
137 void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
138
140
145 void replaceGlobalValues(GlobalOrdinal globalRow,
147 const ArrayView<const Scalar> &vals);
148
150
153 void replaceLocalValues(LocalOrdinal localRow,
155 const ArrayView<const Scalar> &vals);
156
158 virtual void setAllToScalar(const Scalar &alpha);
159
161 void scale(const Scalar &alpha);
162
164
166
167
176 void resumeFill(const RCP<ParameterList> &params = null);
177
189 void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
190
204 // TODO : Get ride of "Tpetra"::OptimizeOption
205 void fillComplete(const RCP<ParameterList> &params = null);
206
208
210
213
215
218
220 size_t getLocalNumRows() const;
221
224
226 size_t getLocalNumEntries() const;
227
229
230 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
231
233
234 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
235
237
239 size_t getGlobalMaxNumRowEntries() const;
240
242
244 size_t getLocalMaxNumRowEntries() const;
245
247 bool isLocallyIndexed() const;
248
250 bool isGloballyIndexed() const;
251
253 bool isFillComplete() const;
254
256
268 void getLocalRowCopy(LocalOrdinal LocalRow,
269 const ArrayView<LocalOrdinal> &Indices,
270 const ArrayView<Scalar> &Values,
271 size_t &NumEntries) const;
272
274
283 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
284
286
295 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
296
298
301
304
307
310
313
316
318 bool haveGlobalConstants() const;
319
321
323
324
326
335 // TODO virtual=0 // TODO: Add default parameters ?
336 // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
337 // matrixData_->multiply(X, Y, trans, alpha, beta);
338 // }
339
341
343
344
346
352 Scalar alpha = ScalarTraits<Scalar>::one(),
353 Scalar beta = ScalarTraits<Scalar>::zero()) const;
354
358 Teuchos::ETransp mode,
359 Scalar alpha,
360 Scalar beta,
361 bool sumInterfaceValues,
362 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter,
363 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const;
364
368
372
375 const RCP<const Map> &getColMap() const;
376
378 const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
379
381
383
385 //{@
386
389
391 void doImport(const Matrix &source,
393
395 void doExport(const Matrix &dest,
397
399 void doImport(const Matrix &source,
401
403 void doExport(const Matrix &dest,
405
406 // @}
407
409
410
412 std::string description() const;
413
416
418
419 void setObjectLabel(const std::string &objectLabel);
421
422#ifdef HAVE_XPETRA_TPETRA
424#if KOKKOS_VERSION >= 40799
425 virtual typename local_matrix_type::host_mirror_type getLocalMatrixHost() const;
426#else
427 virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const;
428#endif
429#else
430#ifdef __GNUC__
431#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."
432#endif
433#endif
434
435 // JG: Added:
436
437 bool hasCrsGraph() const;
438
441
443
445 LocalOrdinal GetStorageBlockSize() const;
446
451
454
456 private:
457 // Default view is created after fillComplete()
458 // Because ColMap might not be available before fillComplete().
459 void CreateDefaultView();
460
461 // 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.
462 // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
463 void updateDefaultView() const;
464 // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
465 // See also CrsMatrixWrap::updateDefaultView()
466 mutable bool finalDefaultView_;
467
468 // The underlying matrix object
470
471}; // class CrsMatrixWrap
472
473template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
476 return Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
477}
478
479template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
482 return Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
483}
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>
499toXpetra(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
501 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
503}
504
505template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
507toXpetra(Teuchos::RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
509 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
511}
512
513template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
515toXpetra(Teuchos::RCP<Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
517 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
519}
520
521template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
523toXpetra(Teuchos::RCP<const Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
525 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
527}
528
529template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
534
535template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
540
541template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
542Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
546
547template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
548Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
552
553template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
556 return Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrix();
557}
558
559template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
562 return Teuchos::rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrixNonConst();
563}
564
565template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
566Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
570
571template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
572const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
576
577template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
582
583template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
588
589template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
590Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
594
595template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
596const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
600
601template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
604 auto crsMat = toCrsMatrix(mat);
605 auto tmp_Crs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
606 if (!tmp_Crs.is_null()) {
607 return tmp_Crs->getTpetra_CrsMatrixNonConst();
608 } else {
609 auto tmp_BlockCrs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
610 if (tmp_BlockCrs.is_null())
611 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
612 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
613 }
614}
615
616template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
619 auto crsMat = toCrsMatrix(mat);
620 auto tmp_Crs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
621 if (!tmp_Crs.is_null()) {
622 return tmp_Crs->getTpetra_CrsMatrix();
623 } else {
624 auto tmp_BlockCrs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
625 if (tmp_BlockCrs.is_null())
626 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
627 return tmp_BlockCrs->getTpetra_BlockCrsMatrix();
628 }
629}
630
631template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
634 auto mat = Teuchos::rcp_dynamic_cast<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
635 auto rmat = Teuchos::rcp_dynamic_cast<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
636 if (!mat.is_null()) {
637 return toTpetraRowMatrix(mat);
638 } else if (!rmat.is_null()) {
639 return rmat->getTpetra_RowMatrixNonConst();
640 } else {
641 auto tpOp = Teuchos::rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true);
642 auto tOp = tpOp->getOperator();
643 auto tRow = Teuchos::rcp_dynamic_cast<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
644 return tRow;
645 }
646}
647
648template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
651 auto mat = Teuchos::rcp_dynamic_cast<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
652 auto rmat = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
653 if (!mat.is_null()) {
654 return toTpetraRowMatrix(mat);
655 } else if (!rmat.is_null()) {
656 return rmat->getTpetra_RowMatrix();
657 } else {
658 auto tpOp = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true);
659 auto tOp = tpOp->getOperator();
660 auto tRow = Teuchos::rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
661 return tRow;
662 }
663}
664
665} // namespace Xpetra
666
667#define XPETRA_CRSMATRIXWRAP_SHORT
668#endif // XPETRA_CRSMATRIXWRAP_DECL_HPP
669
670// NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
671// 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.
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
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
Exception indicating invalid cast attempted.
Xpetra-specific matrix class.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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 CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
Teuchos::RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > toTpetraRowMatrix(const Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mat)
CombineMode
Xpetra::Combine Mode enumerable type.