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

Wraps a Tpetra::RowMatrix in a filter that reorders local rows and columns. More...

#include <Ifpack2_ReorderFilter_decl.hpp>

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

Public Member Functions

Constructor & destructor methods
 ReorderFilter (const Teuchos::RCP< const row_matrix_type > &A, const Teuchos::ArrayRCP< local_ordinal_type > &perm, const Teuchos::ArrayRCP< local_ordinal_type > &reverseperm)
 Constructor. More...
 
virtual ~ReorderFilter ()
 Destructor. More...
 
Matrix query methods
virtual Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 The matrix's 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
 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 getNodeNumRows () const
 Returns the number of rows owned on the calling node. More...
 
virtual size_t getNodeNumCols () 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 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 in this matrix, stored on the calling process, in the row whose global index is globalRow. More...
 
virtual size_t getNumEntriesInLocalRow (local_ordinal_type localRow) const
 The current number of entries in this matrix, stored on the calling process, in the row whose local index is globalRow. More...
 
virtual size_t getGlobalMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on all nodes. More...
 
virtual size_t getNodeMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on this node. 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 (global_ordinal_type GlobalRow, const Teuchos::ArrayView< global_ordinal_type > &Indices, const Teuchos::ArrayView< scalar_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 (local_ordinal_type DropRow, const Teuchos::ArrayView< local_ordinal_type > &Indices, const Teuchos::ArrayView< scalar_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 (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 a copy of the diagonal entries owned by this node, with local row indices. 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
 Returns the Frobenius norm of the 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
 \( Y := \beta Y + \alpha Op(A) X \), where Op(A) is either A, \(A^T\), or \(A^H\). More...
 
virtual bool hasTransposeApply () const
 Whether apply() can apply the transpose or conjugate transpose. More...
 
virtual void permuteOriginalToReordered (const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &originalX, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &reorderedY) const
 Permute multivector: original-to-reordered. More...
 
template<class DomainScalar , class RangeScalar >
void permuteOriginalToReorderedTempl (const Tpetra::MultiVector< DomainScalar, local_ordinal_type, global_ordinal_type, node_type > &originalX, Tpetra::MultiVector< RangeScalar, local_ordinal_type, global_ordinal_type, node_type > &reorderedY) const
 
virtual void permuteReorderedToOriginal (const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &reorderedX, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &originalY) const
 Permute multivector: reordered-to-original. More...
 
template<class DomainScalar , class RangeScalar >
void permuteReorderedToOriginalTempl (const Tpetra::MultiVector< DomainScalar, local_ordinal_type, global_ordinal_type, node_type > &reorderedX, Tpetra::MultiVector< RangeScalar, local_ordinal_type, global_ordinal_type, node_type > &originalY) const
 
- 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::ReorderFilter< MatrixType >

Wraps a Tpetra::RowMatrix in a filter that reorders local rows and columns.

Template Parameters
MatrixTypeA specialization of Tpetra::RowMatrix.

This class is used in AdditiveSchwarz to reorder (if required by the user) the localized matrix. As the localized matrix is defined on a serial communicator only, all maps are trivial (as all elements reside on the same process). This class does not attemp to define properly reordered maps, hence it should not be used for distributed matrices.

To improve the performance of Ifpack2::AdditiveSchwarz, some operations are not performed in the construction phase (like for instance the computation of the 1-norm and infinite-norm, of check whether the reordered matrix is lower/upper triangular or not).

Constructor & Destructor Documentation

template<class MatrixType >
Ifpack2::ReorderFilter< MatrixType >::ReorderFilter ( const Teuchos::RCP< const row_matrix_type > &  A,
const Teuchos::ArrayRCP< local_ordinal_type > &  perm,
const Teuchos::ArrayRCP< local_ordinal_type > &  reverseperm 
)

Constructor.

Parameters
A[in] The matrix to which to apply the filter.
perm[in] Forward permutation of A's rows and columns.
reverseperm[in] Reverse permutation of A's rows and columns.

It must make sense to apply the given permutation to both the rows and columns. This means that the row and column Maps must have the same numbers of entries on all processes, and must have the same order of GIDs on all processes.

perm[i] gives the where OLD index i shows up in the NEW ordering. revperm[i] gives the where NEW index i shows up in the OLD ordering. Note that perm is actually the "inverse permutation," in Zoltan2 terms.

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

Destructor.

Member Function Documentation

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

The matrix's communicator.

template<class MatrixType >
Teuchos::RCP< const typename ReorderFilter< MatrixType >::map_type > Ifpack2::ReorderFilter< MatrixType >::getRowMap ( ) const
virtual

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

template<class MatrixType >
Teuchos::RCP< const typename ReorderFilter< MatrixType >::map_type > Ifpack2::ReorderFilter< MatrixType >::getColMap ( ) const
virtual

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

template<class MatrixType >
Teuchos::RCP< const typename ReorderFilter< MatrixType >::map_type > Ifpack2::ReorderFilter< MatrixType >::getDomainMap ( ) const
virtual

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

template<class MatrixType >
Teuchos::RCP< const typename ReorderFilter< MatrixType >::map_type > Ifpack2::ReorderFilter< 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::ReorderFilter< MatrixType >::getGraph ( ) const
virtual

Returns the RowGraph associated with this matrix.

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

Returns the number of global rows in this matrix.

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

Returns the number of global columns in this matrix.

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

Returns the number of rows owned on the calling node.

template<class MatrixType >
size_t Ifpack2::ReorderFilter< MatrixType >::getNodeNumCols ( ) 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::ReorderFilter< MatrixType >::getIndexBase ( ) const
virtual

Returns the index base for global indices for this matrix.

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

Returns the global number of entries in this matrix.

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

Returns the local number of entries in this matrix.

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

The current number of entries in this matrix, stored on the calling process, in the row whose global index is globalRow.

Returns
The number of entries, or Teuchos::OrdinalTraits<size_t>::invalid() if the specified row is not owned by the calling process.
template<class MatrixType >
size_t Ifpack2::ReorderFilter< MatrixType >::getNumEntriesInLocalRow ( local_ordinal_type  localRow) const
virtual

The current number of entries in this matrix, stored on the calling process, in the row whose local index is globalRow.

Returns
The number of entries, or Teuchos::OrdinalTraits<size_t>::invalid() if the specified row is not owned by the calling process.
template<class MatrixType >
size_t Ifpack2::ReorderFilter< MatrixType >::getGlobalMaxNumRowEntries ( ) const
virtual

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

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

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

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

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

template<class MatrixType >
bool Ifpack2::ReorderFilter< 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::ReorderFilter< 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::ReorderFilter< MatrixType >::isFillComplete ( ) const
virtual

Returns true if fillComplete() has been called.

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

Returns true if RowViews are supported.

template<class MatrixType >
void Ifpack2::ReorderFilter< MatrixType >::getGlobalRowCopy ( global_ordinal_type  GlobalRow,
const Teuchos::ArrayView< global_ordinal_type > &  Indices,
const Teuchos::ArrayView< scalar_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
GlobalRow- (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::ReorderFilter< MatrixType >::getLocalRowCopy ( local_ordinal_type  DropRow,
const Teuchos::ArrayView< local_ordinal_type > &  Indices,
const Teuchos::ArrayView< scalar_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
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().

template<class MatrixType >
void Ifpack2::ReorderFilter< 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
Postcondition
indices.size() == getNumEntriesInGlobalRow(GlobalRow)
Precondition
isLocallyIndexed() == false Note: If GlobalRow does not belong to this node, then indices is set to null.
template<class MatrixType >
void Ifpack2::ReorderFilter< 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) Global column indices corresponding to values.
Values- (Out) Row values
Precondition
isGloballyIndexed() == false
Postcondition
indices.size() == getNumEntriesInDropRow(LocalRow)

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

template<class MatrixType >
void Ifpack2::ReorderFilter< MatrixType >::getLocalDiagCopy ( Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  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::ReorderFilter< 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::ReorderFilter< 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 >
ReorderFilter< MatrixType >::mag_type Ifpack2::ReorderFilter< 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::ReorderFilter< 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

\( Y := \beta Y + \alpha Op(A) X \), where Op(A) is either A, \(A^T\), or \(A^H\).

Apply the reordered version of the matrix (or its transpose or conjugate transpose) to the given multivector X, producing Y.

  • 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 matrix, so that any values in X (including NaNs) are ignored.

This method assumes that X and Y are in the reordered order.

If hasTransposeApply() returns false, then the only valid value of mode is Teuchos::NO_TRANS (the default). Otherwise, it accepts the following values:

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

Whether apply() can apply the transpose or conjugate transpose.

template<class MatrixType >
void Ifpack2::ReorderFilter< MatrixType >::permuteOriginalToReordered ( const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  originalX,
Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  reorderedY 
) const
virtual

Permute multivector: original-to-reordered.

template<class MatrixType >
void Ifpack2::ReorderFilter< MatrixType >::permuteReorderedToOriginal ( const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  reorderedX,
Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  originalY 
) const
virtual

Permute multivector: reordered-to-original.


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