Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Static Public Member Functions | List of all members
Ifpack2::SparseContainer< MatrixType, InverseType > Class Template Reference

Store and solve a local sparse linear problem. More...

#include <Ifpack2_SparseContainer_decl.hpp>

Inheritance diagram for Ifpack2::SparseContainer< MatrixType, InverseType >:
Inheritance graph
[legend]

Public Member Functions

Constructor and destructor
 SparseContainer (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions, const Teuchos::RCP< const import_type > &importer, int OverlapLevel, scalar_type DampingFactor)
 Constructor. More...
 
 SparseContainer (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< local_ordinal_type > &localRows)
 
virtual ~SparseContainer ()
 Destructor (declared virtual for memory safety of derived classes). More...
 
Get and set methods
virtual bool isInitialized () const
 Whether the container has been successfully initialized. More...
 
virtual bool isComputed () const
 Whether the container has been successfully computed. More...
 
virtual void setParameters (const Teuchos::ParameterList &List)
 Set all necessary parameters. More...
 
Mathematical functions
virtual void initialize ()
 Do all set-up operations that only require matrix structure. More...
 
virtual void compute ()
 Initialize and compute all blocks. More...
 
void clearBlocks ()
 
virtual void apply (HostView &X, HostView &Y, int blockIndex, int stride, 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 := alpha * M^{-1} X + beta*Y. More...
 
virtual void weightedApply (HostView &X, HostView &Y, HostView &W, int blockIndex, int stride, 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 := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y. More...
 
Miscellaneous methods
virtual std::ostream & print (std::ostream &os) const
 Print information about this object to the given output stream. More...
 
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 with some verbosity level to the given FancyOStream. More...
 
- Public Member Functions inherited from Ifpack2::Container< MatrixType >
 Container (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions, const Teuchos::RCP< const import_type > &importer, int OverlapLevel, scalar_type DampingFactor)
 Constructor. More...
 
 Container (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< local_ordinal_type > &localRows)
 Constructor for single block. More...
 
virtual ~Container ()
 Destructor. More...
 
Teuchos::ArrayView< const
local_ordinal_type > 
getLocalRows (int blockIndex) const
 Local indices of the rows of the input matrix that belong to this block. More...
 
void setBlockSizes (const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions)
 Initialize arrays with information about block sizes. More...
 
virtual void applyInverseJacobi (const mv_type &, mv_type &, bool, int) const
 Compute Y := (1 - a) Y + a D^{-1} (X - R*Y). More...
 
void applyMV (mv_type &X, mv_type &Y) const
 Wrapper for apply with MVs, used in unit tests (never called by BlockRelaxation) More...
 
void weightedApplyMV (mv_type &X, mv_type &Y, vector_type &W)
 Wrapper for weightedApply with MVs, used in unit tests (never called by BlockRelaxation) More...
 

Static Public Member Functions

static std::string getName ()
 Get the name of this container type for Details::constructContainer() More...
 
- Static Public Member Functions inherited from Ifpack2::Container< MatrixType >
static std::string getName ()
 

Additional Inherited Members

- Protected Types inherited from Ifpack2::Container< MatrixType >
typedef MatrixType::scalar_type scalar_type
 
typedef
MatrixType::local_ordinal_type 
local_ordinal_type
 
typedef
MatrixType::global_ordinal_type 
global_ordinal_type
 
typedef MatrixType::node_type node_type
 
typedef Tpetra::MultiVector
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type > 
mv_type
 
typedef Tpetra::Vector
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type > 
vector_type
 
typedef Tpetra::Map
< local_ordinal_type,
global_ordinal_type, node_type > 
map_type
 
typedef Teuchos::ScalarTraits
< scalar_type > 
STS
 
typedef Tpetra::Import
< local_ordinal_type,
global_ordinal_type, node_type > 
import_type
 
typedef Partitioner
< Tpetra::RowGraph
< local_ordinal_type,
global_ordinal_type, node_type > > 
partitioner_type
 
typedef
Tpetra::Experimental::BlockCrsMatrix
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type > 
block_crs_matrix_type
 
typedef Tpetra::RowMatrix
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type > 
row_matrix_type
 
typedef
Kokkos::Details::ArithTraits
< scalar_type >::val_type 
impl_scalar_type
 Internal representation of Scalar in Kokkos::View. More...
 
- Protected Attributes inherited from Ifpack2::Container< MatrixType >
Teuchos::RCP< const
row_matrix_type > 
inputMatrix_
 The input matrix to the constructor. More...
 
int numBlocks_
 The number of blocks (partitions) in the container. More...
 
Teuchos::Array
< local_ordinal_type > 
partitions_
 Local indices of the rows of the input matrix that belong to this block. More...
 
Teuchos::Array
< local_ordinal_type > 
blockRows_
 Number of rows in each block. More...
 
Teuchos::Array
< local_ordinal_type > 
partitionIndices_
 Starting index in partitions_ of local row indices for each block. More...
 
Teuchos::RCP< vector_type > Diag_
 Diagonal elements. More...
 
bool IsParallel_
 Whether the problem is distributed across multiple MPI processes. More...
 
int OverlapLevel_
 Number of rows of overlap for adjacent blocks. More...
 
scalar_type DampingFactor_
 Damping factor, passed to apply() as alpha. More...
 
Teuchos::RCP< const
Tpetra::Import
< local_ordinal_type,
global_ordinal_type, node_type > > 
Importer_
 Importer for importing off-process elements of MultiVectors. More...
 
local_ordinal_type NumLocalRows_
 Number of local rows in input matrix. More...
 
global_ordinal_type NumGlobalRows_
 Number of global rows in input matrix. More...
 
global_ordinal_type NumGlobalNonzeros_
 Number of nonzeros in input matrix. More...
 
bool hasBlockCrs_
 Whether the input matrix is a BlockCRS matrix. More...
 
int bcrsBlockSize_
 If hasBlockCrs_, the number of DOFs per vertex. Otherwise 1. More...
 

Detailed Description

template<typename MatrixType, typename InverseType>
class Ifpack2::SparseContainer< MatrixType, InverseType >

Store and solve a local sparse linear problem.

Template Parameters
Aspecialization of Tpetra::RowMatrix.

Please refer to the documentation of the Container interface. Currently, Containers are used by BlockRelaxation. Block relaxations need to be able to do two things:

  1. Store the diagonal blocks
  2. Solve linear systems with each diagonal block

These correspond to the two template parameters:

  1. MatrixType, which stores a sparse matrix
  2. InverseType, which solves linear systems with that matrix

This class stores each block as a sparse matrix. Using a sparse matrix for each block is a good idea when the blocks are large and sparse. For small and / or dense blocks, it would probably be better to use an implementation of Container that stores the blocks densely, like DenseContainer. You may also want to consider BandedContainer.

The InverseType template parameter represents the class to use for solving linear systems with a block. In SparseContainer, this template parameter must be a specialization of Preconditioner. Specifically, InverseType must implement the following methods:

We also assume that InverseType has the following typedefs:

MatrixType and InverseType may store values of different types, and may have different template parameters (e.g., local or global ordinal types). You may mix and match so long as implicit conversions are available. The most obvious use case for this are:

SparseContainer currently assumes the following about the column and row Maps of the input matrix:

  1. On all processes, the column and row Maps begin with the same set of on-process entries, in the same order. That is, on-process row and column indices are the same.
  2. On all processes, all off-process indices in the column Map of the input matrix occur after that initial set.

These assumptions may be violated if the input matrix is a Tpetra::CrsMatrix that was constructed with a user-provided column Map. The assumptions are not mathematically necessary and could be relaxed at any time. Implementers who wish to do so will need to modify the extract() method, so that it translates explicitly between local row and column indices, instead of just assuming that they are the same.

Constructor & Destructor Documentation

template<class MatrixType , class InverseType >
Ifpack2::SparseContainer< MatrixType, InverseType >::SparseContainer ( const Teuchos::RCP< const row_matrix_type > &  matrix,
const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &  partitions,
const Teuchos::RCP< const import_type > &  importer,
int  OverlapLevel,
scalar_type  DampingFactor 
)

Constructor.

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

partitioner [in] The BlockRelaxation partitioner.

template<class MatrixType , class InverseType >
Ifpack2::SparseContainer< MatrixType, InverseType >::~SparseContainer ( )
virtual

Destructor (declared virtual for memory safety of derived classes).

Member Function Documentation

template<class MatrixType , class InverseType >
bool Ifpack2::SparseContainer< MatrixType, InverseType >::isInitialized ( ) const
virtual

Whether the container has been successfully initialized.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
bool Ifpack2::SparseContainer< MatrixType, InverseType >::isComputed ( ) const
virtual

Whether the container has been successfully computed.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
void Ifpack2::SparseContainer< MatrixType, InverseType >::setParameters ( const Teuchos::ParameterList List)
virtual

Set all necessary parameters.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
void Ifpack2::SparseContainer< MatrixType, InverseType >::initialize ( )
virtual

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

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
void Ifpack2::SparseContainer< MatrixType, InverseType >::compute ( )
virtual

Initialize and compute all blocks.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
void Ifpack2::SparseContainer< MatrixType, InverseType >::clearBlocks ( )
virtual

Free all per-block resources: Inverses_, and diagBlocks_. Called by BlockRelaxation when the input matrix is changed. Also calls Container::clearBlocks()

Reimplemented from Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
void Ifpack2::SparseContainer< MatrixType, InverseType >::apply ( HostView &  X,
HostView &  Y,
int  blockIndex,
int  stride,
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 := alpha * M^{-1} X + beta*Y.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
void Ifpack2::SparseContainer< MatrixType, InverseType >::weightedApply ( HostView &  X,
HostView &  Y,
HostView &  W,
int  blockIndex,
int  stride,
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 := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
std::ostream & Ifpack2::SparseContainer< MatrixType, InverseType >::print ( std::ostream &  os) const
virtual

Print information about this object to the given output stream.

operator<< uses this method.

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class InverseType >
std::string Ifpack2::SparseContainer< MatrixType, InverseType >::description ( ) const
virtual

A one-line description of this object.

Reimplemented from Teuchos::Describable.

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

Print the object with some verbosity level to the given FancyOStream.

Reimplemented from Teuchos::Describable.

template<typename MatrixType , typename InverseType >
std::string Ifpack2::SparseContainer< MatrixType, InverseType >::getName ( )
static

Get the name of this container type for Details::constructContainer()


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