IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Public Member Functions | Protected Member Functions | Friends | List of all members
Ifpack_CrsRiluk Class Reference

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

#include <Ifpack_CrsRiluk.h>

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

Public Member Functions

 Ifpack_CrsRiluk (const Ifpack_IlukGraph &Graph_in)
 Ifpack_CrsRiluk constuctor with variable number of indices per row. More...
 
 Ifpack_CrsRiluk (const Ifpack_CrsRiluk &Matrix)
 Copy constructor.
 
virtual ~Ifpack_CrsRiluk ()
 Ifpack_CrsRiluk Destructor.
 
int InitValues (const Epetra_CrsMatrix &A)
 Initialize L and U with values from user matrix A. More...
 
int InitValues (const Epetra_VbrMatrix &A)
 Initialize L and U with values from user matrix A. More...
 
bool ValuesInitialized () const
 If values have been initialized, this query returns true, otherwise it returns false.
 
void SetRelaxValue (double RelaxValue)
 Set RILU(k) relaxation parameter.
 
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 Factor ()
 Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxation parameters. More...
 
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_CrsRiluk forward/back solve on a Epetra_MultiVector X in Y (works for Epetra_Vectors also). More...
 
int Multiply (bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of multiplying U, D and L in that order on an Epetra_MultiVector X in Y. More...
 
int Condest (bool Trans, double &ConditionNumberEstimate) const
 Returns the maximum over all the condition number estimate for each local ILU set of factors. More...
 
double GetRelaxValue ()
 Get RILU(k) relaxation parameter.
 
double GetAbsoluteThreshold ()
 Get absolute threshold value.
 
double GetRelativeThreshold ()
 Get relative threshold value.
 
Epetra_CombineMode GetOverlapMode ()
 Get overlap mode type.
 
int NumGlobalRows () const
 Returns the number of global matrix rows.
 
int NumGlobalCols () const
 Returns the number of global matrix columns.
 
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global graph.
 
virtual int NumGlobalBlockDiagonals () const
 Returns the number of diagonal entries found in the global input graph.
 
long long NumGlobalRows64 () const
 
long long NumGlobalCols64 () const
 
long long NumGlobalNonzeros64 () const
 
virtual long long NumGlobalBlockDiagonals64 () const
 
int NumMyRows () const
 Returns the number of local matrix rows.
 
int NumMyCols () const
 Returns the number of local matrix columns.
 
int NumMyNonzeros () const
 Returns the number of nonzero entries in the local graph.
 
virtual int NumMyBlockDiagonals () const
 Returns the number of diagonal entries found in the local input graph.
 
virtual int NumMyDiagonals () const
 Returns the number of nonzero diagonal values found in matrix.
 
int IndexBase () const
 Returns the index base for row and column indices for this graph.
 
long long IndexBase64 () const
 
const Ifpack_IlukGraphGraph () const
 returns the address of the Ifpack_IlukGraph associated with this factored matrix.
 
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 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
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More...
 
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 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)
 
int BlockGraph2PointGraph (const Epetra_CrsGraph &BG, Epetra_CrsGraph &PG, bool Upper)
 

Friends

std::ostream & operator<< (std::ostream &os, const Ifpack_CrsRiluk &A)
 << operator will work for Ifpack_CrsRiluk.
 

Detailed Description

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

The Ifpack_CrsRiluk class computes a "Relaxed" ILU factorization with level k fill
of a given Epetra_CrsMatrix.  The factorization
that is produced is a function of several parameters:
  1. The pattern of the matrix - All fill is derived from the original matrix nonzero structure. Level zero fill is defined as the original matrix pattern (nonzero structure), even if the matrix value at an entry is stored as a zero. (Thus it is possible to add entries to the ILU factors by adding zero entries the original matrix.)

  2. Level of fill - Starting with the original matrix pattern as level fill of zero, the next level of fill is determined by analyzing the graph of the previous level and determining nonzero fill that is a result of combining entries that were from previous level only (not the current level). This rule limits fill to entries that are direct decendents from the previous level graph. Fill for level k is determined by applying this rule recursively. For sufficiently large values of k, the fill would eventually be complete and an exact LU factorization would be computed. Level of fill is defined during the construction of the Ifpack_IlukGraph object.

  3. Level of overlap - All Ifpack preconditioners work on parallel distributed memory computers by using the row partitioning the user input matrix to determine the partitioning for local ILU factors. If the level of overlap is set to zero, the rows of the user matrix that are stored on a given processor are treated as a self-contained local matrix and all column entries that reach to off-processor entries are ignored. Setting the level of overlap to one tells Ifpack to increase the size of the local matrix by adding rows that are reached to by rows owned by this processor. Increasing levels of overlap are defined recursively in the same way. For sufficiently large levels of overlap, the entire matrix would be part of each processor's local ILU factorization process. Level of overlap is defined during the construction of the Ifpack_IlukGraph object.

    Once the factorization is computed, applying the factorization (LUy = x) results in redundant approximations for any elements of y that correspond to rows that are part of more than one local ILU factor. The OverlapMode (changed by calling SetOverlapMode()) defines how these redundancies are handled using the Epetra_CombineMode enum. The default is to zero out all values of y for rows that were not part of the original matrix row distribution.

  4. Fraction of relaxation - Ifpack_CrsRiluk computes the ILU factorization row-by-row. As entries at a given row are computed, some number of them will be dropped because they do match the prescribed sparsity pattern. The relaxation factor determines how these dropped values will be handled. If the RelaxValue (changed by calling SetRelaxValue()) is zero, then these extra entries will by dropped. This is a classical ILU approach. If the RelaxValue is 1, then the sum of the extra entries will be added to the diagonal. This is a classical Modified ILU (MILU) approach. If RelaxValue is between 0 and 1, then RelaxValue times the sum of extra entries will be added to the diagonal.

    For most situations, RelaxValue should be set to zero. For certain kinds of problems, e.g., reservoir modeling, there is a conservation principle involved such that any operator should obey a zero row-sum property. MILU was designed for these cases and you should set the RelaxValue to 1. For other situations, setting RelaxValue to some nonzero value may improve the stability of factorization, and can be used if the computed ILU factors are poorly conditioned.

  5. 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_CrsRiluk objects

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

  1. Create Ifpack_CrsRiluk 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_CrsRiluk 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_CrsRiluk constructor.

Definition at line 210 of file Ifpack_CrsRiluk.h.

Constructor & Destructor Documentation

Ifpack_CrsRiluk::Ifpack_CrsRiluk ( const Ifpack_IlukGraph Graph_in)

Ifpack_CrsRiluk constuctor with variable number of indices per row.

Creates a Ifpack_CrsRiluk object and allocates storage.

Parameters
InGraph_in - Graph generated by Ifpack_IlukGraph.

Definition at line 55 of file Ifpack_CrsRiluk.cpp.

References Epetra_BlockMap::DistributedGlobal(), Ifpack_IlukGraph::DomainMap(), Ifpack_IlukGraph::LevelOverlap(), and Zero.

Member Function Documentation

int Ifpack_CrsRiluk::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

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
InX - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
OutY -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 427 of file Ifpack_CrsRiluk.h.

References Multiply(), and UseTranspose().

int Ifpack_CrsRiluk::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

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
InX - A Epetra_MultiVector of dimension NumVectors to solve for.
OutY -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 444 of file Ifpack_CrsRiluk.h.

References Solve(), and UseTranspose().

int Ifpack_CrsRiluk::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
InTrans -If true, solve transpose problem.
OutConditionNumberEstimate - The maximum across all processors of the infinity-norm estimate of the condition number of the inverse of LDU.

Definition at line 596 of file Ifpack_CrsRiluk.cpp.

References Solve().

int Ifpack_CrsRiluk::Factor ( void  )

Compute ILU factors L and 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 361 of file Ifpack_CrsRiluk.cpp.

References Epetra_BlockMap::Comm(), Factored(), Ifpack_IlukGraph::L_Graph(), NumMyCols(), NumMyRows(), Epetra_CrsGraph::RowMap(), Epetra_Comm::SumAll(), and ValuesInitialized().

int Ifpack_CrsRiluk::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
InA - 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 182 of file Ifpack_CrsRiluk.cpp.

References Copy, Epetra_CrsMatrix::DomainMap(), Insert, Ifpack_IlukGraph::OverlapGraph(), Ifpack_IlukGraph::OverlapImporter(), and Epetra_CrsMatrix::RangeMap().

int Ifpack_CrsRiluk::InitValues ( const Epetra_VbrMatrix 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
InA - 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 212 of file Ifpack_CrsRiluk.cpp.

References Copy, Insert, Ifpack_IlukGraph::OverlapGraph(), and Ifpack_IlukGraph::OverlapImporter().

int Ifpack_CrsRiluk::Multiply ( bool  Trans,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const

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

Parameters
InTrans -If true, multiply by L^T, D and U^T in that order.
InX - A Epetra_MultiVector of dimension NumVectors to solve for.
OutY -A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns
Integer error code, set to 0 if successful.

Definition at line 554 of file Ifpack_CrsRiluk.cpp.

Referenced by Apply().

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

Implements Epetra_Operator.

Definition at line 413 of file Ifpack_CrsRiluk.h.

int Ifpack_CrsRiluk::Solve ( bool  Trans,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const

Returns the result of a Ifpack_CrsRiluk forward/back solve on a Epetra_MultiVector X in Y (works for Epetra_Vectors also).

Parameters
InTrans -If true, solve transpose problem.
InX - A Epetra_MultiVector of dimension NumVectors to solve for.
OutY -A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns
Integer error code, set to 0 if successful.

Definition at line 513 of file Ifpack_CrsRiluk.cpp.

Referenced by ApplyInverse(), and Condest().


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