|
|
virtual std::string | description () const |
| A one-line description of this object. More...
|
|
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
| Print the object to the given output stream. More...
|
|
|
| LocalFilter (const Teuchos::RCP< const row_matrix_type > &A) |
| Constructor. More...
|
|
virtual | ~LocalFilter () |
| Destructor. More...
|
|
|
virtual Teuchos::RCP< const
Teuchos::Comm< int > > | getComm () const |
| Returns the communicator. More...
|
|
virtual Teuchos::RCP< const
map_type > | getRowMap () const |
| Returns the Map that describes the row distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
map_type > | getColMap () const |
| Returns the Map that describes the column distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
map_type > | getDomainMap () const |
| Returns the Map that describes the domain distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
map_type > | getRangeMap () const |
| Returns the Map that describes the range distribution in this matrix. More...
|
|
virtual Teuchos::RCP< const
Tpetra::RowGraph
< local_ordinal_type,
global_ordinal_type, node_type > > | getGraph () const |
| The (locally filtered) matrix's graph. More...
|
|
virtual global_size_t | getGlobalNumRows () const |
| The number of global rows in this matrix. More...
|
|
virtual global_size_t | getGlobalNumCols () const |
| The number of global columns in this matrix. More...
|
|
virtual size_t | getNodeNumRows () const |
| The number of rows owned on the calling process. More...
|
|
virtual size_t | getNodeNumCols () const |
| The number of columns in the (locally filtered) matrix. More...
|
|
virtual global_ordinal_type | 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 (global_ordinal_type globalRow) const |
| The current number of entries on this node in the specified global row. More...
|
|
virtual size_t | getNumEntriesInLocalRow (local_ordinal_type localRow) const |
| The current number of entries on this node in the specified local row. More...
|
|
virtual size_t | getGlobalMaxNumRowEntries () const |
| The maximum number of entries across all rows/columns on all processes. More...
|
|
virtual size_t | getNodeMaxNumRowEntries () const |
| The maximum number of entries across all rows/columns on this process. More...
|
|
virtual bool | hasColMap () const |
| Whether this matrix has a well-defined column Map. More...
|
|
virtual bool | isLocallyIndexed () const |
| Whether the underlying sparse matrix is locally (opposite of globally) indexed. More...
|
|
virtual bool | isGloballyIndexed () const |
| Whether the underlying sparse matrix is globally (opposite of locally) indexed. 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 (global_ordinal_type GlobalRow, const Teuchos::ArrayView< global_ordinal_type > &Indices, const Teuchos::ArrayView< scalar_type > &Values, size_t &NumEntries) const |
| Get the entries in the given row, using global indices. More...
|
|
virtual void | getLocalRowCopy (local_ordinal_type LocalRow, const Teuchos::ArrayView< local_ordinal_type > &Indices, const Teuchos::ArrayView< scalar_type > &Values, size_t &NumEntries) const |
| Get the entries in the given row, using local indices. More...
|
|
virtual void | getGlobalRowView (global_ordinal_type GlobalRow, Teuchos::ArrayView< const global_ordinal_type > &indices, Teuchos::ArrayView< const scalar_type > &values) const |
| Extract a const, non-persisting view of global indices in a specified row of the matrix. More...
|
|
virtual void | getLocalRowView (local_ordinal_type LocalRow, Teuchos::ArrayView< const local_ordinal_type > &indices, Teuchos::ArrayView< const scalar_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_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const |
| Get the diagonal entries of the (locally filtered) matrix. More...
|
|
|
virtual void | leftScale (const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &x) |
| Scales the RowMatrix on the left with the Vector x. More...
|
|
virtual void | rightScale (const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &x) |
| Scales the RowMatrix on the right with the Vector x. More...
|
|
virtual mag_type | getFrobeniusNorm () const |
| The Frobenius norm of the (locally filtered) matrix. More...
|
|
virtual void | apply (const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const |
| Compute Y = beta*Y + alpha*A_local*X. More...
|
|
virtual bool | hasTransposeApply () const |
| Whether this operator supports applying the transpose or conjugate transpose. More...
|
|
virtual Teuchos::RCP< const
row_matrix_type > | getUnderlyingMatrix () const |
| Return matrix that LocalFilter was built on. More...
|
|
virtual | ~RowMatrix ()=default |
| Destructor (virtual for memory safety of derived classes) More...
|
|
template<class MatrixType>
class Ifpack2::LocalFilter< MatrixType >
Access only local rows and columns of a sparse matrix.
- Template Parameters
-
MatrixType | A specialization of Tpetra::RowMatrix. |
Summary
LocalFilter makes a global matrix "appear" local. Solving a linear system with the LocalFilter of a matrix A, if possible, is equivalent to applying one step of nonoverlapping additive Schwarz to A. This class is an implementation detail of Ifpack2. It is not meant for users. It may be useful to Ifpack2 developers who want to implement a "local" preconditioner (i.e., that only works within a single MPI process), and would like the preconditioner to work for a global matrix. For example, ILUT uses LocalFilter to extract and compute the incomplete factorizaton of the "local
matrix." This makes the input to ILUT appear like a square matrix, assuming that the domain and range Maps have the same number of entries on each process. This in turn makes the result of ILUT suitable for solving linear systems.
How does it work?
View of the local rows and columns
LocalFilter provides a view of only the local rows and columns of an existing sparse matrix. "Local rows" are those rows in the row Map that are owned by the range Map on the calling process. "Local columns" are those columns in the column Map that are owned by the domain Map on the calling process. The view's communicator contains only the local process (in MPI terms, MPI_COMM_SELF
), so each process will have its own distinct view of its local part of the matrix.
Square diagonal block of the original matrix
If the following conditions hold on the Maps of a matrix A, then the view resulting from a LocalFilter of A will be that of a square diagonal block of A:
- Domain and range Maps are the same
- On every process, the row Map's entries are the same as the range Map's entries, possibly followed by additional remote entries
- On every process, the column Map's entries are the same as the domain Map's entries, possibly followed by additional remote entries
These conditions commonly hold for a Tpetra::CrsMatrix constructed without a column Map, after fillComplete() has been called on it, with default domain and range Maps.
Remapping of global indices
The global indices of the view's domain and range Maps will be different than those in the original matrix's domain and range Maps. The global indices of the new (domain, range) Map will be remapped to a consecutive sequence, corresponding exactly to the local indices of the original (domain, range) Map. This ensures that the local indices of the old and new Maps match.
Not necessarily a copy
This class does not necessarily copy the entire sparse matrix. It may choose instead just to "filter out" the nonlocal entries. The effect on the LocalFilter of changing the entries or structure of the original matrix is undefined. LocalFilter may even make different implementation choices on different MPI processes. It may do so because it does not invoke MPI communication outside of the calling process.
Examples
Here is an example of how to apply a LocalFilter to an existing Tpetra sparse matrix:
#include "Ifpack2_LocalFilter.hpp"
typedef Tpetra::CrsMatrix<double> crs_matrix_type;
RCP<crs_matrix_type> A = ...;
A->FillComplete ();
Here is an example of how to apply a LocalFilter, using A
and A_local
from the example above.
typedef Tpetra::Vector<double> vec_type;
vec_type x (A.domainMap ());
vec_type y (A.rangeMap ());
A->apply (x, y);
RCP<const vec_type> x_local = x.offsetView (A_local.getDomainMap (), 0);
RCP<vec_type> y_local = y.offsetViewNonConst (A_local.getRangeMap (), 0);
A_local.apply (*x_local, *y_local);