Ifpack Package Browser (Single Doxygen Collection)
Development
|
Ifpack_BlockRelaxation: a class to define block relaxation preconditioners of Epetra_RowMatrix's. More...
#include <Ifpack_BlockRelaxation.h>
Public Member Functions | |
int | NumLocalBlocks () const |
Returns the number local blocks. More... | |
virtual bool | IsInitialized () const |
Returns true if the preconditioner has been successfully computed. More... | |
virtual bool | IsComputed () const |
Returns true if the preconditioner has been successfully computed. More... | |
virtual int | SetParameters (Teuchos::ParameterList &List) |
Sets all the parameters for the preconditioner. More... | |
virtual int | Initialize () |
Initializes the preconditioner. More... | |
virtual int | Compute () |
Computes the preconditioner. More... | |
virtual const Epetra_RowMatrix & | Matrix () const |
Returns a pointer to the matrix to be preconditioned. More... | |
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. More... | |
virtual double | Condest () const |
Returns the computed condition number estimate, or -1.0 if not computed. More... | |
std::ostream & | Print (std::ostream &os) const |
Prints basic information on iostream. This function is used by operator<<. More... | |
virtual int | NumInitialize () const |
Returns the number of calls to Initialize(). More... | |
virtual int | NumCompute () const |
Returns the number of calls to Compute(). More... | |
virtual int | NumApplyInverse () const |
Returns the number of calls to ApplyInverse(). More... | |
virtual double | InitializeTime () const |
Returns the time spent in Initialize(). More... | |
virtual double | ComputeTime () const |
Returns the time spent in Compute(). More... | |
virtual double | ApplyInverseTime () const |
Returns the time spent in ApplyInverse(). More... | |
virtual double | InitializeFlops () const |
Returns the number of flops in the initialization phase. More... | |
virtual double | ComputeFlops () const |
Returns the number of flops in the computation phase. More... | |
virtual double | ApplyInverseFlops () const |
Returns the number of flops in the application of the preconditioner. More... | |
Private Member Functions | |
Ifpack_BlockRelaxation (const Ifpack_BlockRelaxation &rhs) | |
Copy constructor (PRIVATE, should not be used). More... | |
Ifpack_BlockRelaxation & | operator= (const Ifpack_BlockRelaxation &rhs) |
operator= (PRIVATE, should not be used). More... | |
virtual int | ApplyInverseJacobi (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | DoJacobi (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | ApplyInverseGS (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | DoGaussSeidel (Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | ApplyInverseSGS (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | DoSGS (const Epetra_MultiVector &X, Epetra_MultiVector &Xtmp, Epetra_MultiVector &Y) const |
int | ExtractSubmatrices () |
Ifpack_BlockRelaxation (const Epetra_RowMatrix *Matrix) | |
Ifpack_BlockRelaxation constructor with given Epetra_RowMatrix. More... | |
virtual | ~Ifpack_BlockRelaxation () |
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) More... | |
virtual int | SetUseTranspose (bool UseTranspose_in) |
virtual const char * | Label () const |
virtual bool | UseTranspose () const |
Returns the current UseTranspose setting. More... | |
virtual bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. More... | |
virtual const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this operator. More... | |
virtual const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. More... | |
virtual const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. More... | |
bool | IsInitialized_ |
If true, the preconditioner has been successfully initialized. More... | |
bool | IsComputed_ |
If true, the preconditioner has been successfully computed. More... | |
int | NumInitialize_ |
Contains the number of successful calls to Initialize(). More... | |
int | NumCompute_ |
Contains the number of successful call to Compute(). More... | |
int | NumApplyInverse_ |
Contains the number of successful call to ApplyInverse(). More... | |
double | InitializeTime_ |
Contains the time for all successful calls to Initialize(). More... | |
double | ComputeTime_ |
Contains the time for all successful calls to Compute(). More... | |
double | ApplyInverseTime_ |
Contains the time for all successful calls to ApplyInverse(). More... | |
double | InitializeFlops_ |
Contains the number of flops for Initialize(). More... | |
double | ComputeFlops_ |
Contains the number of flops for Compute(). More... | |
double | ApplyInverseFlops_ |
Contain sthe number of flops for ApplyInverse(). More... | |
int | NumSweeps_ |
Number of preconditioning sweeps. More... | |
double | DampingFactor_ |
Damping parameter. More... | |
int | NumLocalBlocks_ |
Number of local blocks. More... | |
Teuchos::ParameterList | List_ |
Parameters list to be used to solve on each subblock. More... | |
Teuchos::RefCountPtr< const Epetra_RowMatrix > | Matrix_ |
Containers_[i] contains all the necessary information to solve on each subblock. More... | |
std::vector < Teuchos::RefCountPtr< T > > | Containers_ |
Epetra_Vector | Diagonal_ |
Teuchos::RefCountPtr < Ifpack_Partitioner > | Partitioner_ |
Contains information about non-overlapping partitions. More... | |
std::string | PartitionerType_ |
int | PrecType_ |
std::string | Label_ |
Label for this object. More... | |
bool | ZeroStartingSolution_ |
If true , starting solution is the zero vector. More... | |
Teuchos::RefCountPtr < Ifpack_Graph > | Graph_ |
Teuchos::RefCountPtr < Epetra_Vector > | W_ |
Weights for overlapping Jacobi only. More... | |
int | OverlapLevel_ |
Epetra_Time | Time_ |
bool | IsParallel_ |
Teuchos::RefCountPtr < Epetra_Import > | Importer_ |
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:
Then, we declare the preconditioner. Note that this is done through the class Ifpack_AdditiveSchwarz (see note below in this section).
The complete list of supported parameters is reported in page ifp_params. For a presentation of basic relaxation schemes, please refer to page Ifpack_PointRelaxation.
Definition at line 139 of file Ifpack_BlockRelaxation.h.
Ifpack_BlockRelaxation< T >::Ifpack_BlockRelaxation | ( | const Epetra_RowMatrix * | Matrix | ) |
Ifpack_BlockRelaxation constructor with given Epetra_RowMatrix.
Creates an Ifpack_Preconditioner preconditioner.
In | Matrix - Pointer to matrix to be preconditioned. |
Definition at line 455 of file Ifpack_BlockRelaxation.h.
|
virtual |
Definition at line 485 of file Ifpack_BlockRelaxation.h.
|
private |
Copy constructor (PRIVATE, should not be used).
|
virtual |
Applies the matrix to an Epetra_MultiVector.
In | X - A Epetra_MultiVector of dimension NumVectors to multiply with matrix. |
Out | Y -A Epetra_MultiVector of dimension NumVectors containing the result. |
Implements Epetra_Operator.
Definition at line 499 of file Ifpack_BlockRelaxation.h.
|
virtual |
Applies the block Jacobi preconditioner to X, returns the result in Y.
In | X - A Epetra_MultiVector of dimension NumVectors to be preconditioned. |
Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. |
Implements Ifpack_Preconditioner.
Definition at line 614 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the infinity norm of the global matrix (not implemented)
Implements Epetra_Operator.
Definition at line 184 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Implements Epetra_Operator.
Definition at line 192 of file Ifpack_BlockRelaxation.h.
|
virtual |
Implements Epetra_Operator.
Definition at line 491 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the current UseTranspose setting.
Implements Epetra_Operator.
Definition at line 202 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Implements Epetra_Operator.
Definition at line 208 of file Ifpack_BlockRelaxation.h.
|
virtual |
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Implements Epetra_Operator.
Definition at line 508 of file Ifpack_BlockRelaxation.h.
|
virtual |
Returns the Epetra_Map object associated with the domain of this operator.
Implements Epetra_Operator.
Definition at line 516 of file Ifpack_BlockRelaxation.h.
|
virtual |
Returns the Epetra_Map object associated with the range of this operator.
Implements Epetra_Operator.
Definition at line 524 of file Ifpack_BlockRelaxation.h.
|
inline |
Returns the number local blocks.
Definition at line 224 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns true
if the preconditioner has been successfully computed.
Implements Ifpack_Preconditioner.
Definition at line 230 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns true
if the preconditioner has been successfully computed.
Implements Ifpack_Preconditioner.
Definition at line 236 of file Ifpack_BlockRelaxation.h.
|
virtual |
Sets all the parameters for the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 1166 of file Ifpack_BlockRelaxation.h.
|
virtual |
Initializes the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 1238 of file Ifpack_BlockRelaxation.h.
|
virtual |
Computes the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 581 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns a pointer to the matrix to be preconditioned.
Implements Ifpack_Preconditioner.
Definition at line 250 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Computes the condition number estimate, returns its value.
Implements Ifpack_Preconditioner.
Definition at line 255 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the computed condition number estimate, or -1.0 if not computed.
Implements Ifpack_Preconditioner.
Definition at line 263 of file Ifpack_BlockRelaxation.h.
|
virtual |
Prints basic information on iostream. This function is used by operator<<.
Implements Ifpack_Preconditioner.
Definition at line 1108 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the number of calls to Initialize().
Implements Ifpack_Preconditioner.
Definition at line 271 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the number of calls to Compute().
Implements Ifpack_Preconditioner.
Definition at line 277 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the number of calls to ApplyInverse().
Implements Ifpack_Preconditioner.
Definition at line 283 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the time spent in Initialize().
Implements Ifpack_Preconditioner.
Definition at line 289 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the time spent in Compute().
Implements Ifpack_Preconditioner.
Definition at line 295 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the time spent in ApplyInverse().
Implements Ifpack_Preconditioner.
Definition at line 301 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the number of flops in the initialization phase.
Implements Ifpack_Preconditioner.
Definition at line 307 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the number of flops in the computation phase.
Implements Ifpack_Preconditioner.
Definition at line 325 of file Ifpack_BlockRelaxation.h.
|
inlinevirtual |
Returns the number of flops in the application of the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 340 of file Ifpack_BlockRelaxation.h.
|
inlineprivate |
operator= (PRIVATE, should not be used).
Definition at line 362 of file Ifpack_BlockRelaxation.h.
|
privatevirtual |
Definition at line 656 of file Ifpack_BlockRelaxation.h.
|
privatevirtual |
Definition at line 687 of file Ifpack_BlockRelaxation.h.
|
privatevirtual |
Definition at line 795 of file Ifpack_BlockRelaxation.h.
|
privatevirtual |
Definition at line 816 of file Ifpack_BlockRelaxation.h.
|
privatevirtual |
Definition at line 924 of file Ifpack_BlockRelaxation.h.
|
privatevirtual |
Definition at line 942 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 531 of file Ifpack_BlockRelaxation.h.
|
private |
If true, the preconditioner has been successfully initialized.
Definition at line 391 of file Ifpack_BlockRelaxation.h.
|
private |
If true, the preconditioner has been successfully computed.
Definition at line 393 of file Ifpack_BlockRelaxation.h.
|
private |
Contains the number of successful calls to Initialize().
Definition at line 395 of file Ifpack_BlockRelaxation.h.
|
private |
Contains the number of successful call to Compute().
Definition at line 397 of file Ifpack_BlockRelaxation.h.
|
mutableprivate |
Contains the number of successful call to ApplyInverse().
Definition at line 399 of file Ifpack_BlockRelaxation.h.
|
private |
Contains the time for all successful calls to Initialize().
Definition at line 401 of file Ifpack_BlockRelaxation.h.
|
private |
Contains the time for all successful calls to Compute().
Definition at line 403 of file Ifpack_BlockRelaxation.h.
|
mutableprivate |
Contains the time for all successful calls to ApplyInverse().
Definition at line 405 of file Ifpack_BlockRelaxation.h.
|
private |
Contains the number of flops for Initialize().
Definition at line 407 of file Ifpack_BlockRelaxation.h.
|
private |
Contains the number of flops for Compute().
Definition at line 409 of file Ifpack_BlockRelaxation.h.
|
mutableprivate |
Contain sthe number of flops for ApplyInverse().
Definition at line 411 of file Ifpack_BlockRelaxation.h.
|
private |
Number of preconditioning sweeps.
Definition at line 416 of file Ifpack_BlockRelaxation.h.
|
private |
Damping parameter.
Definition at line 418 of file Ifpack_BlockRelaxation.h.
|
private |
Number of local blocks.
Definition at line 420 of file Ifpack_BlockRelaxation.h.
|
private |
Parameters list to be used to solve on each subblock.
Definition at line 422 of file Ifpack_BlockRelaxation.h.
|
private |
Containers_[i] contains all the necessary information to solve on each subblock.
Pointers to the matrix to be preconditioned.
Definition at line 428 of file Ifpack_BlockRelaxation.h.
|
mutableprivate |
Definition at line 429 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 430 of file Ifpack_BlockRelaxation.h.
|
private |
Contains information about non-overlapping partitions.
Definition at line 433 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 434 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 435 of file Ifpack_BlockRelaxation.h.
|
private |
Label for this
object.
Definition at line 437 of file Ifpack_BlockRelaxation.h.
|
private |
If true
, starting solution is the zero vector.
Definition at line 439 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 440 of file Ifpack_BlockRelaxation.h.
|
private |
Weights for overlapping Jacobi only.
Definition at line 442 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 444 of file Ifpack_BlockRelaxation.h.
|
mutableprivate |
Definition at line 445 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 446 of file Ifpack_BlockRelaxation.h.
|
private |
Definition at line 447 of file Ifpack_BlockRelaxation.h.