Ifpack2 Templated Preconditioning Package
Version 1.0
|
Wrapper class for direct solvers in Amesos2. More...
#include <Ifpack2_Details_Amesos2Wrapper_decl.hpp>
Public Types | |
Typedefs | |
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 | 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 |
Type of the Tpetra::RowMatrix specialization that this class uses. More... | |
typedef Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > | map_type |
Type of the Tpetra::Map specialization that this class uses. More... | |
typedef Tpetra::CrsMatrix < scalar_type, local_ordinal_type, global_ordinal_type, node_type > | crs_matrix_type |
Type of the Tpetra::CrsMatrix specialization that this class uses. 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 | |
Constructors and destructor | |
Amesos2Wrapper (const Teuchos::RCP< const row_matrix_type > &A) | |
Constructor. More... | |
virtual | ~Amesos2Wrapper () |
Destructor. More... | |
Methods for setting up and computing the factorization | |
void | setParameters (const Teuchos::ParameterList ¶ms) |
Set parameters. More... | |
void | initialize () |
Compute the preordering and symbolic factorization of the matrix. More... | |
bool | isInitialized () const |
Returns true if the preconditioner has been successfully initialized. More... | |
void | compute () |
Compute the numeric factorization of the matrix. More... | |
bool | isComputed () const |
True if compute() completed successfully, else false. 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, resulting in Y. More... | |
Teuchos::RCP< const map_type > | getDomainMap () const |
Tpetra::Map representing the domain of this operator. More... | |
Teuchos::RCP< const map_type > | getRangeMap () const |
Tpetra::Map representing the range of this operator. More... | |
bool | hasTransposeApply () const |
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable). More... | |
Mathematical functions | |
Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const |
The input matrix's communicator. More... | |
Teuchos::RCP< const row_matrix_type > | getMatrix () const |
The input matrix; the matrix to be preconditioned. 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 in seconds spent in successful calls to initialize(). More... | |
double | getComputeTime () const |
The total time in seconds spent in successful calls to compute(). More... | |
double | getApplyTime () const |
The total time in seconds spent in successful calls to apply(). More... | |
Implementation of Teuchos::Describable | |
std::string | description () const |
A 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 the given output stream. 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... | |
virtual void | setZeroStartingSolution (bool zeroStartingSolution) |
Set this preconditioner's parameters. 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... | |
Wrapper class for direct solvers in Amesos2.
MatrixType | A specialization of Tpetra::RowMatrix. |
This class computes a sparse factorization of the input matrix A using Amesos2. The apply() method solves linear system(s) using that factorization. As with all Ifpack2 preconditioners, initialize() computes the symbolic factorization, and compute() computes the numeric factorization.
typedef MatrixType::scalar_type Ifpack2::Details::Amesos2Wrapper< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::Details::Amesos2Wrapper< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::Details::Amesos2Wrapper< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type Ifpack2::Details::Amesos2Wrapper< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::Details::Amesos2Wrapper< 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::Details::Amesos2Wrapper< MatrixType >::row_matrix_type |
Type of the Tpetra::RowMatrix specialization that this class uses.
typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::Details::Amesos2Wrapper< MatrixType >::map_type |
Type of the Tpetra::Map specialization that this class uses.
typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::Details::Amesos2Wrapper< MatrixType >::crs_matrix_type |
Type of the Tpetra::CrsMatrix specialization that this class uses.
|
explicit |
Constructor.
A | [in] The sparse matrix to factor, as a Tpetra::RowMatrix. (Tpetra::CrsMatrix inherits from this, so you may use a Tpetra::CrsMatrix here instead.) |
|
virtual |
Destructor.
|
virtual |
|
virtual |
Compute the preordering and symbolic factorization of the matrix.
You must call this method before you may call compute() or apply(), under the following conditions:
Please also see the documentation of compute() for the conditions under which you must call compute() before calling apply().
|
inlinevirtual |
Returns true
if the preconditioner has been successfully initialized.
|
virtual |
Compute the numeric factorization of the matrix.
You must call this method before you may call apply(), under the following conditions:
Please also see the documentation of initialize() for the conditions under which you must call initialize() before calling compute() or apply().
|
inlinevirtual |
True if compute() completed successfully, else false.
|
virtual |
Change the matrix to be preconditioned.
A | [in] The new matrix. |
! isInitialized ()
! isComputed ()
Calling this method resets the preconditioner's state. After calling this method with a nonnull input, you must first call initialize() and compute() (in that order) before you may call apply().
You may call this method with a null input. If A is null, then you may not call initialize() or compute() until you first call this method again with a nonnull input. This method invalidates any previous factorization whether or not A is null, so calling setMatrix() with a null input is one way to clear the preconditioner's state (and free any memory that it may be using).
The new matrix A need not necessarily have the same Maps or even the same communicator as the original matrix.
|
virtual |
Apply the preconditioner to X, resulting in Y.
If \(M^{-1}\) represents the action of the preconditioner, then this routine computes \(Y := beta \cdot Y + alpha \cdot Op(M^{-1}) X\), where \(Op(M^{-1})\) can be either \(M^{-1}\), \(M^{-T}\), or \(M^{-H}\), depending on mode
.
X | [in] Input multivector; "right-hand side" of the solve. |
Y | [out] Output multivector; result of the solve. |
mode | [in] Whether to solve with the original matrix, its transpose, or its conjugate transpose. |
alpha | [in] Scaling factor for the result of applying the preconditioner. |
alpha | [in] Scaling factor for Y. |
|
virtual |
Tpetra::Map representing the domain of this operator.
|
virtual |
Tpetra::Map representing the range of this operator.
bool Ifpack2::Details::Amesos2Wrapper< MatrixType >::hasTransposeApply | ( | ) | const |
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable).
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::Details::Amesos2Wrapper< MatrixType >::getComm | ( | ) | const |
The input matrix's communicator.
|
virtual |
The input matrix; the matrix to be preconditioned.
|
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 in seconds spent in successful calls to initialize().
|
virtual |
The total time in seconds spent in successful calls to compute().
|
virtual |
The total time in seconds spent in successful calls to apply().
std::string Ifpack2::Details::Amesos2Wrapper< MatrixType >::description | ( | ) | const |
A one-line description of this object.
void Ifpack2::Details::Amesos2Wrapper< MatrixType >::describe | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
) | const |
Print the object with some verbosity level to the given output stream.