Ifpack2 Templated Preconditioning Package
Version 1.0
|
Store and solve a local sparse linear problem. More...
#include <Ifpack2_SparseContainer_decl.hpp>
Public Member Functions | |
Constructor and destructor | |
SparseContainer (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed) | |
Constructor. More... | |
virtual | ~SparseContainer () |
Destructor (declared virtual for memory safety of derived classes). More... | |
Get and set methods | |
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, 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 W, 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... | |
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::ContainerImpl< MatrixType, InverseType::scalar_type > | |
virtual | ~ContainerImpl () |
Destructor. 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... | |
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 () |
Get the name of this container type for Details::constructContainer() More... | |
Static Public Member Functions inherited from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type > | |
static std::string | getName () |
Static Public Member Functions inherited from Ifpack2::Container< MatrixType > | |
static std::string | getName () |
Additional Inherited Members | |
Protected Types inherited from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type > | |
using | local_scalar_type = InverseType::scalar_type |
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 = InverseType::scalar_type |
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 inherited from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type > | |
LO | translateRowToCol (LO row) |
Details::StridedRowView< SC, LO, GO, NO > | getInputRowView (LO row) const |
View a row of the input matrix. 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 inherited from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type > | |
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... | |
Store and solve a local sparse linear problem.
A | specialization 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:
These correspond to the two template parameters:
MatrixType
, which stores a sparse matrix 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:
RCP<const MatrixType>
setParameters(Teuchos::ParameterList&)
initialize()
compute()
apply (const mv_type& X, mv_type& Y, ...)
, where mv_type
is the appropriate specialization of Tpetra::MultiVector We also assume that InverseType
has the following typedefs:
scalar_type
local_ordinal_type
global_ordinal_type
node_type
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:
MatrixType::global_ordinal_type=long long
and InverseType::global_ordinal_type=short
MatrixType::scalar_type=float
and InverseType::scalar_type=double
SparseContainer currently assumes the following about the column and row Maps of the input matrix:
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.
Ifpack2::SparseContainer< MatrixType, InverseType >::SparseContainer | ( | const Teuchos::RCP< const row_matrix_type > & | matrix, |
const Teuchos::Array< Teuchos::Array< LO > > & | partitions, | ||
const Teuchos::RCP< const import_type > & | importer, | ||
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 (declared virtual for memory safety of derived classes).
|
virtual |
Set all necessary parameters.
Reimplemented from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Do all set-up operations that only require matrix structure.
Implements Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Initialize and compute all blocks.
Implements Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
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::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Compute Y := alpha * M^{-1} X + beta*Y
.
Reimplemented from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y
.
Reimplemented from Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
Print information about this object to the given output stream.
operator<< uses this method.
Implements Ifpack2::ContainerImpl< MatrixType, InverseType::scalar_type >.
|
virtual |
A one-line description of this object.
Reimplemented from Teuchos::Describable.
|
virtual |
Print the object with some verbosity level to the given FancyOStream.
Reimplemented from Teuchos::Describable.
|
static |
Get the name of this container type for Details::constructContainer()