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

◆ CsrIlu0()

Sub CsrIlu0 ( N As  Long,
Val() As  Double,
Rowptr() As  Long,
Colind() As  Long,
Val2() As  Double,
D() As  Double,
Optional Info As  Long,
Optional Base As  Long = -1 
)

Initialize incomplete LU decomposition without fill-in (ILU0) preconditioner (CSR)

Purpose
This routine computes the incomplete LU decomposition, without fill-in, 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
This routine outputs the lower triangular matrix L and the upper triangular matrix U in Val2(). And the diagonal elements of U are copied to D(). Val2() and D() will be used by CsrIluSolve().
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]Rowptr()Array Rowptr(LRowptr - 1) (LRowptr >= N + 1)
Row pointers of matrix A.
[in]Colind()Array Colind(LColind - 1) (LColind >= Nnz)
Column indices of matrix A.
[out]Val2()Array Val2(LVal2 - 1) (LVal2 >= Nnz)
Values of non-zero elements of the lower triangular matrix L and the upper triangular matrix U. (Values to be stored in the same location of lower and upper triangular elements of A)
[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 Rowptr() and Colind().
= 0: Zero-based (C style) indexing: Starting index is 0.
= 1: One-based (Fortran style) indexing: Starting index is 1.
(default: Assumes 1 if Rowptr(0) = 1, 0 otherwise)
Example Program
Solve the system of linear equations Ax = B by FGMRES method with DS 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_Ilu0_Csr()
Const N = 3, Nnz = N * N, 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.11: A(2) = -0.93: A(3) = -0.32: A(4) = 0.81: A(5) = 0.37: A(6) = -0.8: A(7) = -0.92: 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(Nnz - 1) As Double, D(N - 1) As Double
Call CsrIlu0(N, A(), Ia(), Ja(), M(), D(), Info)
If Info <> 0 Then Debug.Print "Ilu0 Info =" + Str(Info)
IRev = 0
Do
Call Fgmres_r(N, B(), X(), Info, XX(), YY(), IRev, Iter, Res)
If IRev = 1 Then '- Matvec
Call CsrDusmv("N", N, N, 1, A(), Ia(), Ja(), XX(), 0, YY())
ElseIf IRev = 3 Then '- Psolve
Call CsrIluSolve("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 CsrDusmv(Trans As String, M As Long, N As Long, Alpha As Double, Val() As Double, Rowptr() As Long, Colind() 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 (CSR)
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 CsrIlu0(N As Long, Val() As Double, Rowptr() As Long, Colind() As Long, Val2() As Double, D() As Double, Optional Info As Long, Optional Base As Long=-1)
Initialize incomplete LU decomposition without fill-in (ILU0) preconditioner (CSR)
Sub CsrIluSolve(ByVal Trans As String, N As Long, Val() As Double, Rowptr() As Long, Colind() 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) (CSR)
Example Results
X = 0.86 0.64 0.51
Iter = 1, Res = 2.00597268669339E-16, Info = 0