|
|
| SingletonFilter (const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix) |
| Constructor. More...
|
|
virtual | ~SingletonFilter () |
| Destructor. More...
|
|
|
virtual Teuchos::RCP< const
Teuchos::Comm< int > > | getComm () const |
| Returns the communicator. More...
|
|
virtual Teuchos::RCP< const
Tpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > | getRowMap () const |
| Returns the Map that describes the row distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
Tpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > | getColMap () const |
| Returns the Map that describes the column distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
Tpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > | getDomainMap () const |
| Returns the Map that describes the domain distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
Tpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > | getRangeMap () const |
| Returns the Map that describes the range distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
Tpetra::RowGraph< LocalOrdinal,
GlobalOrdinal, Node > > | getGraph () const |
| Returns the RowGraph associated with this matrix. More...
|
|
virtual global_size_t | getGlobalNumRows () const |
| Returns the number of global rows in this matrix. More...
|
|
virtual global_size_t | getGlobalNumCols () const |
| Returns the number of global columns in this matrix. More...
|
|
virtual size_t | getLocalNumRows () const |
| Returns the number of rows owned on the calling node. More...
|
|
virtual 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. More...
|
|
virtual GlobalOrdinal | getIndexBase () const |
| Returns the index base for global indices for this matrix. More...
|
|
virtual global_size_t | getGlobalNumEntries () const |
| Returns the global number of entries in this matrix. More...
|
|
virtual size_t | getLocalNumEntries () const |
| Returns the local number of entries in this matrix. More...
|
|
virtual size_t | getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const |
| Returns the current number of entries on this node in the specified global row. More...
|
|
virtual size_t | getNumEntriesInLocalRow (LocalOrdinal localRow) const |
| Returns the current number of entries on this node in the specified local row. More...
|
|
virtual LocalOrdinal | getBlockSize () const |
| The number of degrees of freedom per mesh point. More...
|
|
virtual size_t | getGlobalMaxNumRowEntries () const |
| Returns the maximum number of entries across all rows/columns on all nodes. More...
|
|
virtual size_t | getLocalMaxNumRowEntries () const |
| Returns the maximum number of entries across all rows/columns on this node. More...
|
|
virtual bool | hasColMap () const |
| Indicates whether this matrix has a well-defined column map. More...
|
|
virtual bool | isLocallyIndexed () const |
| If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */. More...
|
|
virtual bool | isGloballyIndexed () const |
| If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */. More...
|
|
virtual bool | isFillComplete () const |
| Returns true if fillComplete() has been called. More...
|
|
virtual bool | supportsRowViews () const |
| Returns true if RowViews are supported. More...
|
|
|
virtual void | getGlobalRowCopy (GlobalOrdinal GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const |
| Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage. More...
|
|
virtual void | getLocalRowCopy (LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const |
| Extract a list of entries in a specified local row of the graph. Put into storage allocated by calling routine. More...
|
|
virtual void | getGlobalRowView (GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const |
| Extract a const, non-persisting view of global indices in a specified row of the matrix. More...
|
|
virtual void | getLocalRowView (LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const |
| Extract a const, non-persisting view of local indices in a specified row of the matrix. More...
|
|
virtual void | getLocalDiagCopy (Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const |
| Get a copy of the diagonal entries owned by this node, with local row indices. More...
|
|
|
virtual void | leftScale (const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) |
| Scales the RowMatrix on the left with the Vector x. More...
|
|
virtual void | rightScale (const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) |
| Scales the RowMatrix on the right with the Vector x. More...
|
|
virtual mag_type | getFrobeniusNorm () const |
| Returns the Frobenius norm of the matrix. More...
|
|
virtual void | apply (const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::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. More...
|
|
virtual bool | hasTransposeApply () const |
| Indicates whether this operator supports applying the adjoint operator. More...
|
|
virtual void | SolveSingletons (const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS) |
| Solve the singleton components of the linear system. More...
|
|
template<class DomainScalar , class RangeScalar > |
void | SolveSingletonsTempl (const Tpetra::MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &LHS) |
|
virtual void | CreateReducedRHS (const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS, const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedRHS) |
| Creates a RHS for the reduced singleton-free system. More...
|
|
template<class DomainScalar , class RangeScalar > |
void | CreateReducedRHSTempl (const Tpetra::MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &LHS, const Tpetra::MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedRHS) |
|
virtual void | UpdateLHS (const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedLHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS) |
| Updates a full LHS from a reduces LHS. More...
|
|
template<class DomainScalar , class RangeScalar > |
void | UpdateLHSTempl (const Tpetra::MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedLHS, Tpetra::MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &LHS) |
|
virtual | ~RowMatrix ()=default |
| Destructor (virtual for memory safety of derived classes) More...
|
|
template<class MatrixType>
class Ifpack2::SingletonFilter< MatrixType >
Filter based on matrix entries.
- Template Parameters
-
MatrixType | A specialization of Tpetra::RowMatrix. |
- Warning
- This is an implementation detail of Ifpack2.