Ifpack2 Templated Preconditioning Package
Version 1.0
|
Overlapping Schwarz where redundant patches are not stored explicitly. More...
#include <Ifpack2_DatabaseSchwarz_decl.hpp>
Public Types | |
Typedefs | |
typedef MatrixType | matrix_type |
The template parameter of this class. More... | |
typedef MatrixType::scalar_type | scalar_type |
The type of the entries of the input MatrixType. More... | |
typedef MatrixType::local_ordinal_type | local_ordinal_type |
The type of local indices in the input MatrixType. More... | |
typedef MatrixType::global_ordinal_type | global_ordinal_type |
The type of global indices in the input MatrixType. More... | |
typedef MatrixType::node_type::device_type | device_type |
The Kokkos::Device specialization used by the input MatrixType. More... | |
typedef MatrixType::node_type | node_type |
The Node type used by the input MatrixType. More... | |
typedef Teuchos::ScalarTraits < scalar_type >::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
typedef Tpetra::RowMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | row_matrix_type |
The Tpetra::RowMatrix specialization matching MatrixType. More... | |
typedef Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > | map_type |
The Tpetra::Map specialization matching MatrixType. More... | |
typedef Tpetra::Vector < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | vector_type |
The Tpetra::Vector specialization matching MatrixType. More... | |
Public Types inherited from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > | |
typedef Teuchos::ScalarTraits < MatrixType::scalar_type > ::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
Public Member Functions | |
DatabaseSchwarz (const Teuchos::RCP< const row_matrix_type > &A) | |
Constructor. More... | |
DatabaseSchwarz (const Teuchos::RCP< const row_matrix_type > &A, Teuchos::ParameterList ¶ms) | |
Constructor. More... | |
virtual | ~DatabaseSchwarz () |
Destructor. More... | |
Preconditioner computation methods | |
void | setParameters (const Teuchos::ParameterList ¶ms) |
Set (or reset) parameters. More... | |
bool | supportsZeroStartingSolution () |
void | setZeroStartingSolution (bool zeroStartingSolution) |
Set this preconditioner's parameters. More... | |
void | initialize () |
Initialize the preconditioner. More... | |
bool | isInitialized () const |
void | compute () |
(Re)compute the left scaling, and (if applicable) estimate max and min eigenvalues of D_inv * A. More... | |
bool | isComputed () const |
Whether compute() has been called at least once. More... | |
Implementation of Ifpack2::Details::CanChangeMatrix | |
virtual void | setMatrix (const Teuchos::RCP< const row_matrix_type > &A) |
Change the matrix to be preconditioned. More... | |
Implementation of Tpetra::Operator | |
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 |
Apply the preconditioner to X, returning the result in Y. Y = alpha*Op(A)*X + beta*Y. More... | |
Teuchos::RCP< const map_type > | getDomainMap () const |
The Tpetra::Map representing the domain of this operator. More... | |
Teuchos::RCP< const map_type > | getRangeMap () const |
The Tpetra::Map representing the range of this operator. More... | |
bool | hasTransposeApply () const |
Whether it's possible to apply the transpose of this operator. More... | |
void | applyMat (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) const |
Compute Y = Op(A)*X, where Op(A) is either A, \(A^T\), or \(A^H\). More... | |
Implementation of Teuchos::Describable | |
std::string | description () const |
A simple one-line description of this object. More... | |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
Print the object with some verbosity level to a Teuchos::FancyOStream. More... | |
Public Member Functions inherited from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > | |
virtual | ~Preconditioner () |
Destructor. More... | |
Public Member Functions inherited from Ifpack2::Details::CanChangeMatrix< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > | |
virtual void | setMatrix (const Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > &A)=0 |
Set the new matrix. More... | |
virtual | ~CanChangeMatrix () |
Destructor. More... | |
Attribute accessor methods | |
Teuchos::RCP< const row_matrix_type > | A_ |
The matrix for which this is a preconditioner. More... | |
Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const |
The communicator over which the matrix is distributed. More... | |
Teuchos::RCP< const row_matrix_type > | getMatrix () const |
The matrix for which this is a preconditioner. More... | |
Teuchos::RCP< const Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > > | getCrsMatrix () const |
Attempt to return the matrix A as a Tpetra::CrsMatrix. More... | |
double | getComputeFlops () const |
The total number of floating-point operations taken by all calls to compute(). More... | |
double | getApplyFlops () const |
The total number of floating-point operations taken by all calls to apply(). More... | |
int | getNumInitialize () const |
The total number of successful calls to initialize(). More... | |
int | getNumCompute () const |
The total number of successful calls to compute(). More... | |
int | getNumApply () const |
The total number of successful calls to apply(). More... | |
double | getInitializeTime () const |
The total time spent in all calls to initialize(). More... | |
double | getComputeTime () const |
The total time spent in all calls to compute(). More... | |
double | getApplyTime () const |
The total time spent in all calls to apply(). More... | |
size_t | getNodeSmootherComplexity () const |
Get a rough estimate of cost per iteration. More... | |
Overlapping Schwarz where redundant patches are not stored explicitly.
MatrixType | The type of matrix to use. |
This class implements an overlapping Schwarz assuming it has already been provided overlapping data containers. Additionally, it is assumed there is at least one row with PatchSize_ nonzeros for each patch. This means that a matrix obtained from piecewise linear FEMs in 2d, for example, will not have rows with 4 nonzeros for each patch. However, piecewise quadratic FEMs in 2d will have one row with 9 nonzeros for each cell in the mesh, meaning all the patches will be detected. For a matrix corresponding to a Taylor-Hood discretization in 2d (quadratic velocities and linear pressures), the patch size would be 22.
The general algorithm proceeds as follows:
In general, there is a noticeable speedup when using this method compared to a typical method. This speedup may be further improved by using more advanced linear algebra interfaces such as batched Kokkos solves instead of the current LAPACK approach.
typedef MatrixType Ifpack2::DatabaseSchwarz< MatrixType >::matrix_type |
The template parameter of this class.
typedef MatrixType::scalar_type Ifpack2::DatabaseSchwarz< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::DatabaseSchwarz< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::DatabaseSchwarz< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type::device_type Ifpack2::DatabaseSchwarz< MatrixType >::device_type |
The Kokkos::Device specialization used by the input MatrixType.
typedef MatrixType::node_type Ifpack2::DatabaseSchwarz< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::DatabaseSchwarz< MatrixType >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::DatabaseSchwarz< MatrixType >::row_matrix_type |
The Tpetra::RowMatrix specialization matching MatrixType.
MatrixType must be a Tpetra::RowMatrix specialization. This typedef will always be a Tpetra::RowMatrix specialization.
typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::DatabaseSchwarz< MatrixType >::map_type |
The Tpetra::Map specialization matching MatrixType.
typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::DatabaseSchwarz< MatrixType >::vector_type |
The Tpetra::Vector specialization matching MatrixType.
If you wish to supply setParameters() a precomputed vector of diagonal entries of the matrix, use a pointer to an object of this type.
|
explicit |
Constructor.
[in] | A | The sparse matrix to which to apply DatabaseSchwarz iteration. The matrix A must be square, and its domain Map and range Map must be the same. The latter means that the vectors x and y in the sparse matrix-vector product y = A*x must both have the same distribution over process(es). |
Ifpack2::DatabaseSchwarz< MatrixType >::DatabaseSchwarz | ( | const Teuchos::RCP< const row_matrix_type > & | A, |
Teuchos::ParameterList & | params | ||
) |
Constructor.
[in] | A | The sparse matrix to which to apply DatabaseSchwarz iteration. The matrix A must be square, and its domain Map and range Map must be the same. The latter means that the vectors x and y in the sparse matrix-vector product y = A*x must both have the same distribution over process(es). |
[in] | params | The parameterlist containing settings for the object, such as the patch size to search for. |
|
virtual |
Destructor.
|
virtual |
Set (or reset) parameters.
|
virtual |
Set this preconditioner's parameters.
Reimplemented from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >.
|
virtual |
Initialize the preconditioner.
The compute() method will call initialize() automatically if it has not yet been called, so you do not normally need to call this. However, it is correct to call initialize() yourself, and compute() will not call it again if it already has been called.
|
inlinevirtual |
Whether the preconditioner has been successfully initialized (by calling initialize()).
|
virtual |
(Re)compute the left scaling, and (if applicable) estimate max and min eigenvalues of D_inv * A.
|
inlinevirtual |
Whether compute() has been called at least once.
|
virtual |
Change the matrix to be preconditioned.
|
virtual |
Apply the preconditioner to X, returning the result in Y. Y = alpha*Op(A)*X + beta*Y.
|
virtual |
The Tpetra::Map representing the domain of this operator.
|
virtual |
The Tpetra::Map representing the range of this operator.
bool Ifpack2::DatabaseSchwarz< MatrixType >::hasTransposeApply | ( | ) | const |
Whether it's possible to apply the transpose of this operator.
void Ifpack2::DatabaseSchwarz< MatrixType >::applyMat | ( | 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 |
||
) | const |
Compute Y = Op(A)*X, where Op(A) is either A, \(A^T\), or \(A^H\).
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::DatabaseSchwarz< MatrixType >::getComm | ( | ) | const |
The communicator over which the matrix is distributed.
|
virtual |
The matrix for which this is a preconditioner.
Teuchos::RCP< const Tpetra::CrsMatrix< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::DatabaseSchwarz< MatrixType >::getCrsMatrix | ( | ) | const |
Attempt to return the matrix A as a Tpetra::CrsMatrix.
This class does not require that A be a Tpetra::CrsMatrix. If it is NOT, this method will return Teuchos::null.
double Ifpack2::DatabaseSchwarz< MatrixType >::getComputeFlops | ( | ) | const |
The total number of floating-point operations taken by all calls to compute().
double Ifpack2::DatabaseSchwarz< MatrixType >::getApplyFlops | ( | ) | const |
The total number of floating-point operations taken by all calls to apply().
|
virtual |
The total number of successful calls to initialize().
|
virtual |
The total number of successful calls to compute().
|
virtual |
The total number of successful calls to apply().
|
virtual |
The total time spent in all calls to initialize().
|
virtual |
The total time spent in all calls to compute().
|
virtual |
The total time spent in all calls to apply().
size_t Ifpack2::DatabaseSchwarz< MatrixType >::getNodeSmootherComplexity | ( | ) | const |
Get a rough estimate of cost per iteration.
std::string Ifpack2::DatabaseSchwarz< MatrixType >::description | ( | ) | const |
A simple one-line description of this object.
void Ifpack2::DatabaseSchwarz< MatrixType >::describe | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
) | const |
Print the object with some verbosity level to a Teuchos::FancyOStream.
Teuchos::RCP<const row_matrix_type> Ifpack2::DatabaseSchwarz< MatrixType >::A_ |
The matrix for which this is a preconditioner.