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
34#include <Tpetra_Filter.hpp>
35
40namespace Xpetra {
41
42typedef std::string viewLabel_t;
43
48template <class Scalar,
49 class LocalOrdinal,
50 class GlobalOrdinal,
51 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
52class CrsMatrixWrap : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
57#ifdef HAVE_XPETRA_TPETRA
59#endif
62#ifdef HAVE_XPETRA_TPETRA
66#endif
67
68 public:
70
71
73 CrsMatrixWrap(const RCP<const Map> &rowMap);
74
76 CrsMatrixWrap(const RCP<const Map> &rowMap,
77 size_t maxNumEntriesPerRow);
78
80 CrsMatrixWrap(const RCP<const Map> &rowMap,
81 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
82
84 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow);
85
87 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
88
89#ifdef HAVE_XPETRA_TPETRA
91 CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params = null);
92
94 CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
95 const RCP<const Map> &domainMap = Teuchos::null, const RCP<const Map> &rangeMap = Teuchos::null,
96 const Teuchos::RCP<Teuchos::ParameterList> &params = null);
97#else
98#ifdef __GNUC__
99#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."
100#endif
101#endif
102
104
105 CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList = Teuchos::null);
106
109 const RCP<ParameterList> &paramList = Teuchos::null);
110
112 virtual ~CrsMatrixWrap();
113
115
117
118
120
131 void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
132
134
141 void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
142
144
149 void replaceGlobalValues(GlobalOrdinal globalRow,
151 const ArrayView<const Scalar> &vals);
152
154
157 void replaceLocalValues(LocalOrdinal localRow,
159 const ArrayView<const Scalar> &vals);
160
162 virtual void setAllToScalar(const Scalar &alpha);
163
165 void scale(const Scalar &alpha);
166
168
170
171
180 void resumeFill(const RCP<ParameterList> &params = null);
181
193 void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
194
208 // TODO : Get ride of "Tpetra"::OptimizeOption
209 void fillComplete(const RCP<ParameterList> &params = null);
210
212
214
217
219
222
224 size_t getLocalNumRows() const;
225
228
230 size_t getLocalNumEntries() const;
231
233
234 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
235
237
238 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
239
241
243 size_t getGlobalMaxNumRowEntries() const;
244
246
248 size_t getLocalMaxNumRowEntries() const;
249
251 bool isLocallyIndexed() const;
252
254 bool isGloballyIndexed() const;
255
257 bool isFillComplete() const;
258
260
272 void getLocalRowCopy(LocalOrdinal LocalRow,
273 const ArrayView<LocalOrdinal> &Indices,
274 const ArrayView<Scalar> &Values,
275 size_t &NumEntries) const;
276
278
287 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
288
290
299 void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
300
302
305
308
311
314
317
320
322 bool haveGlobalConstants() const;
323
325
327
328
330
339 // TODO virtual=0 // TODO: Add default parameters ?
340 // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
341 // matrixData_->multiply(X, Y, trans, alpha, beta);
342 // }
343
345
347
348
350
356 Scalar alpha = ScalarTraits<Scalar>::one(),
357 Scalar beta = ScalarTraits<Scalar>::zero()) const;
358
362 Teuchos::ETransp mode,
363 Scalar alpha,
364 Scalar beta,
365 bool sumInterfaceValues,
366 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter,
367 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const;
368
372
376
379 const RCP<const Map> &getColMap() const;
380
382 const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
383
385
387
389 //{@
390
393
395 void doImport(const Matrix &source,
397
399 void doExport(const Matrix &dest,
401
403 void doImport(const Matrix &source,
405
407 void doExport(const Matrix &dest,
409
410 // @}
411
413
414
416 std::string description() const;
417
420
422
423 void setObjectLabel(const std::string &objectLabel);
425
426#ifdef HAVE_XPETRA_TPETRA
428#if KOKKOS_VERSION >= 40799
429 virtual typename local_matrix_type::host_mirror_type getLocalMatrixHost() const;
430#else
431 virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const;
432#endif
433#else
434#ifdef __GNUC__
435#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."
436#endif
437#endif
438
439 // JG: Added:
440
441 bool hasCrsGraph() const;
442
445
447
449 LocalOrdinal GetStorageBlockSize() const;
450
455
458
460 private:
461 // Default view is created after fillComplete()
462 // Because ColMap might not be available before fillComplete().
463 void CreateDefaultView();
464
465 // 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.
466 // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
467 void updateDefaultView() const;
468 // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
469 // See also CrsMatrixWrap::updateDefaultView()
470 mutable bool finalDefaultView_;
471
472 // The underlying matrix object
474
475}; // class CrsMatrixWrap
476
477template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
480 return Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
481}
482
483template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
486 return Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix, true)->getCrsMatrix();
487}
488
489template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
494
495template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
500
501template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
503toXpetra(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
505 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
507}
508
509template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
511toXpetra(Teuchos::RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
513 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
515}
516
517template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
519toXpetra(Teuchos::RCP<Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
521 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
523}
524
525template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
527toXpetra(Teuchos::RCP<const Tpetra::FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &A) {
529 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs, true);
531}
532
533template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
538
539template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
544
545template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
546Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
550
551template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
552Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
556
557template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
560 return Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrix();
561}
562
563template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
566 return Teuchos::rcp_dynamic_cast<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true)->getTpetra_BlockCrsMatrixNonConst();
567}
568
569template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
570Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
574
575template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
576const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
580
581template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
586
587template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
592
593template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
594Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
598
599template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
600const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
604
605template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
608 auto crsMat = toCrsMatrix(mat);
609 auto tmp_Crs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
610 if (!tmp_Crs.is_null()) {
611 return tmp_Crs->getTpetra_CrsMatrixNonConst();
612 } else {
613 auto tmp_BlockCrs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
614 if (tmp_BlockCrs.is_null())
615 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
616 return tmp_BlockCrs->getTpetra_BlockCrsMatrixNonConst();
617 }
618}
619
620template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
623 auto crsMat = toCrsMatrix(mat);
624 auto tmp_Crs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
625 if (!tmp_Crs.is_null()) {
626 return tmp_Crs->getTpetra_CrsMatrix();
627 } else {
628 auto tmp_BlockCrs = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsMat);
629 if (tmp_BlockCrs.is_null())
630 throw Exceptions::BadCast("Cast from Xpetra::CrsMatrix to Xpetra::TpetraCrsMatrix and Xpetra::TpetraBlockCrsMatrix failed");
631 return tmp_BlockCrs->getTpetra_BlockCrsMatrix();
632 }
633}
634
635template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
638 auto mat = Teuchos::rcp_dynamic_cast<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
639 auto rmat = Teuchos::rcp_dynamic_cast<Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
640 if (!mat.is_null()) {
641 return toTpetraRowMatrix(mat);
642 } else if (!rmat.is_null()) {
643 return rmat->getTpetra_RowMatrixNonConst();
644 } else {
645 auto tpOp = Teuchos::rcp_dynamic_cast<Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true);
646 auto tOp = tpOp->getOperator();
647 auto tRow = Teuchos::rcp_dynamic_cast<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
648 return tRow;
649 }
650}
651
652template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
655 auto mat = Teuchos::rcp_dynamic_cast<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
656 auto rmat = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraRowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A);
657 if (!mat.is_null()) {
658 return toTpetraRowMatrix(mat);
659 } else if (!rmat.is_null()) {
660 return rmat->getTpetra_RowMatrix();
661 } else {
662 auto tpOp = Teuchos::rcp_dynamic_cast<const Xpetra::TpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A, true);
663 auto tOp = tpOp->getOperator();
664 auto tRow = Teuchos::rcp_dynamic_cast<const Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tOp, true);
665 return tRow;
666 }
667}
668
669template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
672#ifdef HAVE_XPETRA_TPETRA
673 auto tpCrsA = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsA);
674 if (!tpCrsA.is_null()) {
675 auto tpCrsFiltered = Tpetra::applyFilter_vals(*tpCrsA, filter);
676 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tpCrsFiltered);
677 }
678#endif
679#ifdef HAVE_XPETRA_EPETRA
680 auto epCrsA = Teuchos::rcp_dynamic_cast<Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>(crsA);
681 if (!epCrsA.is_null()) {
682 auto epCrsFiltered = Tpetra::applyFilter_vals(*epCrsA, filter);
683 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(epCrsFiltered);
684 }
685#endif
686 return crsFiltered;
687}
688
689template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
691 auto crsA = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A)->getCrsMatrix();
692 auto crsFiltered = applyFilter_vals(crsA, filter);
693 auto filtered = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsFiltered));
694 return filtered;
695}
696
697template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
700#ifdef HAVE_XPETRA_TPETRA
701 auto tpCrsA = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsA);
702 if (!tpCrsA.is_null()) {
703 auto tpCrsFiltered = Tpetra::applyFilter_LID(*tpCrsA, filter);
704 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tpCrsFiltered);
705 }
706#endif
707#ifdef HAVE_XPETRA_EPETRA
708 auto epCrsA = Teuchos::rcp_dynamic_cast<Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>(crsA);
709 if (!epCrsA.is_null()) {
710 auto epCrsFiltered = Tpetra::applyFilter_LID(*epCrsA, filter);
711 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(epCrsFiltered);
712 }
713#endif
714 return crsFiltered;
715}
716
717template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
719 auto crsA = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A)->getCrsMatrix();
720 auto crsFiltered = applyFilter_LID(crsA, filter);
721 auto filtered = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsFiltered));
722 return filtered;
723}
724
725template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
728#ifdef HAVE_XPETRA_TPETRA
729 auto tpCrsA = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(crsA);
730 if (!tpCrsA.is_null()) {
731 auto tpCrsFiltered = Tpetra::applyFilter_GID(*tpCrsA, filter);
732 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(tpCrsFiltered);
733 }
734#endif
735#ifdef HAVE_XPETRA_EPETRA
736 auto epCrsA = Teuchos::rcp_dynamic_cast<Xpetra::EpetraCrsMatrixT<GlobalOrdinal, Node>>(crsA);
737 if (!epCrsA.is_null()) {
738 auto epCrsFiltered = Tpetra::applyFilter_GID(*epCrsA, filter);
739 crsFiltered = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(epCrsFiltered);
740 }
741#endif
742 return crsFiltered;
743}
744
745template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, class filter_type>
747 auto crsA = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(A)->getCrsMatrix();
748 auto crsFiltered = applyFilter_GID(crsA, filter);
749 auto filtered = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(crsFiltered));
750 return filtered;
751}
752
753} // namespace Xpetra
754
755#define XPETRA_CRSMATRIXWRAP_SHORT
756#endif // XPETRA_CRSMATRIXWRAP_DECL_HPP
757
758// NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
759// 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.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_host_type local_matrix_host_type
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
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_device_type
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.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
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.
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...
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_device_type
typename local_matrix_device_type::HostMirror local_matrix_host_type
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.
RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > applyFilter_vals(const RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &crsA, const filter_type &filter)
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< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > applyFilter_LID(const RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &crsA, const filter_type &filter)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > applyFilter_GID(const RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &crsA, const filter_type &filter)
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.