Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
Public Member Functions
Ifpack2::Container< MatrixType > Class Template Reference

Interface for creating and solving a local linear problem. More...

#include <Ifpack2_Container.hpp>

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

List of all members.

Public Member Functions

 Container (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::ArrayView< const local_ordinal_type > &localRows)
 Constructor.
virtual ~Container ()
 Destructor.
Teuchos::RCP< const
row_matrix_type > 
getMatrix () const
 The input matrix to the constructor.
virtual size_t getNumRows () const =0
 The number of rows in the diagonal block.
Teuchos::ArrayView< const
local_ordinal_type > 
getLocalRows () const
 Local indices of the rows of the input matrix that belong to this block.
virtual void initialize ()=0
 Do all set-up operations that only require matrix structure.
virtual void compute ()=0
 Extract the local diagonal block and prepare the solver.
virtual void setParameters (const Teuchos::ParameterList &List)=0
 Set parameters.
virtual bool isInitialized () const =0
 Return true if the container has been successfully initialized.
virtual bool isComputed () const =0
 Return true if the container has been successfully computed.
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 =0
 Compute Y := alpha * M^{-1} X + beta*Y.
virtual void weightedApply (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, const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &D, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const =0
 Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.
virtual std::ostream & print (std::ostream &os) const =0
 Print basic information about the container to os.

Detailed Description

template<class MatrixType>
class Ifpack2::Container< MatrixType >

Interface for creating and solving a local linear problem.

This class is mainly useful for the implementation of BlockRelaxation, and other preconditioners that need to solve linear systems with diagonal blocks of a sparse matrix.

Users of BlockRelaxation (and any analogous preconditioners) do not normally need to interact with the Container interface. However, they do need to specify a specific Container subclass to use, for example as the second template parameter (ContainerType) of BlockRelaxation. Implementations of Container specify

For example, the SparseContainer subclass uses a sparse matrix (in particular, Tpetra::CrsMatrix) to store each diagonal block, and can use any given Ifpack2 Preconditioner subclass to solve linear systems.

A Container can create, populate, and solve a local linear system. The local linear system matrix, B, is a submatrix of the local components of a distributed matrix, A. The idea of Container is to specify the rows of A that are contained in B, then extract B from A, and compute all it is necessary to solve a linear system in B. Then, set the initial guess (if necessary) and right-hand side for B, and solve the linear system in B.

If you are writing a class (comparable to BlockRelaxation) that uses Container, you should use it in the following way:

  1. Create a Container object, specifying the global matrix A and the indices of the local rows of A that are contained in B. The latter indices come from a Partitioner object.
  2. Optionally, set linear solve parameters using setParameters().
  3. Initialize the container by calling initialize().
  4. Prepare the linear system solver using compute().
  5. Solve the linear system using apply().

For an example of Steps 1-5 above, see the implementation of BlockRelaxation::ExtractSubmatrices() in Ifpack2_BlockRelaxation_def.hpp.


Constructor & Destructor Documentation

template<class MatrixType>
Ifpack2::Container< MatrixType >::Container ( const Teuchos::RCP< const row_matrix_type > &  matrix,
const Teuchos::ArrayView< const local_ordinal_type > &  localRows 
) [inline]

Constructor.

matrix [in] The original input matrix. This Container will construct a local diagonal block from the rows given by localRows.

Parameters:
localRows[in] The set of (local) rows assigned to this container. localRows[i] == j, where i (from 0 to getNumRows() - 1) indicates the Container's row, and j indicates the local row in the calling process. Subclasses must always pass along these indices to the base class.
template<class MatrixType>
virtual Ifpack2::Container< MatrixType >::~Container ( ) [inline, virtual]

Destructor.


Member Function Documentation

template<class MatrixType>
Teuchos::RCP<const row_matrix_type> Ifpack2::Container< MatrixType >::getMatrix ( ) const [inline]

The input matrix to the constructor.

This is not the local diagonal block that the Container extracts; this is the original input matrix, of which that diagonal block is a submatrix. You can't get the local diagonal block through the Container interface, because the subclass of container can store that block in any format it likes (e.g., dense or sparse).

template<class MatrixType>
virtual size_t Ifpack2::Container< MatrixType >::getNumRows ( ) const [pure virtual]
template<class MatrixType>
Teuchos::ArrayView<const local_ordinal_type> Ifpack2::Container< MatrixType >::getLocalRows ( ) const [inline]

Local indices of the rows of the input matrix that belong to this block.

The set of (local) rows assigned to this Container is defined by passing in a set of indices localRows[i] = j to the constructor, where i (from 0 to getNumRows() - 1) indicates the Container's row, and j indicates the local row in the calling process. Subclasses must always pass along these indices to the base class.

The indices are usually used to reorder the local row index (on the calling process) of the i-th row in the Container.

For an example of how to use these indices, see the implementation of BlockRelaxation::ExtractSubmatrices() in Ifpack2_BlockRelaxation_def.hpp.

template<class MatrixType>
virtual void Ifpack2::Container< MatrixType >::initialize ( ) [pure virtual]

Do all set-up operations that only require matrix structure.

If the input matrix's structure changes, you must call this method before you may call compute(). You must then call compute() before you may call apply() or weightedApply().

"Structure" refers to the graph of the matrix: the local and global dimensions, and the populated entries in each row.

Implemented in Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::DenseContainer< MatrixType, LocalScalarType >, and Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.

template<class MatrixType>
virtual void Ifpack2::Container< MatrixType >::compute ( ) [pure virtual]

Extract the local diagonal block and prepare the solver.

If any entries' values in the input matrix have changed, you must call this method before you may call apply() or weightedApply().

Implemented in Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::DenseContainer< MatrixType, LocalScalarType >, and Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.

template<class MatrixType>
virtual void Ifpack2::Container< MatrixType >::setParameters ( const Teuchos::ParameterList List) [pure virtual]
template<class MatrixType>
virtual bool Ifpack2::Container< MatrixType >::isInitialized ( ) const [pure virtual]
template<class MatrixType>
virtual bool Ifpack2::Container< MatrixType >::isComputed ( ) const [pure virtual]
template<class MatrixType>
virtual void Ifpack2::Container< 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 [pure virtual]

Compute Y := alpha * M^{-1} X + beta*Y.

X is in the domain Map of the original matrix (the argument to compute()), and Y is in the range Map of the original matrix. This method only reads resp. modifies the permuted subset of entries of X resp. Y related to the diagonal block M. That permuted subset is defined by the indices passed into the constructor.

This method is marked const for compatibility with Tpetra::Operator's method of the same name. This might require subclasses to mark some of their instance data as mutable.

Implemented in Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::DenseContainer< MatrixType, LocalScalarType >, and Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.

template<class MatrixType>
virtual void Ifpack2::Container< MatrixType >::weightedApply ( 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,
const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &  D,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
scalar_type  alpha = Teuchos::ScalarTraits< scalar_type >::one(),
scalar_type  beta = Teuchos::ScalarTraits< scalar_type >::zero() 
) const [pure virtual]

Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.

X is in the domain Map of the original matrix (the argument to compute()), and Y is in the range Map of the original matrix. This method only reads resp. modifies the permuted subset of entries of X resp. Y related to the diagonal block M. That permuted subset is defined by the indices passed into the constructor. The D scaling vector must have the same number of entries on each process as X and Y, but otherwise need not have the same Map. (For example, D could be locally replicated, or could be a different object on each process with a local (MPI_COMM_SELF) communicator.)

This method supports overlap techniques, such as those used in Schwarz methods.

This method is marked const by analogy with apply(), which itself is marked const for compatibility with Tpetra::Operator's method of the same name. This might require subclasses to mark some of their instance data as mutable.

Implemented in Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::DenseContainer< MatrixType, LocalScalarType >, and Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.

template<class MatrixType>
virtual std::ostream& Ifpack2::Container< MatrixType >::print ( std::ostream &  os) const [pure virtual]

The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends