Ifpack2 Templated Preconditioning Package
Version 1.0
|
Interface for creating and solving a set of local linear problems. More...
#include <Ifpack2_Container_decl.hpp>
Public Member Functions | |
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... | |
virtual void | initialize ()=0 |
Do all set-up operations that only require matrix structure. 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... | |
virtual void | compute ()=0 |
Extract the local diagonal blocks and prepare the solver. More... | |
virtual void | setParameters (const Teuchos::ParameterList &List)=0 |
Set parameters, if any. More... | |
virtual void | apply (ConstHostView X, HostView Y, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const =0 |
Compute Y := alpha * M^{-1} X + beta*Y . More... | |
virtual void | weightedApply (ConstHostView X, HostView Y, ConstHostView D, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const =0 |
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 =0 |
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y) . More... | |
virtual void | applyMV (const mv_type &X, mv_type &Y) const |
Wrapper for apply with MultiVector. More... | |
virtual void | weightedApplyMV (const mv_type &X, mv_type &Y, vector_type &W) const |
Wrapper for weightedApply with MultiVector. More... | |
virtual std::ostream & | print (std::ostream &os) const =0 |
Print basic information about the container to os . More... | |
Static Public Member Functions | |
static std::string | getName () |
Protected Types | |
using | ISC = typename Kokkos::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 | |
virtual void | DoGSBlock (ConstHostView 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... | |
Protected Attributes | |
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... | |
Interface for creating and solving a set of local linear problems.
MatrixType | A specialization of Tpetra::RowMatrix. |
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 sparse matrices (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 local linear systems. A local linear system matrix, B, is a submatrix of the local components of the distributed matrix, A. The idea of Container is to specify the rows of A that are contained in each 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 sides for B, and solve the linear systems in B.
If you are writing a class (comparable to BlockRelaxation) that uses Container, you should use it in the following way:
For an example of Steps 1-5 above, see ExtractSubmatrices() and ApplyInverseGS() in Ifpack2_BlockRelaxation_def.hpp.
|
protected |
Internal representation of Scalar in Kokkos::View.
|
protected |
HostView (the host-space internal representation for Tpetra::Multivector) is the type of the vector arguments of DoJacobi, DoGaussSeidel, and DoSGS.
Ifpack2::Container< MatrixType >::Container | ( | const Teuchos::RCP< const row_matrix_type > & | matrix, |
const Teuchos::Array< Teuchos::Array< LO > > & | partitions, | ||
bool | pointIndexed | ||
) |
Constructor.
matrix | [in] The original input matrix. This Container will construct local diagonal blocks from its rows according to partitions . |
partitioner | [in] The Partitioner object that assigns local rows of the input matrix to blocks. |
pointIndexed | [in] If the input matrix is a Tpetra::BlockCrsMatrix , whether elements of partitions [k] identify rows within blocks (true) or whole blocks (false). |
|
virtual |
Destructor.
Teuchos::ArrayView< const typename MatrixType::local_ordinal_type > Ifpack2::Container< MatrixType >::getBlockRows | ( | int | blockIndex | ) | const |
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 blockRows[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.
|
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::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplNotAvailTag >, Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >, Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >, Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::BandedContainer< MatrixType, LocalScalarType >, Ifpack2::DenseContainer< MatrixType, LocalScalarType >, and Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.
void Ifpack2::Container< MatrixType >::setBlockSizes | ( | const Teuchos::Array< Teuchos::Array< LO > > & | partitions | ) |
Initialize arrays with information about block sizes.
bool Ifpack2::Container< MatrixType >::isInitialized | ( | ) | const |
Whether the container has been successfully initialized.
bool Ifpack2::Container< MatrixType >::isComputed | ( | ) | const |
Whether the container has been successfully computed.
|
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().
Implemented in Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplNotAvailTag >, Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >, Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >, Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::BandedContainer< MatrixType, LocalScalarType >, Ifpack2::DenseContainer< MatrixType, LocalScalarType >, and Ifpack2::TriDiContainer< MatrixType, LocalScalarType >.
|
pure virtual |
Set parameters, if any.
Implemented in Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplNotAvailTag >, Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >, Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >, Ifpack2::SparseContainer< MatrixType, InverseType >, and Ifpack2::BandedContainer< MatrixType, LocalScalarType >.
|
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::ContainerImpl< MatrixType, LocalScalarType >, Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >, and Ifpack2::SparseContainer< MatrixType, InverseType >.
|
pure virtual |
Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y
.
Implemented in Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >, and Ifpack2::SparseContainer< MatrixType, InverseType >.
|
pure virtual |
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y)
.
Implemented in Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, and Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Wrapper for apply with MultiVector.
Reimplemented in Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, and Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Wrapper for weightedApply with MultiVector.
Reimplemented in Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, and Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
pure virtual |
Print basic information about the container to os
.
Implemented in Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >, Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplNotAvailTag >, Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >, Ifpack2::SparseContainer< MatrixType, InverseType >, Ifpack2::BandedContainer< MatrixType, LocalScalarType >, Ifpack2::TriDiContainer< MatrixType, LocalScalarType >, and Ifpack2::DenseContainer< MatrixType, LocalScalarType >.
|
static |
Returns string describing the container. See Details::ContainerFactory
.
|
protectedvirtual |
Do one step of Gauss-Seidel on block i (used by DoGaussSeidel and DoSGS)
Reimplemented in Ifpack2::ContainerImpl< MatrixType, LocalScalarType >, and Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
protected |
The input matrix to the constructor.
|
protected |
The input matrix, dynamic cast to CrsMatrix. May be null.
|
protected |
The input matrix, dynamic cast to BlockCrsMatrix. May be null.
|
protected |
The number of blocks (partitions) in the container.
|
protected |
Local indices of the rows of the input matrix that belong to this block.
|
protected |
Number of rows in each block.
|
protected |
Starting index in blockRows_ of local row indices for each block.
|
mutableprotected |
Diagonal elements.
|
protected |
Whether the problem is distributed across multiple MPI processes.
|
protected |
Number of local rows in input matrix.
|
protected |
Number of global rows in input matrix.
|
protected |
Number of nonzeros in input matrix.
|
protected |
Whether the input matrix is a BlockCRS matrix.
|
protected |
If hasBlockCrs_, the number of DOFs per vertex. Otherwise 1.
|
protected |
(If hasBlockCrs_) Whether the blocks are described using sub-block row indices instead of full block rows.
|
protected |
Number of scalars corresponding to one element of blockRows_ If pointIndexed_, always 1. Otherwise if hasBlockCrs_, then bcrsBlockSize_. Otherwise 1.
|
protected |
If true
, the container has been successfully initialized.
|
protected |
If true
, the container has been successfully computed.