Ifpack_PointRelaxation: a class to define point relaxation preconditioners of for Epetra_RowMatrix's. More...
#include <Ifpack_PointRelaxation.h>
Public Member Functions | |
virtual int | SetUseTranspose (bool UseTranspose_in) |
Ifpack_PointRelaxation (const Epetra_RowMatrix *Matrix) | |
Ifpack_PointRelaxation constructor with given Epetra_RowMatrix. More... | |
virtual | ~Ifpack_PointRelaxation () |
Destructor. | |
virtual int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Applies the matrix to an Epetra_MultiVector. More... | |
virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Applies the preconditioner to X, returns the result in Y. More... | |
virtual double | NormInf () const |
Returns the infinity norm of the global matrix (not implemented) | |
virtual const char * | Label () const |
virtual bool | UseTranspose () const |
Returns the current UseTranspose setting. | |
virtual bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. | |
virtual const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this operator. | |
virtual const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. | |
virtual const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. | |
virtual int | Initialize () |
Computes all it is necessary to initialize the preconditioner. | |
virtual bool | IsInitialized () const |
Returns true if the preconditioner has been successfully initialized, false otherwise. | |
virtual bool | IsComputed () const |
Returns true if the preconditioner has been successfully computed. | |
virtual int | Compute () |
Computes the preconditioners. | |
virtual const Epetra_RowMatrix & | Matrix () const |
Returns a pointer to the matrix to be preconditioned. | |
virtual double | Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix=0) |
Computes the condition number estimates and returns the value. | |
virtual double | Condest () const |
Returns the condition number estimate, or -1.0 if not computed. | |
virtual int | SetParameters (Teuchos::ParameterList &List) |
Sets all the parameters for the preconditioner. | |
virtual std::ostream & | Print (std::ostream &os) const |
Prints object to an output stream. | |
virtual int | NumInitialize () const |
Returns the number of calls to Initialize(). | |
virtual int | NumCompute () const |
Returns the number of calls to Compute(). | |
virtual int | NumApplyInverse () const |
Returns the number of calls to ApplyInverse(). | |
virtual double | InitializeTime () const |
Returns the time spent in Initialize(). | |
virtual double | ComputeTime () const |
Returns the time spent in Compute(). | |
virtual double | ApplyInverseTime () const |
Returns the time spent in ApplyInverse(). | |
virtual double | InitializeFlops () const |
Returns the number of flops in the initialization phase. | |
virtual double | ComputeFlops () const |
Returns the number of flops in the computation phase. | |
virtual double | ApplyInverseFlops () const |
Returns the number of flops for the application of the preconditioner. | |
Ifpack_PointRelaxation: a class to define point relaxation preconditioners of for Epetra_RowMatrix's.
The Ifpack_PointRelaxation class enables the construction of point relaxation preconditioners of an Epetra_RowMatrix. Ifpack_PointRelaxation is derived from the Ifpack_Preconditioner class, which is itself derived from Epetra_Operator. Therefore this object can be used as preconditioner everywhere an ApplyInverse() method is required in the preconditioning step.
This class enables the construction of the following simple preconditioners:
We now briefly describe the main features of the above preconditioners. Consider a linear system of type
\[ A x = b, \]
where \( A\) is a square, real matrix, and \( x, b\) are two real vectors. We begin with the decomposition
\[ A = D - E - F \]
where \( D\) is the diagonal of A, \( -E\) is the strict lower part, and \( -F\) is the strict upper part. It is assumed that the diagonal entries of \( A\) are different from zero.
Given an starting solution \( x_0\), an iteration of the (damped) Jacobi method can be written in matrix form as follows:
\[ x_{k+1} = \omega D^{-1}(E + F) x_k + D_{-1}b, \]
for \( k < k_{max}\), and \(\omega \) a damping parameter.
Using Ifpack_Jacobi, the user can apply the specified number of sweeps ( \( k_{max}\)), and the damping parameter. If only one sweep is used, then the class simply applies the inverse of the diagonal of A to the input vector.
Given an starting solution \( x_0\), an iteration of the (damped) GaussSeidel method can be written in matrix form as follows:
\[ (D - E) x_{k+1} = \omega F x_k + b, \]
for \( k < k_{max}\), and \(\omega \) a damping parameter. Equivalently, the Gauss-Seidel preconditioner can be defined as
\[ P_{GS}^{-1} = (D - E)^{-1}. \]
Clearly, the role of E and F can be interchanged. However, Ifpack_GaussSeidel does not consider backward Gauss-Seidel methods.
For a list of supported parameters, please refer to page List of Supported Parameters.
The complete list of supported parameters is reported in page List of Supported Parameters. For a presentation of basic relaxation schemes, please refer to page Ifpack_PointRelaxation.
Definition at line 130 of file Ifpack_PointRelaxation.h.
Ifpack_PointRelaxation::Ifpack_PointRelaxation | ( | const Epetra_RowMatrix * | Matrix | ) |
Ifpack_PointRelaxation constructor with given Epetra_RowMatrix.
Creates an instance of Ifpack_PointRelaxation class.
Matrix | - (In) Pointer to matrix to precondition. |
Definition at line 63 of file Ifpack_PointRelaxation.cpp.
|
inlinevirtual |
Applies the matrix to an Epetra_MultiVector.
X | - (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix. |
Y | - (Out) A Epetra_MultiVector of dimension NumVectors containing the result. |
Implements Epetra_Operator.
Definition at line 172 of file Ifpack_PointRelaxation.h.
References IsComputed(), and UseTranspose().
|
virtual |
Applies the preconditioner to X, returns the result in Y.
X | - (In) A Epetra_MultiVector of dimension NumVectors to be preconditioned. |
Y | - (InOut) A Epetra_MultiVector of dimension NumVectors containing result. |
Implements Ifpack_Preconditioner.
Definition at line 414 of file Ifpack_PointRelaxation.cpp.
References IsComputed().
|
inlinevirtual |
This flag can be used to apply the preconditioner to the transpose of the input operator.
Implements Epetra_Operator.
Definition at line 153 of file Ifpack_PointRelaxation.h.