Tpetra parallel linear algebra
Version of the Day
|
A read-only, row-oriented interface to a sparse matrix. More...
#include <Tpetra_RowMatrix_decl.hpp>
Public Types | |
Typedefs | |
typedef Scalar | scalar_type |
The type of the entries in the sparse matrix. More... | |
typedef LocalOrdinal | local_ordinal_type |
The type of local indices. More... | |
typedef GlobalOrdinal | global_ordinal_type |
The type of global indices. More... | |
typedef Node | node_type |
The Kokkos Node type. More... | |
typedef MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type | mag_type |
Type of a norm result. More... | |
Public Member Functions | |
Destructor | |
virtual | ~RowMatrix () |
Destructor (virtual for memory safety of derived classes). More... | |
Matrix query methods | |
virtual Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const =0 |
The communicator over which this matrix is distributed. More... | |
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > | getRowMap () const =0 |
The Map that describes the distribution of rows over processes. More... | |
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > | getColMap () const =0 |
The Map that describes the distribution of columns over processes. More... | |
virtual Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > | getGraph () const =0 |
The RowGraph associated with this matrix. More... | |
virtual global_size_t | getGlobalNumRows () const =0 |
The global number of rows of this matrix. More... | |
virtual global_size_t | getGlobalNumCols () const =0 |
The global number of columns of this matrix. More... | |
virtual size_t | getNodeNumRows () const =0 |
The number of rows owned by the calling process. More... | |
virtual size_t | getNodeNumCols () const =0 |
The number of columns needed to apply the forward operator on this node. More... | |
virtual GlobalOrdinal | getIndexBase () const =0 |
The index base for global indices in this matrix. More... | |
virtual global_size_t | getGlobalNumEntries () const =0 |
The global number of stored (structurally nonzero) entries. More... | |
virtual size_t | getNodeNumEntries () const =0 |
The local number of stored (structurally nonzero) entries. More... | |
virtual size_t | getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const =0 |
The current number of entries on the calling process in the specified global row. More... | |
virtual size_t | getNumEntriesInLocalRow (LocalOrdinal localRow) const =0 |
The current number of entries on the calling process in the specified local row. More... | |
virtual size_t | getGlobalMaxNumRowEntries () const =0 |
Maximum number of entries in any row of the matrix, over all processes. More... | |
virtual size_t | getNodeMaxNumRowEntries () const =0 |
Maximum number of entries in any row of the matrix, on this process. More... | |
virtual bool | hasColMap () const =0 |
Whether this matrix has a well-defined column Map. More... | |
virtual bool | isLocallyIndexed () const =0 |
Whether matrix indices are locally indexed. More... | |
virtual bool | isGloballyIndexed () const =0 |
Whether matrix indices are globally indexed. More... | |
virtual bool | isFillComplete () const =0 |
Whether fillComplete() has been called. More... | |
virtual bool | supportsRowViews () const =0 |
Whether this object implements getLocalRowView() and getGlobalRowView(). More... | |
Extraction Methods | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual void | getLocalDiagCopy (Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0 |
Get a copy of the diagonal entries, distributed by the row Map. More... | |
Mathematical methods | |
virtual void | leftScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0 |
Scale the matrix on the left with the given Vector. More... | |
virtual void | rightScale (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0 |
Scale the matrix on the right with the given Vector. More... | |
virtual mag_type | getFrobeniusNorm () const =0 |
The Frobenius norm of the matrix. More... | |
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 . More... | |
Pure virtual functions to be overridden by subclasses. | |
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > | getDomainMap () const =0 |
The Map associated with the domain of this operator, which must be compatible with X.getMap(). More... | |
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > | getRangeMap () const =0 |
The Map associated with the range of this operator, which must be compatible with Y.getMap(). More... | |
virtual 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 =0 |
Computes the operator-multivector application. More... | |
virtual bool | hasTransposeApply () const |
Whether this operator supports applying the transpose or conjugate transpose. More... | |
Implementation of Packable interface | |
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. More... | |
A read-only, row-oriented interface to a sparse matrix.
Scalar | The type of the entries in the sparse matrix. Same as the Scalar typedef of Operator. |
LocalOrdinal | The type of local indices. See the documentation of Map for requirements. |
GlobalOrdinal | The type of global indices. See the documentation of Map for requirements. |
Node | The Kokkos Node type. See the documentation of Map for requirements. |
RowMatrix provides a read-only, row-oriented interface to view the entries of a sparse matrix, which is distributed over one or more (MPI) processes. "Read only" means that you may view the entries, but not modify them. "Row-oriented" means that you may view all the entries owned by the calling process, in any particular row of the matrix that is owned by the calling process.
CrsMatrix implements this interface, but also lets you add or modify entries of the sparse matrix. Ifpack2 provides other implementations of RowMatrix, which do useful things like wrapping an existing matrix to view only certain desired entries.
Definition at line 84 of file Tpetra_RowMatrix_decl.hpp.
typedef Scalar Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type |
The type of the entries in the sparse matrix.
Definition at line 93 of file Tpetra_RowMatrix_decl.hpp.
typedef LocalOrdinal Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type |
The type of local indices.
Definition at line 95 of file Tpetra_RowMatrix_decl.hpp.
typedef GlobalOrdinal Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type |
The type of global indices.
Definition at line 97 of file Tpetra_RowMatrix_decl.hpp.
typedef Node Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type |
The Kokkos Node type.
Definition at line 99 of file Tpetra_RowMatrix_decl.hpp.
typedef MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::mag_type Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type |
Type of a norm result.
This is usually the same as the type of the magnitude (absolute value) of Scalar
, but may differ for certain Scalar
types.
Definition at line 106 of file Tpetra_RowMatrix_decl.hpp.
|
virtual |
Destructor (virtual for memory safety of derived classes).
Definition at line 52 of file Tpetra_RowMatrix_def.hpp.
|
pure virtual |
The communicator over which this matrix is distributed.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The Map that describes the distribution of rows over processes.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The Map that describes the distribution of columns over processes.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The RowGraph associated with this matrix.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The global number of rows of this matrix.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The global number of columns of this matrix.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The number of rows owned by the calling process.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The number of columns needed to apply the forward operator on this node.
This is the same as the number of elements listed in the column Map. It is not necessarily the same as the number of domain Map elements owned by the calling process.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The index base for global indices in this matrix.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The global number of stored (structurally nonzero) entries.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The local number of stored (structurally nonzero) entries.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The current number of entries on the calling process in the specified global row.
Note that if the row Map is overlapping, then the calling process might not necessarily store all the entries in the row. Some other process might have the rest of the entries.
Teuchos::OrdinalTraits<size_t>::invalid()
if the specified global row does not belong to this graph, else the number of entries. Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
The current number of entries on the calling process in the specified local row.
Note that if the row Map is overlapping, then the calling process might not necessarily store all the entries in the row. Some other process might have the rest of the entries.
Teuchos::OrdinalTraits<size_t>::invalid()
if the specified local row is not valid for this graph, else the number of entries. Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Maximum number of entries in any row of the matrix, over all processes.
This method only uses the matrix's graph. Explicitly stored zeros count as "entries."
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Maximum number of entries in any row of the matrix, on this process.
This method only uses the matrix's graph. Explicitly stored zeros count as "entries."
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Whether this matrix has a well-defined column Map.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Whether matrix indices are locally indexed.
A RowMatrix may store column indices either as global indices (of type GlobalOrdinal
), or as local indices (of type LocalOrdinal
). In some cases (for example, if the column Map has not been computed), it is not possible to switch from global to local indices without extra work. Furthermore, some operations only work for one or the other case.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Whether matrix indices are globally indexed.
A RowMatrix may store column indices either as global indices (of type GlobalOrdinal
), or as local indices (of type LocalOrdinal
). In some cases (for example, if the column Map has not been computed), it is not possible to switch from global to local indices without extra work. Furthermore, some operations only work for one or the other case.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Whether fillComplete() has been called.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Whether this object implements getLocalRowView() and getGlobalRowView().
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Get a copy of the given global row's entries.
This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.
GlobalRow | [in] Global index of the row. |
Indices | [out] Global indices of the columns corresponding to values. |
Values | [out] Matrix values. |
NumEntries | [out] Number of stored entries on the calling process; length of Indices and Values. |
This method throws std::runtime_error
if either Indices or Values is not large enough to hold the data associated with row GlobalRow. If GlobalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid()
, and does not modify Indices or Values.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Get a copy of the given local row's entries.
This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.
LocalRow | [in] Local index of the row. |
Indices | [out] Local indices of the columns corresponding to values. |
Values | [out] Matrix values. |
NumEntries | [out] Number of stored entries on the calling process; length of Indices and Values. |
This method throws std::runtime_error
if either Indices or Values is not large enough to hold the data associated with row LocalRow. If LocalRow does not belong to the calling process, then the method sets NumIndices to Teuchos::OrdinalTraits<size_t>::invalid()
, and does not modify Indices or Values.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Get a constant, nonpersisting, globally indexed view of the given row of the matrix.
The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this
. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.
This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.
isGloballyIndexed () && supportsRowViews ()
indices.size () == getNumEntriesInGlobalRow (GlobalRow)
GlobalRow | [in] Global index of the row. |
Indices | [out] Global indices of the columns corresponding to values. |
Values | [out] Matrix values. |
If GlobalRow
does not belong to this node, then indices
is set to null
.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
pure virtual |
Get a constant, nonpersisting, locally indexed view of the given row of the matrix.
The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this
. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.
This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.
isLocallyIndexed () && supportsRowViews ()
indices.size () == getNumEntriesInGlobalRow (LocalRow)
LocalRow | [in] Local index of the row. |
Indices | [out] Local indices of the columns corresponding to values. |
Values | [out] Matrix values. |
If LocalRow
does not belong to this node, then indices
is set to null
.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
virtual |
Get a constant, nonpersisting, locally indexed view of the given row of the matrix, using "raw" pointers instead of Teuchos::ArrayView.
The returned views of the column indices and values are not guaranteed to persist beyond the lifetime of this
. Furthermore, some RowMatrix implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.
This method only gets the entries in the given row that are stored on the calling process. Note that if the matrix has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.
isLocallyIndexed () && supportsRowViews ()
numEnt == getNumEntriesInGlobalRow (LocalRow)
lclRow | [in] Local index of the row. |
numEnt | [out] Number of entries in the row that are stored on the calling process. |
lclColInds | [out] Local indices of the columns corresponding to values. |
vals | [out] Matrix values. |
Reimplemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 598 of file Tpetra_RowMatrix_def.hpp.
|
pure virtual |
Get a copy of the diagonal entries, distributed by the row Map.
On input, the Vector's Map must be the same as the row Map of the matrix. (That is, this->getRowMap ()->isSameAs (* (diag.getMap ())) == true
.)
On return, the entries of diag
are filled with the diagonal entries of the matrix stored on this process. Note that if the row Map is overlapping, multiple processes may own the same diagonal element. You may combine these overlapping diagonal elements by doing an Export from the row Map Vector to a range Map Vector.
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Scale the matrix on the left with the given Vector.
On return, for all entries i,j in the matrix, .
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Scale the matrix on the right with the given Vector.
On return, for all entries i,j in the matrix, .
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The Frobenius norm of the matrix.
This method computes and returns the Frobenius norm of the matrix. The Frobenius norm for the matrix is defined as . It has the same value as the Euclidean norm of a vector made by stacking the columns of .
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >.
|
virtual |
Return a new RowMatrix which is the result of beta*this + alpha*A
.
The new RowMatrix is actually a CrsMatrix (which see). Note that RowMatrix is a read-only interface (not counting the left and right scale methods), so it is impossible to implement an in-place add using just that interface.
For brevity, call this matrix B, and the result matrix C. C's row Map will be identical to B's row Map. It is correct, though less efficient, for A and B not to have the same row Maps. We could make C's row Map the union of the two row Maps in that case. However, we don't want row Maps to grow for a repeated sequence of additions with matrices with different row Maps. Furthermore, the fact that the user called this method on B, rather than on A, suggests a preference for using B's distribution. The most reasonable thing to do, then, is to use B's row Map for C.
A and B must have identical or congruent communicators. This method must be called as a collective over B's communicator.
The parameters are optional and may be null. Here are the parameters that this function accepts:
bool
): If true, call fillComplete on the result matrix C. This is true
by default.It is not strictly necessary that a RowMatrix always have a domain and range Map. For example, a CrsMatrix does not have a domain and range Map until after its first fillComplete call. Neither A nor B need to have a domain and range Map in order to call add(). If at least one of them has a domain and range Map, you need not supply a domain and range Map to this method. If you ask this method to call fillComplete on C (it does by default), it will supply the any missing domain or range Maps from either B's or A's (in that order) domain and range Maps. If neither A nor B have a domain and range Map, and if you ask this method to call fillComplete, then you must supply both a domain Map and a range Map to this method.
This method comes with a default implementation, since the RowMatrix interface suffices for implementing it. Subclasses (like CrsMatrix) may override this implementation, for example to improve its performance, given additional knowledge about the subclass. Subclass implementations may need to do a dynamic cast on A in order to know its type.
Reimplemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 57 of file Tpetra_RowMatrix_def.hpp.
|
virtual |
Pack this object's data for an Import or Export.
Subclasses may override this method to speed up or otherwise improve the implementation by exploiting more specific details of the subclass.
Implements Tpetra::Packable< char, LocalOrdinal >.
Definition at line 298 of file Tpetra_RowMatrix_def.hpp.
|
pure virtualinherited |
The Map associated with the domain of this operator, which must be compatible with X.getMap().
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::CrsMatrixMultiplyOp< Scalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >, Tpetra::details::ApplyOp< Scalar, OperatorType >, Tpetra::RTI::KernelOp< S, LO, GO, Node, Kernel >, and Tpetra::RTI::KernelOp< S, LO, GO, Node, Tpetra::RTI::detail::BinaryFunctorAdapterWithAlphaBeta< Op, S > >.
|
pure virtualinherited |
The Map associated with the range of this operator, which must be compatible with Y.getMap().
Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::CrsMatrixMultiplyOp< Scalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >, Tpetra::details::ApplyOp< Scalar, OperatorType >, Tpetra::RTI::KernelOp< S, LO, GO, Node, Kernel >, and Tpetra::RTI::KernelOp< S, LO, GO, Node, Tpetra::RTI::detail::BinaryFunctorAdapterWithAlphaBeta< Op, S > >.
|
pure virtualinherited |
Computes the operator-multivector application.
Loosely, performs . However, the details of operation vary according to the values of alpha
and beta
. Specifically
beta == 0
, apply() must overwrite Y
, so that any values in Y
(including NaNs) are ignored.alpha == 0
, apply() may short-circuit the operator, so that any values in X
(including NaNs) are ignored. Implemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::CrsMatrixMultiplyOp< Scalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::details::ApplyOp< Scalar, OperatorType >, Tpetra::RTI::KernelOp< S, LO, GO, Node, Kernel >, and Tpetra::RTI::KernelOp< S, LO, GO, Node, Tpetra::RTI::detail::BinaryFunctorAdapterWithAlphaBeta< Op, S > >.
|
virtualinherited |
Whether this operator supports applying the transpose or conjugate transpose.
By default, this returns false. Subclasses must override this method if they can support apply() with mode=Teuchos::TRANS
or mode=Teuchos::CONJ_TRANS
.
Reimplemented in Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::CrsMatrixMultiplyOp< Scalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node >, Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >, and Tpetra::details::ApplyOp< Scalar, OperatorType >.
Definition at line 138 of file Tpetra_Operator.hpp.