10 #ifndef XPETRA_TPETRAROWMATRIX_HPP 
   11 #define XPETRA_TPETRAROWMATRIX_HPP 
   17 #include <Teuchos_Describable.hpp> 
   19 #include "Xpetra_Map.hpp" 
   22 #include "Xpetra_TpetraMap.hpp" 
   24 #include "Tpetra_RowMatrix.hpp" 
   31 template <
class Scalar,
 
   34           class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
 
   36   : 
virtual public RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
 
   50   const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > 
getRowMap()
 const {
 
   56   const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > 
getColMap()
 const {
 
   64     return mtx_->getGlobalNumRows();
 
   70     return mtx_->getGlobalNumCols();
 
   76     return mtx_->getLocalNumRows();
 
   82     return mtx_->getLocalNumCols();
 
   88     return mtx_->getGlobalNumEntries();
 
   94     return mtx_->getNodeNumEntries();
 
  100     return mtx_->getNumEntriesInLocalRow(localRow);
 
  106     return mtx_->getGlobalMaxNumRowEntries();
 
  112     return mtx_->getLocalMaxNumRowEntries();
 
  118     return mtx_->isLocallyIndexed();
 
  124     return mtx_->isGloballyIndexed();
 
  130     return mtx_->isFillComplete();
 
  136     return mtx_->supportsRowViews();
 
  145   void getLocalRowCopy(LocalOrdinal LocalRow, 
const Teuchos::ArrayView<LocalOrdinal> &Indices, 
const Teuchos::ArrayView<Scalar> &Values, 
size_t &NumEntries)
 const {
 
  147     typename Tpetra::RowGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_local_inds_host_view_type indices(
"indices", Indices.size());
 
  148     typename Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values(
"values", Values.size());
 
  150     mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
 
  151     for (
size_t i = 0; i < NumEntries; ++i) {
 
  152       Indices[i] = indices(i);
 
  153       Values[i]  = values(i);
 
  158   void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values)
 const {
 
  160     typename Tpetra::RowGraph<LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type k_indices;
 
  161     typename Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type k_values;
 
  163     mtx_->getGlobalRowView(GlobalRow, k_indices, k_values);
 
  164     indices = ArrayView<const GlobalOrdinal>(k_indices.data(), k_indices.extent(0));
 
  165     values  = ArrayView<const Scalar>(
reinterpret_cast<const Scalar *
>(k_values.data()), k_values.extent(0));
 
  169   void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values)
 const {
 
  171     typename Tpetra::RowGraph<LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type k_indices;
 
  172     typename Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type k_values;
 
  174     mtx_->getLocalRowView(LocalRow, k_indices, k_values);
 
  175     indices = ArrayView<const LocalOrdinal>(k_indices.data(), k_indices.extent(0));
 
  176     values  = ArrayView<const Scalar>(
reinterpret_cast<const Scalar *
>(k_values.data()), k_values.extent(0));
 
  191     return mtx_->getFrobeniusNorm();
 
  200   const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > 
getDomainMap()
 const {
 
  206   const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > 
getRangeMap()
 const {
 
  212   void apply(
const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, 
MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode = Teuchos::NO_TRANS, Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(), Scalar beta = Teuchos::ScalarTraits<Scalar>::zero())
 const { TEUCHOS_TEST_FOR_EXCEPTION(1, 
Xpetra::Exceptions::NotImplemented, 
"TODO"); }
 
  223   TpetraRowMatrix(
const Teuchos::RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx)
 
  227   void setTpetra_RowMatrix(
const Teuchos::RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) { 
mtx_ = mtx; }
 
  238   RCP<Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > 
mtx_;
 
  244 #define XPETRA_TPETRAROWMATRIX_SHORT 
  245 #endif  // XPETRA_TPETRAROWMATRIX_HPP 
Tpetra::global_size_t getGlobalNumCols() const 
Returns the number of global columns in this matrix. 
RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const 
Returns the Map associated with the range of this operator, which must be compatible with Y...
void setTpetra_RowMatrix(const Teuchos::RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
Set the underlying Tpetra matrix. 
Tpetra::global_size_t getGlobalNumEntries() const 
Returns the global number of entries in this matrix. 
TpetraRowMatrix()=default
size_t getLocalNumRows() const 
Returns the number of rows owned on the calling node. 
virtual ~TpetraRowMatrix()=default
Destructor. 
bool isFillComplete() const 
Returns true if fillComplete() has been called. 
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const 
Returns the Map that describes the column distribution in this matrix. 
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const 
Returns the Map that describes the row distribution in this matrix. 
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const 
Returns the Frobenius norm of the matrix. 
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const 
Computes the operator-multivector application. 
Exception throws when you call an unimplemented method of Xpetra. 
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const 
Get a copy of the diagonal entries owned by this node, with local row indices. 
size_t global_size_t
Global size_t object. 
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const 
Returns the current number of entries on this node in the specified local row. 
bool isLocallyIndexed() const 
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const 
Returns the Map associated with the domain of this operator, which must be compatible with X...
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. 
Tpetra::global_size_t getGlobalNumRows() const 
Returns the number of global rows in this matrix. 
TpetraRowMatrix(const Teuchos::RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object. 
size_t getLocalNumCols() const 
Returns the number of columns needed to apply the forward operator on this node, i.e., the number of elements listed in the column map. 
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. 
RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_RowMatrix() const 
Get the underlying Tpetra matrix. 
#define XPETRA_MONITOR(funcName)
RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_RowMatrixNonConst() const 
Get the underlying Tpetra matrix. 
bool supportsRowViews() const 
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class. 
size_t getLocalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on this node. 
size_t getGlobalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on all nodes. 
size_t getNodeNumEntries() const 
Returns the local number of entries in this matrix. 
void getLocalRowCopy(LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const 
Extract a list of entries in a specified local row of the graph. Put into storage allocated by callin...
bool isGloballyIndexed() const 
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)