Ifpack2 Templated Preconditioning Package
Version 1.0
|
Access only local rows and columns of a sparse matrix. More...
#include <Ifpack2_LocalFilter_decl.hpp>
Public Types | |
Typedefs | |
typedef MatrixType::scalar_type | scalar_type |
The type of the entries of the input MatrixType. More... | |
typedef MatrixType::local_ordinal_type | local_ordinal_type |
The type of local indices in the input MatrixType. More... | |
typedef MatrixType::global_ordinal_type | global_ordinal_type |
The type of global indices in the input MatrixType. More... | |
typedef MatrixType::node_type | node_type |
The Node type used by the input MatrixType. More... | |
typedef MatrixType::global_inds_host_view_type | global_inds_host_view_type |
typedef MatrixType::local_inds_host_view_type | local_inds_host_view_type |
typedef MatrixType::values_host_view_type | values_host_view_type |
typedef MatrixType::nonconst_global_inds_host_view_type | nonconst_global_inds_host_view_type |
typedef MatrixType::nonconst_local_inds_host_view_type | nonconst_local_inds_host_view_type |
typedef MatrixType::nonconst_values_host_view_type | nonconst_values_host_view_type |
typedef Teuchos::ScalarTraits < scalar_type >::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
typedef Tpetra::RowMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | row_matrix_type |
Type of the Tpetra::RowMatrix specialization that this class uses. More... | |
typedef Tpetra::RowGraph < local_ordinal_type, global_ordinal_type, node_type > | row_graph_type |
Type of the Tpetra::RowGraph specialization that this class uses. More... | |
typedef Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > | map_type |
Type of the Tpetra::Map specialization that this class uses. More... | |
typedef row_matrix_type::mag_type | mag_type |
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 |
Public Member Functions | |
Implementation of Teuchos::Describable | |
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... | |
Constructor and destructor | |
LocalFilter (const Teuchos::RCP< const row_matrix_type > &A) | |
Constructor. More... | |
virtual | ~LocalFilter () |
Destructor. More... | |
Matrix Query Methods | |
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 | getLocalNumRows () const |
The number of rows owned on the calling process. More... | |
virtual size_t | getLocalNumCols () 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 | getLocalNumEntries () const |
Returns the local number of entries in this matrix. More... | |
virtual local_ordinal_type | getBlockSize () const |
The number of degrees of freedom per mesh point. 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 | getLocalMaxNumRowEntries () 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... | |
Methods for getting the entries in a row. | |
virtual void | getGlobalRowCopy (global_ordinal_type GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const |
Get the entries in the given row, using global indices. More... | |
virtual void | getLocalRowCopy (local_ordinal_type LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const |
Get the entries in the given row, using local indices. More... | |
virtual void | getGlobalRowView (global_ordinal_type 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 (local_ordinal_type 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_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const |
Get the diagonal entries of the (locally filtered) matrix. More... | |
Mathematical methods | |
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... | |
Public Member Functions inherited from Ifpack2::Details::RowMatrix< MatrixType > | |
virtual | ~RowMatrix ()=default |
Destructor (virtual for memory safety of derived classes) More... | |
Access only local rows and columns of a sparse matrix.
MatrixType | A specialization of Tpetra::RowMatrix. |
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.
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.
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:
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.
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.
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.
Here is an example of how to apply a LocalFilter to an existing Tpetra sparse matrix:
Here is an example of how to apply a LocalFilter, using A
and A_local
from the example above.
typedef MatrixType::scalar_type Ifpack2::LocalFilter< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::LocalFilter< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::LocalFilter< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type Ifpack2::LocalFilter< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::LocalFilter< MatrixType >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::LocalFilter< MatrixType >::row_matrix_type |
Type of the Tpetra::RowMatrix specialization that this class uses.
typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::LocalFilter< MatrixType >::row_graph_type |
Type of the Tpetra::RowGraph specialization that this class uses.
typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::LocalFilter< MatrixType >::map_type |
Type of the Tpetra::Map specialization that this class uses.
|
explicit |
Constructor.
A | [in] The sparse matrix to which to apply the local filter. |
This class will not modify the input matrix.
|
virtual |
Destructor.
|
virtual |
A one-line description of this object.
Reimplemented from Teuchos::Describable.
|
virtual |
Print the object to the given output stream.
Reimplemented from Teuchos::Describable.
|
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 |
The (locally filtered) matrix's graph.
|
virtual |
The number of global rows in this matrix.
|
virtual |
The number of global columns in this matrix.
|
virtual |
The number of rows owned on the calling process.
|
virtual |
The number of columns in the (locally filtered) matrix.
|
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 |
The number of degrees of freedom per mesh point.
|
virtual |
The current number of entries on this node in the specified global row.
Teuchos::OrdinalTraits<size_t>::invalid()
if the specified row is not owned by this process, otherwise the number of entries in that row on this process.
|
virtual |
The current number of entries on this node in the specified local row.
Teuchos::OrdinalTraits<size_t>::invalid()
if the specified local row is not valid on this process, otherwise the number of entries in that row on this process.
|
virtual |
The maximum number of entries across all rows/columns on all processes.
|
virtual |
The maximum number of entries across all rows/columns on this process.
|
virtual |
Whether this matrix has a well-defined column Map.
|
virtual |
Whether the underlying sparse matrix is locally (opposite of globally) indexed.
|
virtual |
Whether the underlying sparse matrix is globally (opposite of locally) indexed.
|
virtual |
Returns true
if fillComplete() has been called.
|
virtual |
Returns true
if RowViews are supported.
|
virtual |
Get the entries in the given row, using global indices.
GlobalRow | [i] Global index of the row for which indices are desired. |
Indices | [out] Global indices of the entries in the given row. |
Values | [out] Values of the entries in the given row. |
NumEntries | [out] Number of entries in the row. |
This method throws an exception if either output array is not large enough to hold the data associated with row GlobalRow
. If GlobalRow
does not belong to the calling process, then Indices
and Values
are unchanged and NumIndices
is Teuchos::OrdinalTraits<size_t>::invalid()
on output.
|
virtual |
Get the entries in the given row, using local indices.
LocalRow | [i] Local index of the row for which indices are desired. |
Indices | [out] Local indices of the entries in the given row. |
Values | [out] Values of the entries in the given row. |
NumEntries | [out] Number of entries in the row. |
This method throws an exception if either output array is not large enough to hold the data associated with row LocalRow
. If LocalRow
does not belong to the calling process, then Indices
and Values
are unchanged and NumIndices
is Teuchos::OrdinalTraits<size_t>::invalid()
on output.
|
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] Local 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 the diagonal entries of the (locally filtered) matrix.
diag | [in/out] On input: a Tpetra::Vector whose Map is the same as the local filter's row Map, which in turn is the same as the original matrix's row Map. On output: filled with the diagonal entries owned by the calling process. |
|
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 |
The Frobenius norm of the (locally filtered) matrix.
This method may return a different value on each process, because this method computes the norm of the locally filtered matrix, which may be different on each process.
The Frobenius norm of a matrix \(A\) is defined as \(\|A\|_F = \sqrt{\sum_{i,j} \|A_{ij}\|^2}\).
|
virtual |
|
virtual |
Whether this operator supports applying the transpose or conjugate transpose.
|
virtual |
Return matrix that LocalFilter was built on.