Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
Ifpack2::LocalFilter< MatrixType > Class Template Reference

Access only local rows and columns of a sparse matrix. More...

#include <Ifpack2_LocalFilter_decl.hpp>

Inheritance diagram for Ifpack2::LocalFilter< MatrixType >:
Inheritance graph
[legend]

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 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::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 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...
 
Methods for getting the entries in a row.
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...
 
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...
 

Detailed Description

template<class MatrixType>
class Ifpack2::LocalFilter< MatrixType >

Access only local rows and columns of a sparse matrix.

Template Parameters
MatrixTypeA 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:

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::RowMatrix<double> row_matrix_type;
typedef Tpetra::CrsMatrix<double> crs_matrix_type;
RCP<crs_matrix_type> A = ...;
// ... fill the entries of 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;
// Make vectors x and y suitable for A->apply()
vec_type x (A.domainMap ());
vec_type y (A.rangeMap ());
// ... fill x ...
A->apply (x, y);
// Reinterpret x and y as local views suitable for A_local->apply()
RCP<const vec_type> x_local = x.offsetView (A_local.getDomainMap (), 0);
RCP<vec_type> y_local = y.offsetViewNonConst (A_local.getRangeMap (), 0);
// Apply the locally filtered version of A
A_local.apply (*x_local, *y_local);

Member Typedef Documentation

template<class MatrixType>
typedef MatrixType::scalar_type Ifpack2::LocalFilter< MatrixType >::scalar_type

The type of the entries of the input MatrixType.

template<class MatrixType>
typedef MatrixType::local_ordinal_type Ifpack2::LocalFilter< MatrixType >::local_ordinal_type

The type of local indices in the input MatrixType.

template<class MatrixType>
typedef MatrixType::global_ordinal_type Ifpack2::LocalFilter< MatrixType >::global_ordinal_type

The type of global indices in the input MatrixType.

template<class MatrixType>
typedef MatrixType::node_type Ifpack2::LocalFilter< MatrixType >::node_type

The Node type used by the input MatrixType.

template<class MatrixType>
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::LocalFilter< MatrixType >::magnitude_type

The type of the magnitude (absolute value) of a matrix entry.

template<class MatrixType>
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.

template<class MatrixType>
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.

Constructor & Destructor Documentation

template<class MatrixType >
Ifpack2::LocalFilter< MatrixType >::LocalFilter ( const Teuchos::RCP< const row_matrix_type > &  A)
explicit

Constructor.

Parameters
A[in] The sparse matrix to which to apply the local filter.

This class will not modify the input matrix.

template<class MatrixType >
Ifpack2::LocalFilter< MatrixType >::~LocalFilter ( )
virtual

Destructor.

Member Function Documentation

template<class MatrixType >
std::string Ifpack2::LocalFilter< MatrixType >::description ( ) const
virtual

A one-line description of this object.

Reimplemented from Teuchos::Describable.

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
virtual

Print the object to the given output stream.

Reimplemented from Teuchos::Describable.

template<class MatrixType >
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::LocalFilter< MatrixType >::getComm ( ) const
virtual

Returns the communicator.

template<class MatrixType >
Teuchos::RCP< const Tpetra::Map< typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::LocalFilter< MatrixType >::getRowMap ( ) const
virtual

Returns the Map that describes the row distribution in this matrix.

template<class MatrixType >
Teuchos::RCP< const Tpetra::Map< typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::LocalFilter< MatrixType >::getColMap ( ) const
virtual

Returns the Map that describes the column distribution in this matrix.

template<class MatrixType >
Teuchos::RCP< const Tpetra::Map< typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::LocalFilter< MatrixType >::getDomainMap ( ) const
virtual

Returns the Map that describes the domain distribution in this matrix.

template<class MatrixType >
Teuchos::RCP< const Tpetra::Map< typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::LocalFilter< MatrixType >::getRangeMap ( ) const
virtual

Returns the Map that describes the range distribution in this matrix.

template<class MatrixType >
Teuchos::RCP< const Tpetra::RowGraph< typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::LocalFilter< MatrixType >::getGraph ( ) const
virtual

The (locally filtered) matrix's graph.

template<class MatrixType >
global_size_t Ifpack2::LocalFilter< MatrixType >::getGlobalNumRows ( ) const
virtual

The number of global rows in this matrix.

template<class MatrixType >
global_size_t Ifpack2::LocalFilter< MatrixType >::getGlobalNumCols ( ) const
virtual

The number of global columns in this matrix.

template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getNodeNumRows ( ) const
virtual

The number of rows owned on the calling process.

template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getNodeNumCols ( ) const
virtual

The number of columns in the (locally filtered) matrix.

template<class MatrixType >
MatrixType::global_ordinal_type Ifpack2::LocalFilter< MatrixType >::getIndexBase ( ) const
virtual

Returns the index base for global indices for this matrix.

template<class MatrixType >
global_size_t Ifpack2::LocalFilter< MatrixType >::getGlobalNumEntries ( ) const
virtual

Returns the global number of entries in this matrix.

template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getNodeNumEntries ( ) const
virtual

Returns the local number of entries in this matrix.

template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getNumEntriesInGlobalRow ( global_ordinal_type  globalRow) const
virtual

The current number of entries on this node in the specified global row.

Returns
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.
template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getNumEntriesInLocalRow ( local_ordinal_type  localRow) const
virtual

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 on this process, otherwise the number of entries in that row on this process.
template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getGlobalMaxNumRowEntries ( ) const
virtual

The maximum number of entries across all rows/columns on all processes.

template<class MatrixType >
size_t Ifpack2::LocalFilter< MatrixType >::getNodeMaxNumRowEntries ( ) const
virtual

The maximum number of entries across all rows/columns on this process.

template<class MatrixType >
bool Ifpack2::LocalFilter< MatrixType >::hasColMap ( ) const
virtual

Whether this matrix has a well-defined column Map.

template<class MatrixType >
bool Ifpack2::LocalFilter< MatrixType >::isLocallyIndexed ( ) const
virtual

Whether the underlying sparse matrix is locally (opposite of globally) indexed.

template<class MatrixType >
bool Ifpack2::LocalFilter< MatrixType >::isGloballyIndexed ( ) const
virtual

Whether the underlying sparse matrix is globally (opposite of locally) indexed.

template<class MatrixType >
bool Ifpack2::LocalFilter< MatrixType >::isFillComplete ( ) const
virtual

Returns true if fillComplete() has been called.

template<class MatrixType >
bool Ifpack2::LocalFilter< MatrixType >::supportsRowViews ( ) const
virtual

Returns true if RowViews are supported.

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::getGlobalRowCopy ( global_ordinal_type  GlobalRow,
const Teuchos::ArrayView< global_ordinal_type > &  Indices,
const Teuchos::ArrayView< scalar_type > &  Values,
size_t &  NumEntries 
) const
virtual

Get the entries in the given row, using global indices.

Parameters
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.

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::getLocalRowCopy ( local_ordinal_type  LocalRow,
const Teuchos::ArrayView< local_ordinal_type > &  Indices,
const Teuchos::ArrayView< scalar_type > &  Values,
size_t &  NumEntries 
) const
virtual

Get the entries in the given row, using local indices.

Parameters
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.

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::getGlobalRowView ( global_ordinal_type  GlobalRow,
Teuchos::ArrayView< const global_ordinal_type > &  indices,
Teuchos::ArrayView< const scalar_type > &  values 
) const
virtual

Extract a const, non-persisting view of global indices in a specified row of the matrix.

Parameters
GlobalRow[in] Global row number for which indices are desired.
Indices[out] Global column indices corresponding to values.
Values[out] Row values
Precondition
isLocallyIndexed() == false
Postcondition
indices.size() == getNumEntriesInGlobalRow(GlobalRow)

Note: If GlobalRow does not belong to this node, then indices is set to null.

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::getLocalRowView ( local_ordinal_type  LocalRow,
Teuchos::ArrayView< const local_ordinal_type > &  indices,
Teuchos::ArrayView< const scalar_type > &  values 
) const
virtual

Extract a const, non-persisting view of local indices in a specified row of the matrix.

Parameters
LocalRow[in] Local row number for which indices are desired.
Indices[out] Local column indices corresponding to values.
Values[out] Row values
Precondition
isGloballyIndexed() == false
Postcondition
indices.size() == getNumEntriesInLocalRow(LocalRow)

Note: If LocalRow does not belong to this node, then indices is set to null.

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::getLocalDiagCopy ( Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  diag) const
virtual

Get the diagonal entries of the (locally filtered) matrix.

Parameters
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.
template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::leftScale ( const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  x)
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.

Parameters
xA vector to left scale this matrix.
template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::rightScale ( const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  x)
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.

Parameters
xA vector to right scale this matrix.
template<class MatrixType >
LocalFilter< MatrixType >::mag_type Ifpack2::LocalFilter< MatrixType >::getFrobeniusNorm ( ) const
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}\).

template<class MatrixType >
void Ifpack2::LocalFilter< MatrixType >::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
virtual

Compute Y = beta*Y + alpha*A_local*X.

If beta == 0, apply() must overwrite Y, so that any values in Y (including NaNs) are ignored. If alpha == 0, apply() must short-circuit the operator, so that any values in X (including NaNs) are ignored.

template<class MatrixType >
bool Ifpack2::LocalFilter< MatrixType >::hasTransposeApply ( ) const
virtual

Whether this operator supports applying the transpose or conjugate transpose.

template<class MatrixType >
Teuchos::RCP< const Tpetra::RowMatrix< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::LocalFilter< MatrixType >::getUnderlyingMatrix ( ) const
virtual

Return matrix that LocalFilter was built on.


The documentation for this class was generated from the following files: