IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
List of all members
Ifpack_ILU Class Reference

Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a given Epetra_RowMatrix. More...

#include <Ifpack_ILU.h>

Inheritance diagram for Ifpack_ILU:
Inheritance graph
[legend]
Collaboration diagram for Ifpack_ILU:
Collaboration graph
[legend]

Public Member Functions

 Ifpack_ILU (Epetra_RowMatrix *A)
 Constructor.
 
 ~Ifpack_ILU ()
 Destructor.
 
int Initialize ()
 Initialize the preconditioner, does not touch matrix values.
 
bool IsInitialized () const
 Returns true if the preconditioner has been successfully initialized.
 
int Compute ()
 Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxation parameters. More...
 
bool IsComputed () const
 If factor is completed, this query returns true, otherwise it returns false.
 
int SetParameters (Teuchos::ParameterList &parameterlist)
 Set parameters using a Teuchos::ParameterList object.
 
int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied. More...
 
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
int Multiply (bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More...
 
double Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix_in=0)
 Computes the estimated condition number and returns the value.
 
double Condest () const
 Returns the computed estimated condition number, or -1.0 if not computed.
 
const Epetra_CrsMatrixL () const
 Returns the address of the L factor associated with this factored matrix.
 
const Epetra_VectorD () const
 Returns the address of the D factor associated with this factored matrix.
 
const Epetra_CrsMatrixU () const
 Returns the address of the L factor associated with this factored matrix.
 
const char * Label () const
 Returns a character string describing the operator.
 
int SetLabel (const char *Label_in)
 Sets label for this object.
 
double NormInf () const
 Returns 0.0 because this class cannot compute Inf-norm.
 
bool HasNormInf () const
 Returns false because this class cannot compute an Inf-norm.
 
bool UseTranspose () const
 Returns the current UseTranspose setting.
 
const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator.
 
const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator.
 
const Epetra_CommComm () const
 Returns the Epetra_BlockMap object associated with the range of this matrix operator.
 
const Epetra_RowMatrixMatrix () const
 Returns a reference to the matrix to be preconditioned.
 
virtual std::ostream & Print (std::ostream &os) const
 Prints on stream basic information about this object.
 
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.
 

Detailed Description

Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a given Epetra_RowMatrix.

The Ifpack_ILU class computes a "relaxed" ILU factorization with level k fill of a given Epetra_RowMatrix. The notion of relaxation is the same as described in "Experimental study of ILU preconditioners for indefinite matrices" by Chow and Saad.

Please refer to General description of incomplete factorizations for a general description of the ILU algorithm.

The complete list of supported parameters is reported in page List of Supported Parameters.

Author
Mike Heroux, Marzio Sala, SNL 9214.
Date
Last modified on 22-Jan-05.

Definition at line 83 of file Ifpack_ILU.h.

Member Function Documentation

int Ifpack_ILU::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.

In this implementation, we use several existing attributes to determine how virtual method ApplyInverse() should call the concrete method Solve(). We pass in the UpperTriangular(), the Epetra_CrsMatrix::UseTranspose(), and NoDiagonal() methods. The most notable warning is that if a matrix has no diagonal values we assume that there is an implicit unit diagonal that should be accounted for when doing a triangular solve.

Parameters
X- (In) A Epetra_MultiVector of dimension NumVectors to solve for.
OutY - (Out) A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implements Ifpack_Preconditioner.

Definition at line 575 of file Ifpack_ILU.cpp.

References Epetra_Time::ElapsedTime(), IsComputed(), Epetra_Time::ResetStartTime(), and UseTranspose().

int Ifpack_ILU::Compute ( )
virtual

Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxation parameters.

This function computes the ILU(k) factors L and U using the current:

  1. Ifpack_IlukGraph specifying the structure of L and U.
  2. Value for the ILU(k) relaxation parameter.
  3. Value for the a priori diagonal threshold values.

InitValues() must be called before the factorization can proceed.

Implements Ifpack_Preconditioner.

Definition at line 334 of file Ifpack_ILU.cpp.

References Comm(), Epetra_Time::ElapsedTime(), Initialize(), IsInitialized(), and Epetra_Time::ResetStartTime().

int Ifpack_ILU::SetUseTranspose ( bool  UseTranspose_in)
inlinevirtual

If set true, transpose of this operator will be applied.

This flag allows the transpose of the given operator to be used implicitly. Setting this flag affects only the Apply() and ApplyInverse() methods. If the implementation of this interface does not support transpose use, this method should return a value of -1.

Parameters
UseTranspose_in- (In) If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Always returns 0.

Implements Epetra_Operator.

Definition at line 145 of file Ifpack_ILU.h.


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