|
|
◆ csc_ilu()
| void csc_ilu |
( |
int |
n, |
|
|
const double |
val[], |
|
|
const int |
colptr[], |
|
|
const int |
rowind[], |
|
|
int |
base, |
|
|
int |
p, |
|
|
int |
nnz2, |
|
|
double |
val2[], |
|
|
int |
colptr2[], |
|
|
int |
rowind2[], |
|
|
int |
base2, |
|
|
double |
d[], |
|
|
double |
work[], |
|
|
int |
iwork[], |
|
|
int * |
info |
|
) |
| |
Incomplete LU decomposition with level (ILU(p)) (CSC)
- Purpose
- This routine computes the incomplete LU decomposition of the coefficient matrix A of the sparse linear equations. where R is the difference from the complete LU decomposition. Assuming that R is small, the following preconditioner matrix is obtained by solving the equations using this decomposition. In order to suppress fill-ins, new non-zero elements generated during decomposition with level values larger than the specified value p are dropped. This is called as ILU(p).
Level values are initialized as follows. lev(a(i,j)) = 0 (if a(i,j) != 0 or i = j),
= ∞ (if a(i,j) = 0)
In the elimination process of LU decomposition, the operation a(i,j) = a(i,j) - a(i,k)*a(k,j) is required. At the time, the level value of a(i,j) is updated as follows. lev(a(i,j)) = min{ lev(a(i,j)), lev(a(i,k)) + lev(a(k,j)) + 1 }
Due to this algorithm, the level of the position of the non-zero element of A before decomposition remains 0 until the last. When decomposition is completed, elements with level values greater than p are dropped.
This routine outputs the lower triangular matrix L and the upper triangular matrix U in val2[], rowind2[] and colptr2[]. The diagonal elements of U are copied to d[]. val2[], rowind2[], colptr2[] and d[] will be used by csc_ilu_solve().
- Parameters
-
| [in] | n | Dimension of matrix A. (n >= 0) (If n = 0, returns without computation) |
| [in] | val[] | Array val[lval] (lval >= nnz) (nnz is number of non-zero elements of matrix A)
Values of non-zero elements of matrix A. |
| [in] | colptr[] | Array colptr[lcolptr] (lcolptr >= n + 1)
Column pointers of matrix A. |
| [in] | rowind[] | Array rowind[lrowind] (lrowind >= nnz)
Row indices of matrix A. |
| [in] | base | Indexing of colptr[] and rowind[].
= 0: Zero-based (C style) indexing: Starting index is 0.
= 1: One-based (Fortran style) indexing: Starting index is 1. |
| [in] | p | The value of level p of ILU(p) algorithm. (p >= 0)
Note: it is possible to set p = 0. However, it is recommended to use csc_ilu0() since it is faster. |
| [in] | nnz2 | Maximum number of non-zero elements of matrix L*U. If the number of non-zero elements exceeds this during the factorization, the routine will terminate with error (info = -7). |
| [out] | val2[] | Array val2[lval2] (lval2 >= nnz2)
Values of non-zero elements of the lower triangular matrix L and the upper triangular matrix U. |
| [out] | colptr2[] | Array colptr2[lcolptr2] (lcolptr2 >= n + 1)
Column pointers of matrix L*U. |
| [out] | rowind2[] | Array rowind2[lrowind2] (lrowind2 >= nnz2)
Row indices of matrix L*U. (nnz is number of non-zero elements) |
| [in] | base2 | Indexing of colptr2[] and rowind2[].
= 0: Zero-based (C style) indexing: Starting index is 0.
= 1: One-based (Fortran style) indexing: Starting index is 1. |
| [out] | d[] | Array d[ld] (ld >= n)
The diagonal elements of upper triangular matrix U. |
| [out] | work[] | Array work[lwork] (lwork >= n)
Work array. |
| [out] | iwork[] | Array iwork[liwork] (liwork >= nnz2 + 2*n)
Integer work array. |
| [out] | info | = 0: Successful exit.
= i < 0: The (-i)-th argument is invalid.
= j > 0: Matrix is singular (j-th diagonal element is zero). |
|