#include <Ifpack_CrsRiluk.h>
Public Member Functions | |
Ifpack_CrsRiluk (const Ifpack_IlukGraph &Graph_in) | |
Ifpack_CrsRiluk constuctor with variable number of indices per row. | |
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. | |
int | InitValues (const Epetra_VbrMatrix &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. | |
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 ¶meterlist, 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. | |
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). | |
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. | |
int | Condest (bool Trans, double &ConditionNumberEstimate) const |
Returns the maximum over all the condition number estimate for each local ILU set of factors. | |
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. | |
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. | |
const Ifpack_IlukGraph & | Graph () const |
Returns the address of the Ifpack_IlukGraph associated with this factored matrix. | |
const Epetra_CrsMatrix & | L () const |
Returns the address of the L factor associated with this factored matrix. | |
const Epetra_Vector & | D () const |
Returns the address of the D factor associated with this factored matrix. | |
const Epetra_CrsMatrix & | U () 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. | |
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_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. | |
const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. | |
const Epetra_Comm & | Comm () 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 | |
ostream & | operator<< (ostream &os, const Ifpack_CrsRiluk &A) |
<< operator will work for Ifpack_CrsRiluk. |
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:
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.
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.
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 .
The condition of a matrix , called
, is defined as
in some appropriate norm
.
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
or
may be meaningless.
The -norm of a vector
is defined as the maximum of the absolute values of the vector entries, and the
-norm of a matrix C is defined as
. A crude lower bound for the
is
where
. It is a lower bound because
.
For our purposes, we want to estimate , where
and
are our incomplete factors. Edmond in his Ph.D. thesis demonstrates that
provides an effective estimate for
. Furthermore, since finding
such that
is a basic kernel for applying the preconditioner, computing this estimate of
is performed by setting
, calling the solve kernel to compute
and then computing
.
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 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
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 with
,
, where
is the matrix dimension and
returns the sign of the diagonal entry. This has the effect of forcing the diagonal values to have minimal magnitude of
and to increase each by an amount proportional to
, 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:
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.
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.
In | Graph_in - Graph generated by Ifpack_IlukGraph. |
Ifpack_CrsRiluk::Ifpack_CrsRiluk | ( | const Ifpack_CrsRiluk & | Matrix | ) |
Copy constructor.
virtual Ifpack_CrsRiluk::~Ifpack_CrsRiluk | ( | ) | [virtual] |
Ifpack_CrsRiluk Destructor.
bool Ifpack_CrsRiluk::Allocated | ( | ) | const [inline, protected] |
int Ifpack_CrsRiluk::Apply | ( | const Epetra_MultiVector & | X, | |
Epetra_MultiVector & | Y | |||
) | const [inline] |
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.
In | X - A Epetra_MultiVector of dimension NumVectors to multiply with matrix. | |
Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. |
int Ifpack_CrsRiluk::ApplyInverse | ( | const Epetra_MultiVector & | X, | |
Epetra_MultiVector & | Y | |||
) | const [inline] |
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.
In | X - A Epetra_MultiVector of dimension NumVectors to solve for. | |
Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. |
int Ifpack_CrsRiluk::BlockGraph2PointGraph | ( | const Epetra_CrsGraph & | BG, | |
Epetra_CrsGraph & | PG, | |||
bool | Upper | |||
) | [protected] |
const Epetra_Comm& Ifpack_CrsRiluk::Comm | ( | ) | const [inline] |
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
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.
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. |
const Epetra_Vector& Ifpack_CrsRiluk::D | ( | ) | const [inline] |
Returns the address of the D factor associated with this factored matrix.
int Ifpack_CrsRiluk::Factor | ( | ) |
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:
bool Ifpack_CrsRiluk::Factored | ( | ) | const [inline] |
If factor is completed, this query returns true, otherwise it returns false.
double Ifpack_CrsRiluk::GetAbsoluteThreshold | ( | ) | [inline] |
Get absolute threshold value.
Epetra_CombineMode Ifpack_CrsRiluk::GetOverlapMode | ( | ) | [inline] |
Get overlap mode type.
double Ifpack_CrsRiluk::GetRelativeThreshold | ( | ) | [inline] |
Get relative threshold value.
double Ifpack_CrsRiluk::GetRelaxValue | ( | ) | [inline] |
Get RILU(k) relaxation parameter.
const Ifpack_IlukGraph& Ifpack_CrsRiluk::Graph | ( | ) | const [inline] |
Returns the address of the Ifpack_IlukGraph associated with this factored matrix.
bool Ifpack_CrsRiluk::HasNormInf | ( | ) | const [inline] |
Returns false because this class cannot compute an Inf-norm.
int Ifpack_CrsRiluk::IndexBase | ( | ) | const [inline] |
Returns the index base for row and column indices for this graph.
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.
In | A - User matrix to be factored. |
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.
In | A - User matrix to be factored. |
const Epetra_CrsMatrix& Ifpack_CrsRiluk::L | ( | ) | const [inline] |
Returns the address of the L factor associated with this factored matrix.
const char* Ifpack_CrsRiluk::Label | ( | ) | const [inline] |
Returns a character string describing the operator.
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.
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. |
double Ifpack_CrsRiluk::NormInf | ( | ) | const [inline] |
Returns 0.0 because this class cannot compute Inf-norm.
virtual int Ifpack_CrsRiluk::NumGlobalBlockDiagonals | ( | ) | const [inline, virtual] |
Returns the number of diagonal entries found in the global input graph.
int Ifpack_CrsRiluk::NumGlobalCols | ( | ) | const [inline] |
Returns the number of global matrix columns.
int Ifpack_CrsRiluk::NumGlobalNonzeros | ( | ) | const [inline] |
Returns the number of nonzero entries in the global graph.
int Ifpack_CrsRiluk::NumGlobalRows | ( | ) | const [inline] |
Returns the number of global matrix rows.
virtual int Ifpack_CrsRiluk::NumMyBlockDiagonals | ( | ) | const [inline, virtual] |
Returns the number of diagonal entries found in the local input graph.
int Ifpack_CrsRiluk::NumMyCols | ( | ) | const [inline] |
Returns the number of local matrix columns.
virtual int Ifpack_CrsRiluk::NumMyDiagonals | ( | ) | const [inline, virtual] |
Returns the number of nonzero diagonal values found in matrix.
int Ifpack_CrsRiluk::NumMyNonzeros | ( | ) | const [inline] |
Returns the number of nonzero entries in the local graph.
int Ifpack_CrsRiluk::NumMyRows | ( | ) | const [inline] |
Returns the number of local matrix rows.
const Epetra_Map& Ifpack_CrsRiluk::OperatorDomainMap | ( | ) | const [inline] |
Returns the Epetra_Map object associated with the domain of this operator.
const Epetra_Map& Ifpack_CrsRiluk::OperatorRangeMap | ( | ) | const [inline] |
Returns the Epetra_Map object associated with the range of this operator.
void Ifpack_CrsRiluk::SetAbsoluteThreshold | ( | double | Athresh | ) | [inline] |
Set absolute threshold value.
int Ifpack_CrsRiluk::SetAllocated | ( | bool | Flag | ) | [inline, protected] |
void Ifpack_CrsRiluk::SetFactored | ( | bool | Flag | ) | [inline, protected] |
void Ifpack_CrsRiluk::SetOverlapMode | ( | Epetra_CombineMode | OverlapMode | ) | [inline] |
Set overlap mode type.
int Ifpack_CrsRiluk::SetParameters | ( | const Teuchos::ParameterList & | parameterlist, | |
bool | cerr_warning_if_unused = false | |||
) |
Set parameters using a Teuchos::ParameterList object.
void Ifpack_CrsRiluk::SetRelativeThreshold | ( | double | Rthresh | ) | [inline] |
Set relative threshold value.
void Ifpack_CrsRiluk::SetRelaxValue | ( | double | RelaxValue | ) | [inline] |
Set RILU(k) relaxation parameter.
int Ifpack_CrsRiluk::SetUseTranspose | ( | bool | UseTranspose_in | ) | [inline] |
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.
In | UseTranspose_in -If true, multiply by the transpose of operator, otherwise just use operator. |
void Ifpack_CrsRiluk::SetValuesInitialized | ( | bool | Flag | ) | [inline, protected] |
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).
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. |
const Epetra_CrsMatrix& Ifpack_CrsRiluk::U | ( | ) | const [inline] |
Returns the address of the L factor associated with this factored matrix.
bool Ifpack_CrsRiluk::UseTranspose | ( | ) | const [inline] |
Returns the current UseTranspose setting.
bool Ifpack_CrsRiluk::ValuesInitialized | ( | ) | const [inline] |
If values have been initialized, this query returns true, otherwise it returns false.
ostream& operator<< | ( | ostream & | os, | |
const Ifpack_CrsRiluk & | A | |||
) | [friend] |
<< operator will work for Ifpack_CrsRiluk.