12#ifndef XPETRA_CRSMATRIXWRAP_DECL_HPP
13#define XPETRA_CRSMATRIXWRAP_DECL_HPP
15#include <Tpetra_KokkosCompat_DefaultNode.hpp>
25#include "Xpetra_Matrix.hpp"
27#include "Xpetra_TpetraBlockCrsMatrix.hpp"
29#include "Xpetra_TpetraCrsMatrix.hpp"
34#include <Tpetra_Filter.hpp>
48template <
class Scalar,
51 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
57#ifdef HAVE_XPETRA_TPETRA
62#ifdef HAVE_XPETRA_TPETRA
77 size_t maxNumEntriesPerRow);
89#ifdef HAVE_XPETRA_TPETRA
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."
165 void scale(
const Scalar &alpha);
275 size_t &NumEntries)
const;
365 bool sumInterfaceValues,
426#ifdef HAVE_XPETRA_TPETRA
428#if KOKKOS_VERSION >= 40799
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."
477template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
480 return Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(matrix,
true)->getCrsMatrix();
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();
489template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
495template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
501template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
505 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs,
true);
509template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
513 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs,
true);
517template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
521 auto xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs,
true);
525template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
529 auto xpCrs = Teuchos::rcp_dynamic_cast<const Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(xpTpCrs,
true);
533template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
539template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
545template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
546Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
551template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
552Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
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();
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();
569template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
570Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
575template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
576const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
581template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
587template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
593template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
594Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
599template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
600const Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
605template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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();
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();
620template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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();
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();
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()) {
642 }
else if (!rmat.is_null()) {
643 return rmat->getTpetra_RowMatrixNonConst();
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);
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()) {
659 }
else if (!rmat.is_null()) {
660 return rmat->getTpetra_RowMatrix();
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);
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);
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);
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();
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);
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);
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();
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);
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);
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();
755#define XPETRA_CRSMATRIXWRAP_SHORT
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 > ¶ms=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.
void updateDefaultView() const
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.
RCP< CrsMatrix > matrixData_
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 > ¶ms=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.