|
Ifpack2 Templated Preconditioning Package Version 1.0
|
Ifpack2_DiagonalFilter: Filter to modify the diagonal entries of a given Tpetra_RowMatrix. More...
#include <Ifpack2_DiagonalFilter_decl.hpp>

Public Member Functions | |
Constructor & destructor methods | |
| DiagonalFilter (const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix, magnitudeType AbsoluteThreshold, magnitudeType RelativeThreshold) | |
| Constructor. | |
| virtual | ~DiagonalFilter () |
| Destructor. | |
Matrix Query Methods | |
| virtual Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const |
| Returns the communicator. | |
| virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > | getRowMap () const |
| Returns the Map that describes the row distribution in this matrix. | |
| virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > | getColMap () const |
| Returns the Map that describes the column distribution in this matrix. | |
| virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > | getDomainMap () const |
| Returns the Map that describes the domain distribution in this matrix. | |
| virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > | getRangeMap () const |
| Returns the Map that describes the range distribution in this matrix. | |
| virtual Teuchos::RCP< const Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node > > | getGraph () const |
| Returns the RowGraph associated with this matrix. | |
| virtual global_size_t | getGlobalNumRows () const |
| Returns the number of global rows in this matrix. | |
| virtual global_size_t | getGlobalNumCols () const |
| Returns the number of global columns in this matrix. | |
| virtual size_t | getLocalNumRows () const |
| Returns the number of rows owned on the calling node. | |
| 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. | |
| virtual GlobalOrdinal | getIndexBase () const |
| Returns the index base for global indices for this matrix. | |
| virtual global_size_t | getGlobalNumEntries () const |
| Returns the global number of entries in this matrix. | |
| virtual size_t | getLocalNumEntries () const |
| Returns the local number of entries in this matrix. | |
| virtual size_t | getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const |
| Returns the current number of entries on this node in the specified global row. | |
| virtual size_t | getNumEntriesInLocalRow (LocalOrdinal localRow) const |
| Returns the current number of entries on this node in the specified local row. | |
| virtual size_t | getGlobalMaxNumRowEntries () const |
| Returns the maximum number of entries across all rows/columns on all nodes. | |
| virtual size_t | getLocalMaxNumRowEntries () const |
| Returns the maximum number of entries across all rows/columns on this node. | |
| virtual LocalOrdinal | getBlockSize () const |
| The number of degrees of freedom per mesh point. | |
| virtual bool | hasColMap () const |
| Indicates whether this matrix has a well-defined column map. | |
| virtual bool | isLocallyIndexed () const |
| If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */. | |
| virtual bool | isGloballyIndexed () const |
| If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */. | |
| virtual bool | isFillComplete () const |
Returns true if fillComplete() has been called. | |
| virtual bool | supportsRowViews () const |
Returns true if RowViews are supported. | |
Extraction Methods | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
Public Member Functions inherited from Ifpack2::Details::RowMatrix< MatrixType > | |
| virtual | ~RowMatrix ()=default |
| Destructor (virtual for memory safety of derived classes) | |
Mathematical Methods | |
| virtual void | leftScale (const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) |
| Scales the RowMatrix on the left with the Vector x. | |
| virtual void | rightScale (const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x) |
| Scales the RowMatrix on the right with the Vector x. | |
| virtual mag_type | getFrobeniusNorm () const |
| Returns the Frobenius norm of the matrix. | |
| 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. | |
| virtual bool | hasTransposeApply () const |
| Indicates whether this operator supports applying the adjoint operator. | |
Additional Inherited Members |
Ifpack2_DiagonalFilter: Filter to modify the diagonal entries of a given Tpetra_RowMatrix.
Ifpack2_DiagonalFilter modifies the elements on the diagonal.
A typical use is as follows:
Note: This operation only really makes sense if the Thresholds are not complex.
\data Last modified on 31-Aug-12.
|
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 |
The number of degrees of freedom per mesh point.
|
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.
| LocalRow | - (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.
| LocalRow | - (In) Local row number for which indices are desired. |
| Indices | - (Out) Local 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 LocalRow. If LocalRow 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 |
isLocallyIndexed() == 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.
| LocalRow | - (In) Local row number for which indices are desired. |
| Indices | - (Out) Global column indices corresponding to values. |
| Values | - (Out) Row values |
isGloballyIndexed() == false indices.size() == getNumEntriesInLocalRow(LocalRow)Note: If LocalRow 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
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.This is analagous to the Multiply function in Ifpack, not the Apply
|
virtual |
Indicates whether this operator supports applying the adjoint operator.