Tpetra parallel linear algebra
Version of the Day
|
Wrap a CrsMatrix instance's triangular solve in an Operator. More...
#include <Tpetra_CrsMatrixSolveOp.hpp>
Public Types | |
typedef CrsMatrix< MatScalar, LocalOrdinal, GlobalOrdinal, Node > | crs_matrix_type |
The specialization of CrsMatrix which this class wraps. More... | |
typedef Map< LocalOrdinal, GlobalOrdinal, Node > | map_type |
The specialization of Map which this class uses. More... | |
Typedefs that give access to the template parameters. | |
typedef Scalar | scalar_type |
The type of the entries of the input and output multivectors. More... | |
typedef LocalOrdinal | local_ordinal_type |
The local index type. More... | |
typedef GlobalOrdinal | global_ordinal_type |
The global index type. More... | |
typedef Node | node_type |
The Kokkos Node type. More... | |
Public Member Functions | |
Constructor and destructor | |
CrsMatrixSolveOp (const Teuchos::RCP< const crs_matrix_type > &A) | |
Constructor; takes a CrsMatrix to use for local triangular solves. More... | |
virtual | ~CrsMatrixSolveOp () |
Destructor (virtual for memory safety of derived classes). More... | |
Implementation of Operator | |
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 |
Compute ![]() ![]() | |
bool | hasTransposeApply () const |
Whether apply() can solve with the (conjugate) transpose of the matrix. More... | |
Teuchos::RCP< const map_type > | getDomainMap () const |
The domain Map of this operator. This is the range map of the underlying CrsMatrix. More... | |
Teuchos::RCP< const map_type > | getRangeMap () const |
The range Map of this operator. This is the domain Map of the underlying CrsMatrix. More... | |
Protected Member Functions | |
void | applyNonTranspose (const MV &X_in, MV &Y_in) const |
Do the non-transpose solve. More... | |
void | applyTranspose (const MV &X_in, MV &Y_in, const Teuchos::ETransp mode) const |
Do the transpose or conjugate transpose solve. More... | |
Protected Attributes | |
const Teuchos::RCP< const crs_matrix_type > | matrix_ |
The underlying CrsMatrix. More... | |
Teuchos::RCP< MV > | importMV_ |
Cached temporary destination of Import operation in apply(). More... | |
Teuchos::RCP< MV > | exportMV_ |
Cached temporary source of Export operation in apply(). More... | |
Related Functions | |
(Note that these are not member functions.) | |
template<class OpScalar , class MatScalar , class LocalOrdinal , class GlobalOrdinal , class Node > | |
Teuchos::RCP< CrsMatrixSolveOp < OpScalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node > > | createCrsMatrixSolveOp (const Teuchos::RCP< const CrsMatrix< MatScalar, LocalOrdinal, GlobalOrdinal, Node > > &A) |
Nonmember function that wraps a CrsMatrix in a CrsMatrixSolveOp. More... | |
Wrap a CrsMatrix instance's triangular solve in an Operator.
Scalar | Same as the first template parameter of Operator. The type of the entries of the MultiVector input and output of apply(). Not necessarily the same as the first template parameter of the CrsMatrix used to create this object. |
MatScalar | Same as the first template parameter of CrsMatrix. The type of the entries of the sparse matrix. Not necessarily the same as the type of the entries of the MultiVector input and output of apply(). |
LocalOrdinal | Same as the second template parameter of CrsMatrix and Operator. |
GlobalOrdinal | Same as the third template parameter of CrsMatrix and Operator. |
Node | Same as the fourth template parameter of CrsMatrix and Operator. |
This class' apply() method does a "local" triangular solve. "Local" is in quotes because apply() does the same communication (Import and Export) operations that CrsMatrix's apply() method would do for a sparse matrix-vector multiply, but the triangular solve is restricted to each process' part of the data. Thus, it is not a triangular solve of a fully distributed triangular matrix.
Here are some situations where this operation is useful:
MatScalar
of entries in the matrix differs from the type Scalar
of entries in the MultiVector input and output of apply(). Definition at line 93 of file Tpetra_CrsMatrixSolveOp.hpp.
typedef CrsMatrix<MatScalar, LocalOrdinal, GlobalOrdinal, Node> Tpetra::CrsMatrixSolveOp< Scalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node >::crs_matrix_type |
The specialization of CrsMatrix which this class wraps.
Definition at line 97 of file Tpetra_CrsMatrixSolveOp.hpp.
typedef Map<LocalOrdinal, GlobalOrdinal, Node> Tpetra::CrsMatrixSolveOp< Scalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node >::map_type |
The specialization of Map which this class uses.
Definition at line 99 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inherited |
The type of the entries of the input and output multivectors.
Definition at line 89 of file Tpetra_Operator.hpp.
|
inherited |
The local index type.
Definition at line 92 of file Tpetra_Operator.hpp.
|
inherited |
The global index type.
Definition at line 95 of file Tpetra_Operator.hpp.
|
inherited |
The Kokkos Node type.
Definition at line 98 of file Tpetra_Operator.hpp.
|
inline |
Constructor; takes a CrsMatrix to use for local triangular solves.
Definition at line 105 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlinevirtual |
Destructor (virtual for memory safety of derived classes).
Definition at line 110 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlinevirtual |
Compute , where
represents the result of the local triangular solve.
Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 119 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlinevirtual |
Whether apply() can solve with the (conjugate) transpose of the matrix.
Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 159 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlinevirtual |
The domain Map of this operator. This is the range map of the underlying CrsMatrix.
Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 165 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlinevirtual |
The range Map of this operator. This is the domain Map of the underlying CrsMatrix.
Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 172 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlineprotected |
Do the non-transpose solve.
Definition at line 189 of file Tpetra_CrsMatrixSolveOp.hpp.
|
inlineprotected |
Do the transpose or conjugate transpose solve.
Definition at line 276 of file Tpetra_CrsMatrixSolveOp.hpp.
|
related |
Nonmember function that wraps a CrsMatrix in a CrsMatrixSolveOp.
The function has the same template parameters of CrsMatrixSolveOp.
A | [in] The CrsMatrix instance to wrap in an CrsMatrixSolveOp. |
Definition at line 379 of file Tpetra_CrsMatrixSolveOp.hpp.
|
protected |
The underlying CrsMatrix.
Definition at line 181 of file Tpetra_CrsMatrixSolveOp.hpp.
|
mutableprotected |
Cached temporary destination of Import operation in apply().
Definition at line 184 of file Tpetra_CrsMatrixSolveOp.hpp.
|
mutableprotected |
Cached temporary source of Export operation in apply().
Definition at line 186 of file Tpetra_CrsMatrixSolveOp.hpp.