IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Public Member Functions | List of all members
Ifpack_BlockRelaxation< T > Class Template Reference

Ifpack_BlockRelaxation: a class to define block relaxation preconditioners of Epetra_RowMatrix's. More...

#include <Ifpack_BlockRelaxation.h>

Inheritance diagram for Ifpack_BlockRelaxation< T >:
Inheritance graph
[legend]
Collaboration diagram for Ifpack_BlockRelaxation< T >:
Collaboration graph
[legend]

Public Member Functions

int NumLocalBlocks () const
 Returns the number local blocks.
 
virtual bool IsInitialized () const
 Returns true if the preconditioner has been successfully computed.
 
virtual bool IsComputed () const
 Returns true if the preconditioner has been successfully computed.
 
virtual int SetParameters (Teuchos::ParameterList &List)
 Sets all the parameters for the preconditioner.
 
virtual int Initialize ()
 Initializes the preconditioner.
 
virtual int Compute ()
 Computes the preconditioner.
 
virtual const Epetra_RowMatrixMatrix () const
 Returns a pointer to the matrix to be preconditioned.
 
virtual double Condest (const Ifpack_CondestType=Ifpack_Cheap, const int=1550, const double=1e-9, Epetra_RowMatrix *=0)
 Computes the condition number estimate, returns its value.
 
virtual double Condest () const
 Returns the computed condition number estimate, or -1.0 if not computed.
 
std::ostream & Print (std::ostream &os) const
 Prints basic information on iostream. This function is used by operator<<.
 
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 in the application of the preconditioner.
 
 Ifpack_BlockRelaxation (const Epetra_RowMatrix *Matrix)
 Ifpack_BlockRelaxation constructor with given Epetra_RowMatrix. More...
 
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 block Jacobi preconditioner to X, returns the result in Y. More...
 
virtual double NormInf () const
 Returns the infinity norm of the global matrix (not implemented)
 
virtual int SetUseTranspose (bool UseTranspose_in)
 
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_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator.
 
virtual const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator.
 
virtual const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator.
 

Detailed Description

template<typename T>
class Ifpack_BlockRelaxation< T >

Ifpack_BlockRelaxation: a class to define block relaxation preconditioners of Epetra_RowMatrix's.

The Ifpack_BlockRelaxation class enables the construction of

block relaxation preconditioners of an Epetra_RowMatrix. Ifpack_PointRelaxation is derived from the Ifpack_Preconditioner class, which is derived from Epetra_Operator. Therefore this object can be used as preconditioner everywhere an ApplyInverse() method is required in the preconditioning step.

The class currently support:

The idea of block relaxation method is to extend their point relaxation counterpart (implemented in Ifpack_PointRelaxation), by working on a group of equation simulteneously. Generally, larger blocks result in better convergence and increased robusteness.

The user can decide:

The following is an example of usage of this preconditioner with dense containers. First, we include the header files:

#include "Ifpack_AdditiveSchwarz.h"
#include "Ifpack_BlockPreconditioner.h"
#include "Ifpack_DenseContainer.h"

Then, we declare the preconditioner. Note that this is done through the class Ifpack_AdditiveSchwarz (see note below in this section).

// A is an Epetra_RowMatrix
// List is a Teuchos::ParameterList
IFPACK_CHK_ERR(Prec.SetParameters(List));
IFPACK_CHK_ERR(Prec.Initialize());
IFPACK_CHK_ERR(Prec.Compute());
// action of the preconditioner is given by ApplyInverse()
// Now use it in AztecOO, solver is an AztecOO object
solver.SetPrecOperator(&Prec);

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.

Author
Marzio Sala, SNL 9214.
Date
Last modified on 25-Jan-05.

Definition at line 139 of file Ifpack_BlockRelaxation.h.

Constructor & Destructor Documentation

template<typename T >
Ifpack_BlockRelaxation< T >::Ifpack_BlockRelaxation ( const Epetra_RowMatrix Matrix)

Ifpack_BlockRelaxation constructor with given Epetra_RowMatrix.

Creates an Ifpack_Preconditioner preconditioner.

Parameters
InMatrix - Pointer to matrix to be preconditioned.

Definition at line 455 of file Ifpack_BlockRelaxation.h.

References Epetra_Operator::Comm(), and Epetra_Comm::NumProc().

Member Function Documentation

template<typename T >
int Ifpack_BlockRelaxation< T >::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Applies the matrix to an Epetra_MultiVector.

Parameters
InX - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
OutY -A Epetra_MultiVector of dimension NumVectors containing the result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 499 of file Ifpack_BlockRelaxation.h.

References Epetra_Operator::Apply().

template<typename T >
int Ifpack_BlockRelaxation< T >::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Applies the block Jacobi preconditioner to X, returns the result in Y.

Parameters
InX - A Epetra_MultiVector of dimension NumVectors to be preconditioned.
OutY -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implements Ifpack_Preconditioner.

Definition at line 614 of file Ifpack_BlockRelaxation.h.


The documentation for this class was generated from the following file: