XLPack 7.0
XLPack Numerical Library (Excel VBA) Reference Manual
Loading...
Searching...
No Matches

◆ CscIlu()

Sub CscIlu ( N As  Long,
Val() As  Double,
Colptr() As  Long,
Rowind() As  Long,
P As  Long,
Val2() As  Double,
Colptr2() As  Long,
Rowind2() As  Long,
D() As  Double,
Optional Info As  Long,
Optional Base As  Long = -1,
Optional Base2 As  Long = 0 
)

Initialize incomplete LU decomposition with level (ILU(p)) preconditioner (CSC)

Purpose
This routine computes the incomplete LU decomposition of the coefficient matrix A of the sparse linear equations.
A = L * U + R
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.
M = L * U
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_ILUSolve().

Parameters
[in]NDimension of the matrix A. (N >= 0) (if N = 0, returns without computation)
[in]Val()Array Val(LVal - 1) (LVal >= Nnz) (Nnz is number of non-zero elements)
Values of non-zero elements of matrix A.
[in]Colptr()Array Colptr(LColptr - 1) (LColptr >= N + 1)
Column pointers of matrix A.
[in]Rowind()Array Rowind(LRowind - 1) (LRowind >= Nnz)
Row indices of matrix A.
[in]PThe 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.
[out]Val2()Array Val2(LVal2 - 1) (LVal2 > 0) Values of non-zero elements of the lower triangular matrix L and the upper triangular matrix U (matrix L*U).
Nnz2 = min(LVal2, LRowind2) is the maximum number of non-zero elements of matrix L*U. If the number of non-zero elements exceeds Nnz2 during the factorization, the routine will terminate with error number -6.
[out]Colptr2()Array Colptr2(LColptr2 - 1) (LColptr2 >= N + 1)
Column pointers of matrix L*U.
[out]Rowind2()Array Rowind2(LRowind2 - 1) (LRowind2 > 0)
Row indices of matrix L*U.
[out]D()Array D(LD) (LD >= N)
The diagonal elements of upper triangular matrix U.
[out]Info(Optional)
= 0: Successful exit.
= i < 0: The (-i)-th argument is invalid.
= j > 0: Matrix is singular (j-th diagonal element is zero).
[in]Base(Optional)
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.
(default: Assumes 1 if Colptr(0) = 1, 0 otherwise)
[in]Base2(Optional)
Indexing of Colptr2() and Rowind2(). (default = 0)
= 0: Zero-based (C style) indexing: Starting index is 0.
= 1: One-based (Fortran style) indexing: Starting index is 1.
Example Program
Solve the system of linear equations Ax = B by FGMRES method with ILU(p) preconditioner, where
( 0.2 -0.11 -0.93 ) ( -0.3727 )
A = ( -0.32 0.81 0.37 ), B = ( 0.4319 )
( -0.8 -0.92 -0.29 ) ( -1.4247 )
Sub Ex_Fgmres_Ilu_Csc()
Const N = 3, Nnz = N * N, NnzM = Nnz, Tol = 0.0000000001 '1.0e-10
Dim A(Nnz - 1) As Double, Ia(N) As Long, Ja(Nnz - 1) As Long
Dim B(N - 1) As Double, X(N - 1) As Double
Dim XX(N - 1) As Double, YY(N - 1) As Double
Dim Iter As Long, Res As Double, IRev As Long, Info As Long
A(0) = 0.2: A(1) = -0.32: A(2) = -0.8: A(3) = -0.11: A(4) = 0.81: A(5) = -0.92: A(6) = -0.93: A(7) = 0.37: A(8) = -0.29
Ia(0) = 0: Ia(1) = 3: Ia(2) = 6: Ia(3) = 9
Ja(0) = 0: Ja(1) = 1: Ja(2) = 2: Ja(3) = 0: Ja(4) = 1: Ja(5) = 2: Ja(6) = 0: Ja(7) = 1: Ja(8) = 2
B(0) = -0.3727: B(1) = 0.4319: B(2) = -1.4247
Dim M(NnzM - 1) As Double, Im(N) As Long, Jm(NnzM - 1) As Long
Dim P As Long, D(N - 1) As Double
P = 3
Call CscIlu(N, A(), Ia(), Ja(), P, M(), Im(), Jm(), D(), Info)
If Info <> 0 Then Debug.Print "Ilu Info =" + Str(Info)
IRev = 0
Do
Call Fgmres_r(N, B(), X(), Info, XX(), YY(), IRev, Iter, Res)
If IRev = 1 Then '- Matvec
Call CscDusmv("N", N, N, 1, A(), Ia(), Ja(), XX(), 0, YY())
ElseIf IRev = 3 Then '- Psolve
Call CscIluSolve("N", N, M(), Ia(), Ja(), D(), YY(), XX(), Info)
If Info <> 0 Then Debug.Print "IluSolve Info =" + Str(Info)
ElseIf IRev = 10 Then '- Check convergence
If Res < Tol Then IRev = 11
End If
Loop While IRev <> 0
Debug.Print "X =", X(0), X(1), X(2)
Debug.Print "Iter = " + CStr(Iter) + ", Res = " + CStr(Res) + ", Info = " + CStr(Info)
End Sub
Sub CscDusmv(Trans As String, M As Long, N As Long, Alpha As Double, Val() As Double, Colptr() As Long, Rowind() As Long, X() As Double, Beta As Double, Y() As Double, Optional Info As Long, Optional Base As Long=-1, Optional IncX As Long=1, Optional IncY As Long=1)
y <- αAx + βy or y <- αATx + βy (CSC)
Sub Fgmres_r(N As Long, B() As Double, X() As Double, Info As Long, XX() As Double, YY() As Double, IRev As Long, Optional Iter As Long, Optional Res As Double, Optional M As Long=0, Optional MaxIter As Long=500)
Solution of linear system Ax = b using generalized minimum residual (FGMRES) method (Reverse communic...
Sub CscIluSolve(ByVal Trans As String, N As Long, Val() As Double, Colptr() As Long, Rowind() As Long, D() As Double, B() As Double, X() As Double, Optional Info As Long, Optional Base As Long=-1)
Incomplete LU decomposition preconditioner (ILU) (CSC)
Sub CscIlu(N As Long, Val() As Double, Colptr() As Long, Rowind() As Long, P As Long, Val2() As Double, Colptr2() As Long, Rowind2() As Long, D() As Double, Optional Info As Long, Optional Base As Long=-1, Optional Base2 As Long=0)
Initialize incomplete LU decomposition with level (ILU(p)) preconditioner (CSC)
Example Results
X = 0.86 0.64 0.51
Iter = 1, Res = 3.40151402401104E-16, Info = 0