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

The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS). This class allows a custom scalar type (LocalScalarType) to be used for storing blocks and solving the block systems. Hiding this template parameter from the Container interface simplifies the BlockRelaxation and ContainerFactory classes. More...

#include <Ifpack2_Container_decl.hpp>

Inheritance diagram for Ifpack2::ContainerImpl< MatrixType, LocalScalarType >:
Inheritance graph
[legend]

Public Member Functions

virtual ~ContainerImpl ()
 Destructor. More...
 
virtual void initialize ()=0
 Do all set-up operations that only require matrix structure. More...
 
virtual void compute ()=0
 Extract the local diagonal blocks and prepare the solver. More...
 
virtual void setParameters (const Teuchos::ParameterList &List)
 Set parameters, if any. More...
 
virtual void apply (HostView X, HostView Y, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const
 Compute Y := alpha * M^{-1} X + beta*Y. More...
 
virtual void weightedApply (HostView X, HostView Y, HostView D, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const
 Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y. More...
 
virtual void applyInverseJacobi (const mv_type &, mv_type &, SC dampingFactor, 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) const
 Wrapper for weightedApply with MVs, used in unit tests (never called by BlockRelaxation) More...
 
virtual std::ostream & print (std::ostream &os) const =0
 Print basic information about the container to os. More...
 
- Public Member Functions inherited from Ifpack2::Container< MatrixType >
 Container (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, bool pointIndexed)
 Constructor. More...
 
virtual ~Container ()
 Destructor. More...
 
Teuchos::ArrayView< const LO > getBlockRows (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< LO > > &partitions)
 Initialize arrays with information about block sizes. More...
 
bool isInitialized () const
 Whether the container has been successfully initialized. More...
 
bool isComputed () const
 Whether the container has been successfully computed. More...
 

Static Public Member Functions

static std::string getName ()
 
- Static Public Member Functions inherited from Ifpack2::Container< MatrixType >
static std::string getName ()
 

Protected Types

Internal typedefs (protected)
using local_scalar_type = LocalScalarType
 
using SC = typename Container< MatrixType >::scalar_type
 
using LO = typename Container< MatrixType >::local_ordinal_type
 
using GO = typename Container< MatrixType >::global_ordinal_type
 
using NO = typename Container< MatrixType >::node_type
 
using StridedRowView = Details::StridedRowView< SC, LO, GO, NO >
 
using LSC = LocalScalarType
 The internal representation of LocalScalarType in Kokkos::View. More...
 
using LISC = typename Kokkos::Details::ArithTraits< LSC >::val_type
 
using local_mv_type = Tpetra::MultiVector< LSC, LO, GO, NO >
 
using HostViewLocal = typename local_mv_type::dual_view_type::t_host
 
using HostSubviewLocal = Kokkos::View< LISC **, Kokkos::LayoutStride, typename HostViewLocal::memory_space >
 
- Protected Types inherited from Ifpack2::Container< MatrixType >
using ISC = typename Kokkos::Details::ArithTraits< SC >::val_type
 Internal representation of Scalar in Kokkos::View. More...
 
using HostView = typename mv_type::dual_view_type::t_host
 

Protected Member Functions

LO translateRowToCol (LO row)
 
Details::StridedRowView< SC,
LO, GO, NO > 
getInputRowView (LO row) const
 View a row of the input matrix. More...
 
virtual void extract ()=0
 Extract the submatrices identified by the local indices set by the constructor. BlockMatrix may be any type that supports direct entry access: Scalar& operator()(size_t row, size_t col). More...
 
void DoGSBlock (HostView X, HostView Y, HostView Y2, HostView Resid, SC dampingFactor, LO i) const
 Do one step of Gauss-Seidel on block i (used by DoGaussSeidel and DoSGS) More...
 
virtual void solveBlock (HostSubviewLocal X, HostSubviewLocal Y, int blockIndex, Teuchos::ETransp mode, const LSC alpha, const LSC beta) const
 

Protected Attributes

HostViewLocal X_local_
 Scratch vectors used in apply(). More...
 
HostViewLocal weightedApplyScratch_
 
std::vector< HostSubviewLocal > X_localBlocks_
 Views for holding pieces of X corresponding to each block. More...
 
std::vector< HostSubviewLocal > Y_localBlocks_
 Views for holding pieces of Y corresponding to each block. More...
 
- Protected Attributes inherited from Ifpack2::Container< MatrixType >
Teuchos::RCP< const
row_matrix_type > 
inputMatrix_
 The input matrix to the constructor. More...
 
Teuchos::RCP< const
crs_matrix_type > 
inputCrsMatrix_
 The input matrix, dynamic cast to CrsMatrix. May be null. More...
 
Teuchos::RCP< const
block_crs_matrix_type > 
inputBlockMatrix_
 The input matrix, dynamic cast to BlockCrsMatrix. May be null. More...
 
int numBlocks_
 The number of blocks (partitions) in the container. More...
 
Teuchos::Array< LO > blockRows_
 Local indices of the rows of the input matrix that belong to this block. More...
 
Teuchos::Array< LO > blockSizes_
 Number of rows in each block. More...
 
Teuchos::Array< LO > blockOffsets_
 Starting index in blockRows_ 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...
 
LO NumLocalRows_
 Number of local rows in input matrix. More...
 
GO NumGlobalRows_
 Number of global rows in input matrix. More...
 
GO 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...
 
bool pointIndexed_
 (If hasBlockCrs_) Whether the blocks are described using sub-block row indices instead of full block rows. More...
 
LO scalarsPerRow_
 
bool IsInitialized_
 If true, the container has been successfully initialized. More...
 
bool IsComputed_
 If true, the container has been successfully computed. More...
 

Detailed Description

template<class MatrixType, class LocalScalarType>
class Ifpack2::ContainerImpl< MatrixType, LocalScalarType >

The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS). This class allows a custom scalar type (LocalScalarType) to be used for storing blocks and solving the block systems. Hiding this template parameter from the Container interface simplifies the BlockRelaxation and ContainerFactory classes.

Member Typedef Documentation

template<class MatrixType, class LocalScalarType>
using Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::LSC = LocalScalarType
protected

The internal representation of LocalScalarType in Kokkos::View.

Constructor & Destructor Documentation

template<class MatrixType , class LocalScalarType >
Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::~ContainerImpl ( )
virtual

Destructor.

Member Function Documentation

template<class MatrixType , class LocalScalarType >
ContainerImpl< MatrixType, LocalScalarType >::LO Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::translateRowToCol ( LO  row)
protected

Translate local row (if pointIndexed_, then a row within a block row/node) to the corresponding column, so that the intersection of the row and column is on the global diagonal. Translates to global ID first using the row map, and then back to a local column ID using the col map. If the corresponding column is off-process, then returns OrdinalTraits<LO>::invalid() and does not throw an exception.

template<class MatrixType , class LocalScalarType >
Details::StridedRowView< typename ContainerImpl< MatrixType, LocalScalarType >::SC, typename ContainerImpl< MatrixType, LocalScalarType >::LO, typename ContainerImpl< MatrixType, LocalScalarType >::GO, typename ContainerImpl< MatrixType, LocalScalarType >::NO > Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::getInputRowView ( LO  row) const
protected

View a row of the input matrix.

template<class MatrixType, class LocalScalarType>
virtual void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::extract ( )
protectedpure virtual

Extract the submatrices identified by the local indices set by the constructor. BlockMatrix may be any type that supports direct entry access: Scalar& operator()(size_t row, size_t col).

Parameters
diagBlocks[in/out] The diagonal block matrices. Its size must be this->numBlocks_. Each BlockMatrix must be ready for entries to be assigned.
template<class MatrixType, class LocalScalarType>
virtual void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::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.

Implements Ifpack2::Container< MatrixType >.

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

template<class MatrixType, class LocalScalarType>
virtual void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::compute ( )
pure virtual

Extract the local diagonal blocks 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().

If DOF decoupling is to be used, it must be enabled with enableDecoupling() before calling compute().

Implements Ifpack2::Container< MatrixType >.

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

template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::setParameters ( const Teuchos::ParameterList List)
virtual
template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::apply ( HostView  X,
HostView  Y,
int  blockIndex,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
SC  alpha = Teuchos::ScalarTraits<SC>::one(),
SC  beta = Teuchos::ScalarTraits<SC>::zero() 
) const
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.

Implements Ifpack2::Container< MatrixType >.

Reimplemented in Ifpack2::SparseContainer< MatrixType, InverseType >.

template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::weightedApply ( HostView  X,
HostView  Y,
HostView  D,
int  blockIndex,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
SC  alpha = Teuchos::ScalarTraits<SC>::one(),
SC  beta = Teuchos::ScalarTraits<SC>::zero() 
) const
virtual

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

Implements Ifpack2::Container< MatrixType >.

Reimplemented in Ifpack2::SparseContainer< MatrixType, InverseType >.

template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::applyInverseJacobi ( const mv_type &  ,
mv_type &  ,
SC  dampingFactor,
bool  ,
int   
) const
virtual

Compute Y := (1 - a) Y + a D^{-1} (X - R*Y).

Implements Ifpack2::Container< MatrixType >.

template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::applyMV ( mv_type &  X,
mv_type &  Y 
) const
virtual

Wrapper for apply with MVs, used in unit tests (never called by BlockRelaxation)

Reimplemented from Ifpack2::Container< MatrixType >.

template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::weightedApplyMV ( mv_type &  X,
mv_type &  Y,
vector_type &  W 
) const
virtual

Wrapper for weightedApply with MVs, used in unit tests (never called by BlockRelaxation)

Reimplemented from Ifpack2::Container< MatrixType >.

template<class MatrixType, class LocalScalarType>
virtual std::ostream& Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::print ( std::ostream &  os) const
pure virtual
template<class MatrixType , class LocalScalarType >
std::string Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::getName ( )
static

Returns string describing the container. See Details::ContainerFactory.

template<class MatrixType , typename LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::DoGSBlock ( HostView  X,
HostView  Y,
HostView  Y2,
HostView  Resid,
SC  dampingFactor,
LO  i 
) const
protectedvirtual

Do one step of Gauss-Seidel on block i (used by DoGaussSeidel and DoSGS)

Reimplemented from Ifpack2::Container< MatrixType >.

template<class MatrixType , class LocalScalarType >
void Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::solveBlock ( HostSubviewLocal  X,
HostSubviewLocal  Y,
int  blockIndex,
Teuchos::ETransp  mode,
const LSC  alpha,
const LSC  beta 
) const
protectedvirtual

Exactly solves the linear system By = x, where B is a diagonal block matrix (blockIndex), and X, Y are multivector subviews with the same length as B's dimensions.

The Dense, Banded and TriDi containers all implement this and it is used in ContainerImpl::apply(). The Sparse and BlockTriDi containers have their own implementation of apply() and do not use solveBlock.

Reimplemented in Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.

Member Data Documentation

template<class MatrixType, class LocalScalarType>
HostViewLocal Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::X_local_
mutableprotected

Scratch vectors used in apply().

template<class MatrixType, class LocalScalarType>
HostViewLocal Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::weightedApplyScratch_
mutableprotected

applyScratch provides space for temporary block-sized vectors in weightedApply(), so that full Kokkos::Views don't need to be created for every block apply (slow).

Layout of applyScratch_:

Name Range
D_local 0 to maxBlockSize_
X_scaled maxBlockSize_ to 2 * maxBlockSize_
Y_temp 2 * maxBlockSize_ to 3 * maxBlockSize_
template<class MatrixType, class LocalScalarType>
std::vector<HostSubviewLocal> Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::X_localBlocks_
mutableprotected

Views for holding pieces of X corresponding to each block.

template<class MatrixType, class LocalScalarType>
std::vector<HostSubviewLocal> Ifpack2::ContainerImpl< MatrixType, LocalScalarType >::Y_localBlocks_
mutableprotected

Views for holding pieces of Y corresponding to each block.


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