42 #ifndef TPETRA_ROWMATRIX_DECL_HPP
43 #define TPETRA_ROWMATRIX_DECL_HPP
45 #include "Tpetra_ConfigDefs.hpp"
48 #include "Tpetra_Operator.hpp"
52 #include "Teuchos_Describable.hpp"
80 template <
class Scalar,
85 virtual public Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
87 public Packable<char, LocalOrdinal> {
120 virtual Teuchos::RCP<const Teuchos::Comm<int> >
getComm()
const = 0;
122 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
123 virtual TPETRA_DEPRECATED Teuchos::RCP<Node> getNode()
const = 0;
125 #endif // TPETRA_ENABLE_DEPRECATED_CODE
128 virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >
getRowMap()
const = 0;
131 virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >
getColMap()
const = 0;
134 virtual Teuchos::RCP<const RowGraph<LocalOrdinal,GlobalOrdinal,Node> >
getGraph()
const = 0;
234 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
253 virtual size_t TPETRA_DEPRECATED getNodeNumDiags ()
const = 0;
265 virtual bool TPETRA_DEPRECATED isLowerTriangular ()
const = 0;
277 virtual bool TPETRA_DEPRECATED isUpperTriangular ()
const = 0;
278 #endif // TPETRA_ENABLE_DEPRECATED_CODE
306 const Teuchos::ArrayView<GlobalOrdinal> &Indices,
307 const Teuchos::ArrayView<Scalar> &Values,
308 size_t &NumEntries)
const = 0;
332 const Teuchos::ArrayView<LocalOrdinal> &Indices,
333 const Teuchos::ArrayView<Scalar> &Values,
334 size_t &NumEntries)
const = 0;
362 Teuchos::ArrayView<const GlobalOrdinal> &indices,
363 Teuchos::ArrayView<const Scalar> &values)
const = 0;
391 Teuchos::ArrayView<const LocalOrdinal>& indices,
392 Teuchos::ArrayView<const Scalar>& values)
const = 0;
422 LocalOrdinal& numEnt,
423 const LocalOrdinal*& lclColInds,
424 const Scalar*& vals)
const;
516 virtual Teuchos::RCP<RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
517 add (
const Scalar& alpha,
522 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null)
const;
528 packRow (
char*
const numEntOut,
532 const LocalOrdinal lclRow)
const;
537 allocatePackSpace (Teuchos::Array<char>& exports,
538 size_t& totalNumEntries,
539 const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs)
const;
546 packImpl (
const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
547 Teuchos::Array<char>& exports,
548 const Teuchos::ArrayView<size_t>& numPacketsPerLID,
549 size_t& constantNumPackets,
563 pack (
const Teuchos::ArrayView<const LocalOrdinal>& exportLIDs,
564 Teuchos::Array<char>& exports,
565 const Teuchos::ArrayView<size_t>& numPacketsPerLID,
566 size_t& constantNumPackets,
572 #endif // TPETRA_ROWMATRIX_DECL_HPP
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const =0
The communicator over which this matrix is distributed.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
The Map that describes the distribution of columns over processes.
virtual bool isLocallyIndexed() const =0
Whether matrix indices are locally indexed.
virtual size_t getNodeMaxNumRowEntries() const =0
Maximum number of entries in any row of the matrix, on this process.
virtual global_size_t getGlobalNumCols() const =0
The global number of columns of this matrix.
virtual GlobalOrdinal getIndexBase() const =0
The index base for global indices in this matrix.
Node node_type
The Kokkos Node type.
virtual size_t getNodeNumRows() const =0
The number of rows owned by the calling process.
virtual void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Scale the matrix on the right with the given Vector.
virtual size_t getNodeNumCols() const =0
The number of columns needed to apply the forward operator on this node.
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Get a copy of the given local row's entries.
Forward declaration of Tpetra::RowMatrix.
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
The current number of entries on the calling process in the specified local row.
virtual void getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
Get a constant, nonpersisting, locally indexed view of the given row of the matrix.
virtual void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Scale the matrix on the left with the given Vector.
Abstract base class for sources of an Import or Export.
virtual Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const =0
The RowGraph associated with this matrix.
LocalOrdinal local_ordinal_type
The type of local indices.
virtual LocalOrdinal getLocalRowViewRaw(const LocalOrdinal lclRow, LocalOrdinal &numEnt, const LocalOrdinal *&lclColInds, const Scalar *&vals) const
Get a constant, nonpersisting, locally indexed view of the given row of the matrix, using "raw" pointers instead of Teuchos::ArrayView.
virtual ~RowMatrix()
Destructor (virtual for memory safety of derived classes).
Forward declaration of Tpetra::RowGraph.
size_t global_size_t
Global size_t object.
virtual void getGlobalRowView(GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const =0
Get a constant, nonpersisting, globally indexed view of the given row of the matrix.
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Get a copy of the given global row's entries.
virtual size_t getGlobalMaxNumRowEntries() const =0
Maximum number of entries in any row of the matrix, over all processes.
Abstract interface for operators (e.g., matrices and preconditioners).
GlobalOrdinal global_ordinal_type
The type of global indices.
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
The current number of entries on the calling process in the specified global row. ...
Sets up and executes a communication plan for a Tpetra DistObject.
virtual void pack(const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< char > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor) const
Pack this object's data for an Import or Export.
CrsGraphType::global_ordinal_type getGlobalNumDiags(const CrsGraphType &G)
Number of populated diagonal entries in the given sparse graph, over all processes in the graph's (MP...
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
Abstract base class for objects that can be the source of an Import or Export operation.
Forward declaration of Tpetra::Vector.
virtual void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0
Get a copy of the diagonal entries, distributed by the row Map.
virtual bool isFillComplete() const =0
Whether fillComplete() has been called.
virtual global_size_t getGlobalNumEntries() const =0
The global number of stored (structurally nonzero) entries.
virtual Teuchos::RCP< RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > add(const Scalar &alpha, const RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null) const
Return a new RowMatrix which is the result of beta*this + alpha*A.
Scalar scalar_type
The type of the entries in the sparse matrix.
A parallel distribution of indices over processes.
virtual bool supportsRowViews() const =0
Whether this object implements getLocalRowView() and getGlobalRowView().
A read-only, row-oriented interface to a sparse matrix.
A distributed dense vector.
virtual size_t getNodeNumEntries() const =0
The local number of stored (structurally nonzero) entries.
virtual bool isGloballyIndexed() const =0
Whether matrix indices are globally indexed.
Declaration of Tpetra::Packable.
Abstract base class for objects that can be the source of an Import or Export operation, and that also know how to pack their data to send to the target object.
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type mag_type
Type of a norm result.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
The Map that describes the distribution of rows over processes.
virtual bool hasColMap() const =0
Whether this matrix has a well-defined column Map.
virtual global_size_t getGlobalNumRows() const =0
The global number of rows of this matrix.
virtual mag_type getFrobeniusNorm() const =0
The Frobenius norm of the matrix.