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

Drop entries of a matrix, based on the sparsity pattern. More...

#include <Ifpack2_SparsityFilter_decl.hpp>

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

Public Member Functions

Constructor & destructor methods
 SparsityFilter (const Teuchos::RCP< const row_matrix_type > &Matrix, size_t AllowedNumEntries, LocalOrdinal AllowedBandwidth=-Teuchos::ScalarTraits< LocalOrdinal >::one())
 Constructor. More...
 
virtual ~SparsityFilter ()
 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 getLocalNumRows () const
 Returns the number of rows owned on the calling node. More...
 
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. 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 getLocalNumEntries () 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 getLocalMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on this node. More...
 
virtual LocalOrdinal getBlockSize () const
 The number of degrees of freedom per mesh point. 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 true if fillComplete() has been called. More...
 
virtual bool supportsRowViews () const
 Returns true if RowViews are supported. More...
 
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. More...
 
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. More...
 
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. More...
 
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. 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
 

Detailed Description

template<class MatrixType>
class Ifpack2::SparsityFilter< MatrixType >

Drop entries of a matrix, based on the sparsity pattern.

Template Parameters
MatrixTypeA specialization of Tpetra::RowMatrix.

This class takes an existing Tpetra::RowMatrix, and wraps it in a matrix interface that drops entries using the following criteria:

  1. Drop elements beyond a certain distance from the diagonal, and / or
  2. Keep only the N largest entries in each row.

Here is a typical use case:

// First filter out all entries in columns
// which are not in A's domain Map.
// Drop all but the largest 20 elements in each row.
const size_t maxEntriesPerRow = 20;
// Exclude elements in each row whose local column index
// is more than 10 greater or less than the local column
// index of the diagonal element in that row.
const int maxBw = 10;
// Now create the sparsity filter. Elements dropped are
// not included in calls to getLocalRowCopy() and apply().
Ifpack2::SparsityFilter<Tpetra::RowMatrix<> > A_drop (A_local, maxEntriesPerRow, maxBw);

This class currently only works if the matrix is the result of a LocalFilter, that is, if the row and column Maps are the same.

Constructor & Destructor Documentation

template<class MatrixType>
Ifpack2::SparsityFilter< MatrixType >::SparsityFilter ( const Teuchos::RCP< const row_matrix_type > &  Matrix,
size_t  AllowedNumEntries,
LocalOrdinal  AllowedBandwidth = -Teuchos::ScalarTraits<LocalOrdinal>::one() 
)
explicit

Constructor.

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

Destructor.

Member Function Documentation

template<class MatrixType >
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::SparsityFilter< 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::SparsityFilter< 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::SparsityFilter< 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::SparsityFilter< 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::SparsityFilter< 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::SparsityFilter< MatrixType >::getGraph ( ) const
virtual

Returns the RowGraph associated with this matrix.

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

Returns the number of global rows in this matrix.

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

Returns the number of global columns in this matrix.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getLocalNumRows ( ) const
virtual

Returns the number of rows owned on the calling node.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getLocalNumCols ( ) const
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.

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

Returns the index base for global indices for this matrix.

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

Returns the global number of entries in this matrix.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getLocalNumEntries ( ) const
virtual

Returns the local number of entries in this matrix.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getNumEntriesInGlobalRow ( GlobalOrdinal  globalRow) const
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.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getNumEntriesInLocalRow ( LocalOrdinal  localRow) const
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.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getGlobalMaxNumRowEntries ( ) const
virtual

Returns the maximum number of entries across all rows/columns on all nodes.

template<class MatrixType >
size_t Ifpack2::SparsityFilter< MatrixType >::getLocalMaxNumRowEntries ( ) const
virtual

Returns the maximum number of entries across all rows/columns on this node.

template<class MatrixType >
MatrixType::local_ordinal_type Ifpack2::SparsityFilter< MatrixType >::getBlockSize ( ) const
virtual

The number of degrees of freedom per mesh point.

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

Indicates whether this matrix has a well-defined column map.

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

If matrix indices are in the local range, this function returns true. Otherwise, this function returns false. */.

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

If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. */.

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

Returns true if fillComplete() has been called.

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

Returns true if RowViews are supported.

template<class MatrixType >
void Ifpack2::SparsityFilter< MatrixType >::getGlobalRowCopy ( GlobalOrdinal  GlobalRow,
nonconst_global_inds_host_view_type &  Indices,
nonconst_values_host_view_type &  Values,
size_t &  NumEntries 
) const
virtual

Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage.

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

template<class MatrixType >
void Ifpack2::SparsityFilter< MatrixType >::getLocalRowCopy ( LocalOrdinal  LocalRow,
nonconst_local_inds_host_view_type &  Indices,
nonconst_values_host_view_type &  Values,
size_t &  NumEntries 
) const
virtual

Extract a list of entries in a specified local row of the graph. Put into storage allocated by calling routine.

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

template<class MatrixType >
void Ifpack2::SparsityFilter< MatrixType >::getGlobalRowView ( GlobalOrdinal  GlobalRow,
global_inds_host_view_type &  indices,
values_host_view_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
isDroplyIndexed() == 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::SparsityFilter< MatrixType >::getLocalRowView ( LocalOrdinal  LocalRow,
local_inds_host_view_type &  indices,
values_host_view_type &  values 
) const
virtual

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

Parameters
DropRow- (In) Drop row number for which indices are desired.
Indices- (Out) Global column indices corresponding to values.
Values- (Out) Row values
Precondition
isGloballyIndexed() == false
Postcondition
indices.size() == getNumEntriesInDropRow(DropRow)

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

template<class MatrixType >
void Ifpack2::SparsityFilter< MatrixType >::getLocalDiagCopy ( Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  diag) const
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.

template<class MatrixType >
void Ifpack2::SparsityFilter< MatrixType >::leftScale ( const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  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::SparsityFilter< MatrixType >::rightScale ( const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  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 >
SparsityFilter< MatrixType >::mag_type Ifpack2::SparsityFilter< MatrixType >::getFrobeniusNorm ( ) const
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} \)

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

  • if beta == 0, apply() must overwrite Y, so that any values in Y (including NaNs) are ignored.
  • if alpha == 0, apply() may short-circuit the operator, so that any values in X (including NaNs) are ignored.
template<class MatrixType >
bool Ifpack2::SparsityFilter< MatrixType >::hasTransposeApply ( ) const
virtual

Indicates whether this operator supports applying the adjoint operator.


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