| Ifpack2 Templated Preconditioning Package
    Version 1.0
    | 
Filter based on matrix entries. More...
#include <Ifpack2_DropFilter_decl.hpp>

| Public Member Functions | |
| Constructor & destructor methods | |
| DropFilter (const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix, magnitudeType DropTol) | |
| Constructor.  More... | |
| virtual | ~DropFilter () | 
| Destructor.  More... | |
| Matrix Query Methods | |
| 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 | getNodeNumRows () const | 
| Returns the number of rows owned on the calling node.  More... | |
| virtual size_t | getNodeNumCols () 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 | getNodeNumEntries () 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 size_t | getGlobalMaxNumRowEntries () const | 
| Returns the maximum number of entries across all rows/columns on all nodes.  More... | |
| virtual size_t | getNodeMaxNumRowEntries () 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 trueif fillComplete() has been called.  More... | |
| virtual bool | supportsRowViews () const | 
| Returns trueif RowViews are supported.  More... | |
| Extraction Methods | |
| virtual void | getGlobalRowCopy (GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &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 DropRow, 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 calling routine.  More... | |
| virtual void | getGlobalRowView (GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const | 
| Extract a const, non-persisting view of global indices in a specified row of the matrix.  More... | |
| virtual void | getLocalRowView (LocalOrdinal DropRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &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... | |
| Mathematical Methods | |
| 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... | |
|  Public Member Functions inherited from Ifpack2::Details::RowMatrix< MatrixType > | |
| virtual | ~RowMatrix ()=default | 
| Destructor (virtual for memory safety of derived classes)  More... | |
| Additional Inherited Members | |
|  Public Types inherited from Ifpack2::Details::RowMatrix< MatrixType > | |
| using | scalar_type = typename MatrixType::scalar_type | 
| using | local_ordinal_type = typename MatrixType::local_ordinal_type | 
| using | global_ordinal_type = typename MatrixType::global_ordinal_type | 
| using | node_type = typename MatrixType::node_type | 
Filter based on matrix entries.
| MatrixType | Tpetra::RowMatrix specialization | 
This class filters a matrix, by dropping all elements whose absolute value is below a specified threshold. The matrix should be "localized," which means that it has no remote entries (e.g., is a LocalFilter).
A typical use is as follows:
| 
 | explicit | 
Constructor.
| 
 | virtual | 
Destructor.
| 
 | virtual | 
Returns the communicator.
| 
 | virtual | 
Returns the Map that describes the row distribution in this matrix.
| 
 | virtual | 
Returns the Map that describes the column distribution in this matrix.
| 
 | virtual | 
Returns the Map that describes the domain distribution in this matrix.
| 
 | virtual | 
Returns the Map that describes the range distribution in this matrix.
| 
 | virtual | 
Returns the RowGraph associated with this matrix.
| 
 | virtual | 
Returns the number of global rows in this matrix.
| 
 | virtual | 
Returns the number of global columns in this matrix.
| 
 | virtual | 
Returns the number of rows owned on the calling node.
| 
 | virtual | 
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.
| 
 | virtual | 
Returns the index base for global indices for this matrix.
| 
 | virtual | 
Returns the global number of entries in this matrix.
| 
 | virtual | 
Returns the local number of entries in this matrix.
| 
 | virtual | 
Returns the current number of entries on this node in the specified global row.
Returns Teuchos::OrdinalTraits<size_t>::invalid() if the specified global row does not belong to this graph.
| 
 | virtual | 
Returns the current number of entries on this node in the specified local row.
Returns Teuchos::OrdinalTraits<size_t>::invalid() if the specified local row is not valid for this graph.
| 
 | virtual | 
Returns the maximum number of entries across all rows/columns on all nodes.
| 
 | virtual | 
Returns the maximum number of entries across all rows/columns on this node.
| 
 | virtual | 
Indicates whether this matrix has a well-defined column map.
| 
 | virtual | 
If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */.
| 
 | virtual | 
If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */.
| 
 | virtual | 
Returns true if fillComplete() has been called. 
| 
 | virtual | 
Returns true if RowViews are supported. 
| 
 | virtual | 
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage.
| DropRow | - (In) Global row number for which indices are desired. | 
| Indices | - (Out) Global column indices corresponding to values. | 
| Values | - (Out) Matrix values. | 
| NumEntries | - (Out) Number of indices. | 
Note: A std::runtime_error exception is thrown if either Indices or Values is not large enough to hold the data associated with row GlobalRow. If GlobalRow does not belong to this node, then Indices and Values are unchanged and NumIndices is returned as Teuchos::OrdinalTraits<size_t>::invalid(). 
| 
 | virtual | 
Extract a list of entries in a specified local row of the graph. Put into storage allocated by calling routine.
| DropRow | - (In) Drop row number for which indices are desired. | 
| Indices | - (Out) Drop column indices corresponding to values. | 
| Values | - (Out) Matrix values. | 
| NumIndices | - (Out) Number of indices. | 
Note: A std::runtime_error exception is thrown if either Indices or Values is not large enough to hold the data associated with row DropRow. If DropRow is not valid for this node, then Indices and Values are unchanged and NumIndices is returned as Teuchos::OrdinalTraits<size_t>::invalid(). 
| 
 | virtual | 
Extract a const, non-persisting view of global indices in a specified row of the matrix.
| GlobalRow | - (In) Global row number for which indices are desired. | 
| Indices | - (Out) Global column indices corresponding to values. | 
| Values | - (Out) Row values | 
isDroplyIndexed() == false indices.size() == getNumEntriesInGlobalRow(GlobalRow)Note: If GlobalRow does not belong to this node, then indices is set to null. 
| 
 | virtual | 
Extract a const, non-persisting view of local indices in a specified row of the matrix.
| DropRow | - (In) Drop row number for which indices are desired. | 
| Indices | - (Out) Global column indices corresponding to values. | 
| Values | - (Out) Row values | 
isGloballyIndexed() == false indices.size() == getNumEntriesInDropRow(DropRow)Note: If DropRow does not belong to this node, then indices is set to null. 
| 
 | virtual | 
Get a copy of the diagonal entries owned by this node, with local row indices.
Returns a distributed Vector object partitioned according to this matrix's row map, containing the the zero and non-zero diagonals owned by this node.
| 
 | virtual | 
Scales the RowMatrix on the left with the Vector x.
This matrix will be scaled such that A(i,j) = x(i)*A(i,j) where i denoes the global row number of A and j denotes the global column number of A.
| x | A vector to left scale this matrix. | 
| 
 | virtual | 
Scales the RowMatrix on the right with the Vector x.
This matrix will be scaled such that A(i,j) = x(j)*A(i,j) where i denoes the global row number of A and j denotes the global column number of A.
| x | A vector to right scale this matrix. | 
| 
 | virtual | 
Returns the Frobenius norm of the matrix.
Computes and returns the Frobenius norm of the matrix, defined as: \( \|A\|_F = \sqrt{\sum_{i,j} \|\a_{ij}\|^2} \)
| 
 | virtual | 
Computes the operator-multivector application.
Loosely, performs \(Y = \alpha \cdot A^{\textrm{mode}} \cdot X + \beta \cdot Y\). However, the details of operation vary according to the values of alpha and beta. Specifically
| 
 | virtual | 
Indicates whether this operator supports applying the adjoint operator.
 1.8.5
 1.8.5