Ifpack_CrsIct Class Reference

Ifpack_CrsIct: A class for constructing and using an incomplete Cholesky factorization of a given Epetra_CrsMatrix. More...

#include <Ifpack_CrsIct.h>

Inheritance diagram for Ifpack_CrsIct:

[legend]
Collaboration diagram for Ifpack_CrsIct:
[legend]
List of all members.

Public Member Functions

 Ifpack_CrsIct (const Epetra_CrsMatrix &A, double Droptol=1.0E-4, int Lfil=20)
 Ifpack_CrsIct constuctor with variable number of indices per row.
 Ifpack_CrsIct (const Ifpack_CrsIct &IctOperator)
 Copy constructor.
virtual ~Ifpack_CrsIct ()
 Ifpack_CrsIct Destructor.
void SetAbsoluteThreshold (double Athresh)
 Set absolute threshold value.
void SetRelativeThreshold (double Rthresh)
 Set relative threshold value.
void SetOverlapMode (Epetra_CombineMode OverlapMode)
 Set overlap mode type.
int SetParameters (const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
 Set parameters using a Teuchos::ParameterList object.
int InitValues (const Epetra_CrsMatrix &A)
 Initialize L and U with values from user matrix A.
bool ValuesInitialized () const
 If values have been initialized, this query returns true, otherwise it returns false.
int Factor ()
 Compute IC factor U using the specified graph, diagonal perturbation thresholds and relaxation parameters.
bool Factored () const
 If factor is completed, this query returns true, otherwise it returns false.
int Solve (bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Ifpack_CrsIct forward/back solve on a Epetra_MultiVector X in Y.
int Multiply (bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of multiplying U, D and U^T in that order on an Epetra_MultiVector X in Y.
int Condest (bool Trans, double &ConditionNumberEstimate) const
 Returns the maximum over all the condition number estimate for each local ILU set of factors.
double GetAbsoluteThreshold ()
 Get absolute threshold value.
double GetRelativeThreshold ()
 Get relative threshold value.
Epetra_CombineMode GetOverlapMode ()
 Get overlap mode type.
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global graph.
int NumMyNonzeros () const
 Returns the number of nonzero entries in the local graph.
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 U factor associated with this factored matrix.
Additional methods required to support the Epetra_Operator interface.
const char * Label () const
 Returns a character string describing the operator.
int SetUseTranspose (bool UseTranspose)
 If set true, transpose of this operator will be applied.
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
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.
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.

Protected Member Functions

void SetFactored (bool Flag)
void SetValuesInitialized (bool Flag)
bool Allocated () const
int SetAllocated (bool Flag)

Friends

ostream & operator<< (ostream &os, const Ifpack_CrsIct &A)
 << operator will work for Ifpack_CrsIct.

Detailed Description

Ifpack_CrsIct: A class for constructing and using an incomplete Cholesky factorization of a given Epetra_CrsMatrix.

The Ifpack_CrsIct class computes a threshold based incomplete LDL^T factorization of a given Epetra_CrsMatrix. The factorization that is produced is a function of several parameters:

  1. Maximum number of entries per row/column in factor - The factorization will contain at most this number of nonzero terms in each row/column of the factorization.

  2. Diagonal perturbation - Prior to computing the factorization, it is possible to modify the diagonal entries of the matrix for which the factorization will be computing. If the absolute and relative perturbation values are zero and one, respectively, the factorization will be compute for the original user matrix A. Otherwise, the factorization will computed for a matrix that differs from the original user matrix in the diagonal values only. Below we discuss the details of diagonal perturbations. The absolute and relative threshold values are set by calling SetAbsoluteThreshold() and SetRelativeThreshold(), respectively.

Estimating Preconditioner Condition Numbers

For ill-conditioned matrices, we often have difficulty computing usable incomplete factorizations. The most common source of problems is that the factorization may encounter a small or zero pivot, in which case the factorization can fail, or even if the factorization succeeds, the factors may be so poorly conditioned that use of them in the iterative phase produces meaningless results. Before we can fix this problem, we must be able to detect it. To this end, we use a simple but effective condition number estimate for $(LU)^{-1}$.

The condition of a matrix $B$, called $cond_p(B)$, is defined as $cond_p(B) = \|B\|_p\|B^{-1}\|_p$ in some appropriate norm $p$. $cond_p(B)$ gives some indication of how many accurate floating point digits can be expected from operations involving the matrix and its inverse. A condition number approaching the accuracy of a given floating point number system, about 15 decimal digits in IEEE double precision, means that any results involving $B$ or $B^{-1}$ may be meaningless.

The $\infty$-norm of a vector $y$ is defined as the maximum of the absolute values of the vector entries, and the $\infty$-norm of a matrix C is defined as $\|C\|_\infty = \max_{\|y\|_\infty = 1} \|Cy\|_\infty$. A crude lower bound for the $cond_\infty(C)$ is $\|C^{-1}e\|_\infty$ where $e = (1, 1, \ldots, 1)^T$. It is a lower bound because $cond_\infty(C) = \|C\|_\infty\|C^{-1}\|_\infty \ge \|C^{-1}\|_\infty \ge |C^{-1}e\|_\infty$.

For our purposes, we want to estimate $cond_\infty(LU)$, where $L$ and $U$ are our incomplete factors. Edmond in his Ph.D. thesis demonstrates that $\|(LU)^{-1}e\|_\infty$ provides an effective estimate for $cond_\infty(LU)$. Furthermore, since finding $z$ such that $LUz = y$ is a basic kernel for applying the preconditioner, computing this estimate of $cond_\infty(LU)$ is performed by setting $y = e$, calling the solve kernel to compute $z$ and then computing $\|z\|_\infty$.

A priori Diagonal Perturbations

Given the above method to estimate the conditioning of the incomplete factors, if we detect that our factorization is too ill-conditioned we can improve the conditioning by perturbing the matrix diagonal and restarting the factorization using this more diagonally dominant matrix. In order to apply perturbation, prior to starting the factorization, we compute a diagonal perturbation of our matrix $A$ and perform the factorization on this perturbed matrix. The overhead cost of perturbing the diagonal is minimal since the first step in computing the incomplete factors is to copy the matrix $A$ into the memory space for the incomplete factors. We simply compute the perturbed diagonal at this point.

The actual perturbation values we use are the diagonal values $(d_1, d_2, \ldots, d_n)$ with $d_i = sgn(d_i)\alpha + d_i\rho$, $i=1, 2, \ldots, n$, where $n$ is the matrix dimension and $sgn(d_i)$ returns the sign of the diagonal entry. This has the effect of forcing the diagonal values to have minimal magnitude of $\alpha$ and to increase each by an amount proportional to $\rho$, and still keep the sign of the original diagonal entry.

Constructing Ifpack_CrsIct objects

Constructing Ifpack_CrsIct objects is a multi-step process. The basic steps are as follows:

  1. Create Ifpack_CrsIct instance, including storage, via constructor.
  2. Enter values via one or more Put or SumInto functions.
  3. Complete construction via FillComplete call.

Note that, even after a matrix is constructed, it is possible to update existing matrix entries. It is not possible to create new entries.

Counting Floating Point Operations

Each Ifpack_CrsIct object keep track of the number of serial floating point operations performed using the specified object as the this argument to the function. The Flops() function returns this number as a double precision number. Using this information, in conjunction with the Epetra_Time class, one can get accurate parallel performance numbers. The ResetFlops() function resets the floating point counter.

Warning:
A Epetra_Map is required for the Ifpack_CrsIct constructor.

Definition at line 157 of file Ifpack_CrsIct.h.


Constructor & Destructor Documentation

Ifpack_CrsIct::Ifpack_CrsIct const Epetra_CrsMatrix A,
double  Droptol = 1.0E-4,
int  Lfil = 20
 

Ifpack_CrsIct constuctor with variable number of indices per row.

Creates a Ifpack_CrsIct object and allocates storage.

Parameters:
In A - User matrix to be factored.
In Graph - Graph generated by Ifpack_IlukGraph.

Definition at line 46 of file Ifpack_CrsIct.cpp.

References Zero.


Member Function Documentation

int Ifpack_CrsIct::Apply const Epetra_MultiVector X,
Epetra_MultiVector Y
const [inline, virtual]
 

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

Note that this implementation of Apply does NOT perform a forward back solve with the LDU factorization. Instead it applies these operators via multiplication with U, D and L respectively. The ApplyInverse() method performs a solve.

Parameters:
In X - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Out Y -A Epetra_MultiVector of dimension NumVectors containing result.
Returns:
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 319 of file Ifpack_CrsIct.h.

References Multiply(), and UseTranspose().

int Ifpack_CrsIct::ApplyInverse const Epetra_MultiVector X,
Epetra_MultiVector Y
const [inline, 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:
In X - A Epetra_MultiVector of dimension NumVectors to solve for.
Out Y -A Epetra_MultiVector of dimension NumVectors containing result.
Returns:
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 336 of file Ifpack_CrsIct.h.

References Solve(), and UseTranspose().

int Ifpack_CrsIct::Condest bool  Trans,
double &  ConditionNumberEstimate
const
 

Returns the maximum over all the condition number estimate for each local ILU set of factors.

This functions computes a local condition number estimate on each processor and return the maximum over all processor of the estimate.

Parameters:
In Trans -If true, solve transpose problem.
Out ConditionNumberEstimate - The maximum across all processors of the infinity-norm estimate of the condition number of the inverse of LDU.

Definition at line 387 of file Ifpack_CrsIct.cpp.

References Epetra_MultiVector::Abs(), Epetra_MultiVector::MaxValue(), Epetra_MultiVector::PutScalar(), Epetra_CrsMatrix::RowMap(), and Solve().

int Ifpack_CrsIct::Factor  ) 
 

Compute IC factor U using the specified graph, diagonal perturbation thresholds and relaxation parameters.

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

  1. Ifpack_IlukGraph specifying the structure of L and U.
  2. Value for the RILU(k) relaxation parameter.
  3. Value for the a priori diagonal threshold values.
InitValues() must be called before the factorization can proceed.

Definition at line 253 of file Ifpack_CrsIct.cpp.

References Epetra_CrsMatrix::Comm(), Epetra_Vector::ExtractView(), Factored(), Epetra_CrsMatrix::FillComplete(), Epetra_MultiVector::GlobalLength(), Epetra_CrsMatrix::InsertMyValues(), Epetra_CrsMatrix::NumGlobalNonzeros(), Epetra_CrsMatrix::OperatorDomainMap(), Epetra_CrsMatrix::OperatorRangeMap(), Epetra_MultiVector::Reciprocal(), Epetra_CrsMatrix::RowMatrixRowMap(), Epetra_Comm::SumAll(), Epetra_CompObject::UpdateFlops(), and View.

int Ifpack_CrsIct::InitValues const Epetra_CrsMatrix A  ) 
 

Initialize L and U with values from user matrix A.

Copies values from the user's matrix into the nonzero pattern of L and U.

Parameters:
In A - User matrix to be factored.
Warning:
The graph of A must be identical to the graph passed in to Ifpack_IlukGraph constructor.

Definition at line 166 of file Ifpack_CrsIct.cpp.

References Epetra_CrsMatrix::Comm(), Epetra_DistObject::DistributedGlobal(), Epetra_CrsMatrix::ExtractMyRowCopy(), Epetra_Vector::ExtractView(), Epetra_CrsMatrix::FillComplete(), Epetra_CrsMatrix::InsertMyValues(), Epetra_Comm::MaxAll(), Epetra_CrsMatrix::MaxNumEntries(), Epetra_CrsMatrix::NumMyRows(), Epetra_CrsMatrix::OperatorDomainMap(), Epetra_CrsMatrix::OperatorRangeMap(), and U().

int Ifpack_CrsIct::Multiply bool  Trans,
const Epetra_MultiVector X,
Epetra_MultiVector Y
const
 

Returns the result of multiplying U, D and U^T in that order on an Epetra_MultiVector X in Y.

Parameters:
In Trans -If true, multiply by L^T, D and U^T in that order.
In X - A Epetra_MultiVector of dimension NumVectors to solve for.
Out Y -A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns:
Integer error code, set to 0 if successful.

Definition at line 363 of file Ifpack_CrsIct.cpp.

References Epetra_CrsMatrix::Multiply(), Epetra_MultiVector::NumVectors(), Epetra_MultiVector::ReciprocalMultiply(), and Epetra_MultiVector::Update().

Referenced by Apply().

int Ifpack_CrsIct::SetUseTranspose bool  UseTranspose  )  [inline, virtual]
 

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:
In UseTranspose -If true, multiply by the transpose of operator, otherwise just use operator.
Returns:
Always returns 0.

Implements Epetra_Operator.

Definition at line 305 of file Ifpack_CrsIct.h.

int Ifpack_CrsIct::Solve bool  Trans,
const Epetra_MultiVector X,
Epetra_MultiVector Y
const
 

Returns the result of a Ifpack_CrsIct forward/back solve on a Epetra_MultiVector X in Y.

Parameters:
In Trans -If true, solve transpose problem.
In X - A Epetra_MultiVector of dimension NumVectors to solve for.
Out Y -A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns:
Integer error code, set to 0 if successful.

Definition at line 342 of file Ifpack_CrsIct.cpp.

References Epetra_MultiVector::Multiply(), Epetra_MultiVector::NumVectors(), and Epetra_CrsMatrix::Solve().

Referenced by ApplyInverse(), and Condest().


The documentation for this class was generated from the following files:
Generated on Thu Sep 18 12:37:08 2008 for IFPACK by doxygen 1.3.9.1