Ifpack2 Templated Preconditioning Package
Version 1.0
|
Interface for all Ifpack2 preconditioners. More...
#include <Ifpack2_Preconditioner.hpp>
Public Types | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
Public Member Functions | |
virtual | ~Preconditioner () |
Destructor. More... | |
virtual void | setParameters (const Teuchos::ParameterList &List)=0 |
Set this preconditioner's parameters. More... | |
virtual void | setZeroStartingSolution (bool zeroStartingSolution) |
Set this preconditioner's parameters. More... | |
virtual void | initialize ()=0 |
Set up the graph structure of this preconditioner. More... | |
virtual bool | isInitialized () const =0 |
True if the preconditioner has been successfully initialized, else false. More... | |
virtual void | compute ()=0 |
Set up the numerical values in this preconditioner. More... | |
virtual bool | isComputed () const =0 |
True if the preconditioner has been successfully computed, else false. More... | |
virtual Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > | getMatrix () const =0 |
The input matrix given to the constructor. More... | |
virtual int | getNumInitialize () const =0 |
The number of calls to initialize(). More... | |
virtual int | getNumCompute () const =0 |
The number of calls to compute(). More... | |
virtual int | getNumApply () const =0 |
The number of calls to apply(). More... | |
virtual double | getInitializeTime () const =0 |
The time (in seconds) spent in initialize(). More... | |
virtual double | getComputeTime () const =0 |
The time (in seconds) spent in compute(). More... | |
virtual double | getApplyTime () const =0 |
The time (in seconds) spent in apply(). More... | |
Methods implementing Tpetra::Operator. | |
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > | getDomainMap () const =0 |
The domain Map of this operator. More... | |
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > | getRangeMap () const =0 |
The range Map of this operator. More... | |
virtual void | apply (const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const =0 |
Apply the preconditioner to X, putting the result in Y. More... | |
Interface for all Ifpack2 preconditioners.
Scalar | Type of the matrix's entries; same as the first template parameter of Tpetra::RowMatrix |
LocalOrdinal | Type of the matrix's local indices; same as the second template parameter of Tpetra::RowMatrix |
GlobalOrdinal | Type of the matrix's global indices; same as the third template parameter of Tpetra::RowMatrix |
Node | The matrix's Node type; same as the fourth template parameter of Tpetra::RowMatrix |
The Preconditioner class defines the interface that all Ifpack2 preconditioners must implement. Preconditioner inherits from Tpetra::Operator. Its apply() method applies the preconditioner.
If you are familiar with the IFPACK package, please be aware that this is different from IFPACK. In IFPACK, the ApplyInverse() method applies or "solves with" the preconditioner \(M^{-1}\), and the Apply() method "applies" the preconditioner \(M\). In Ifpack2, the apply() method applies or "solves with" the preconditioner \(M^{-1}\). Ifpack2 has no method comparable to IFPACK's Apply().
Preconditioner provides the following methods
Implementations of compute() must internally call initialize() if isInitialized() returns false. The preconditioner is applied by apply() (which returns if isComputed() is false). Every time that initialize() is called, the object destroys all the previously allocated information, and reinitializes the preconditioner. Every time compute() is called, the object recomputes the actual values of the preconditioner.
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
|
inlinevirtual |
Destructor.
|
pure virtual |
The domain Map of this operator.
The domain Map describes the distribution of valid input vectors X to the apply() method.
Implemented in Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::MDF< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The range Map of this operator.
The range Map describes the distribution of valid output vectors Y to the apply() method.
Implemented in Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::MDF< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Apply the preconditioner to X, putting the result in Y.
If the result of applying this preconditioner to a vector X is \(F \cdot X$\), then this method computes \(\beta Y + \alpha F \cdot X\). The typical case is \(\beta = 0\) and \(\alpha = 1\).
Implemented in Ifpack2::Chebyshev< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::MDF< MatrixType >, Ifpack2::ILUT< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
Set this preconditioner's parameters.
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Ifpack2::IdentitySolver< MatrixType >.
|
inlinevirtual |
Set this preconditioner's parameters.
Reimplemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::DatabaseSchwarz< MatrixType >, and Ifpack2::Hiptmair< MatrixType >.
|
pure virtual |
Set up the graph structure of this preconditioner.
If the graph structure of the constructor's input matrix has changed, or if you have not yet called initialize(), you must call initialize() before you may call compute() or apply().
Thus, initialize() corresponds to the "symbolic factorization" step of a sparse factorization, whether or not the specific preconditioner actually does a sparse factorization.
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::Experimental::RBILUK< MatrixType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Ifpack2::IdentitySolver< MatrixType >.
|
pure virtual |
True if the preconditioner has been successfully initialized, else false.
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Ifpack2::IdentitySolver< MatrixType >.
|
pure virtual |
Set up the numerical values in this preconditioner.
If the values of the constructor's input matrix have changed, or if you have not yet called compute(), you must call compute() before you may call apply().
Thus, compute() corresponds to the "numeric factorization" step of a sparse factorization, whether or not the specific preconditioner actually does a sparse factorization.
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::Experimental::RBILUK< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Ifpack2::IdentitySolver< MatrixType >.
|
pure virtual |
True if the preconditioner has been successfully computed, else false.
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Ifpack2::IdentitySolver< MatrixType >.
|
pure virtual |
The input matrix given to the constructor.
Implemented in Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::MDF< MatrixType >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The number of calls to initialize().
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The number of calls to compute().
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The number of calls to apply().
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The time (in seconds) spent in initialize().
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The time (in seconds) spent in compute().
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
|
pure virtual |
The time (in seconds) spent in apply().
Implemented in Ifpack2::AdditiveSchwarz< MatrixType, LocalInverseType >, Ifpack2::Chebyshev< MatrixType >, Ifpack2::Relaxation< MatrixType >, Ifpack2::Details::DenseSolver< MatrixType, false >, Ifpack2::Details::TriDiSolver< MatrixType, false >, Ifpack2::RILUK< MatrixType >, Ifpack2::RILUK< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >, Ifpack2::ILUT< MatrixType >, Ifpack2::BlockRelaxation< MatrixType, ContainerType >, Ifpack2::LocalSparseTriangularSolver< MatrixType >, Ifpack2::Details::Amesos2Wrapper< MatrixType >, Ifpack2::DatabaseSchwarz< MatrixType >, Ifpack2::MDF< MatrixType >, Ifpack2::Hiptmair< MatrixType >, Ifpack2::IdentitySolver< MatrixType >, and Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >.