|
|
◆ ZCscIlu()
| Sub ZCscIlu |
( |
N As |
Long, |
|
|
Val() As |
Complex, |
|
|
Colptr() As |
Long, |
|
|
Rowind() As |
Long, |
|
|
P As |
Long, |
|
|
Val2() As |
Complex, |
|
|
Colptr2() As |
Long, |
|
|
Rowind2() As |
Long, |
|
|
D() As |
Complex, |
|
|
Optional Info As |
Long, |
|
|
Optional Base As |
Long = -1, |
|
|
Optional Base2 As |
Long = 0 |
|
) |
| |
Initialize incomplete LU decomposition with level (ILU(p)) preconditioner (Complex matrices) (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_ILUSolve().
- Parameters
-
| [in] | N | Dimension 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] | 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. |
| [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.78+0.16i -0.9-1.46i 0.48-1.08i )
A = ( 0.73+0.63i 1.58-1.24 -0.41-0.91i )
( 0.23-1.37i 0.79+0.64i -0.73-1.5i )
( 0.2126-0.2904i )
B = ( -0.3028+0.3346i )
( -1.2905-1.0346i )
Sub Ex_ZFgmres_Ilu_Csc()
Const N = 3, Nnz = N * N, NnzM = Nnz, Tol = 0.0000000001 '1.0e-10
Dim A(Nnz - 1) As Complex, Ia(N) As Long, Ja(Nnz - 1) As Long
Dim B(N - 1) As Complex, X(N - 1) As Complex
Dim XX(N - 1) As Complex, YY(N - 1) As Complex
Dim Iter As Long, Res As Double, IRev As Long, Info As Long, I As Long
A(0) = Cmplx(0.78, 0.16): A(1) = Cmplx(0.73, 0.63): A(2) = Cmplx(0.23, -1.37): A(3) = Cmplx(-0.9, -1.46): A(4) = Cmplx(1.58, -1.24): A(5) = Cmplx(0.79, 0.64): A(6) = Cmplx(0.48, -1.08): A(7) = Cmplx(-0.41, -0.91): A(8) = Cmplx(-0.73, -1.5)
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) = Cmplx(0.2126, -0.2904): B(1) = Cmplx(-0.3028, 0.3346): B(2) = Cmplx(-1.2905, -1.0346)
Dim M(NnzM - 1) As Complex, Im(N) As Long, Jm(NnzM - 1) As Long
Dim P As Long, D(N - 1) As Complex
P = 3
Call ZCscIlu(N, A(), Ia(), Ja(), P, M(), Im(), Jm(), D(), Info)
If Info <> 0 Then Debug.Print "Ilu Info =" + Str(Info)
IRev = 0
Do
Call ZFgmres_r(N, B(), X(), Info, XX(), YY(), IRev, Iter, Res)
If IRev = 1 Then '- Matvec
ElseIf IRev = 3 Then '- Psolve
Call ZCscIluSolve("N", N, M(), Im(), Jm(), 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 ="
Debug.Print "(" + CStr( Creal(X(0))) + "," + CStr( Cimag(X(0))) + ")"
Debug.Print "(" + CStr( Creal(X(1))) + "," + CStr( Cimag(X(1))) + ")"
Debug.Print "(" + CStr( Creal(X(2))) + "," + CStr( Cimag(X(2))) + ")"
Debug.Print "Iter =" + Str(Iter) + ", Res =" + Str(Res) + ", Info =" + Str(Info)
End Sub
Function Cmplx(R As Double, Optional I As Double=0) As Complex Building complex number
Function Cimag(A As Complex) As Double Imaginary part of complex number
Function Creal(A As Complex) As Double Real part of complex number
Sub CscZusmv(Trans As String, M As Long, N As Long, Alpha As Complex, Val() As Complex, Colptr() As Long, Rowind() As Long, X() As Complex, Beta As Complex, Y() As Complex, Optional Info As Long, Optional Base As Long=-1, Optional IncX As Long=1, Optional IncY As Long=1) y <- αAx + βy, y <- αATx + βy or y <- αAHx + βy (Complex matrices) (CSC)
Sub ZFgmres_r(N As Long, B() As Complex, X() As Complex, Info As Long, XX() As Complex, YY() As Complex, 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 (Complex matrices...
Sub ZCscIluSolve(ByVal Trans As String, N As Long, Val() As Complex, Colptr() As Long, Rowind() As Long, D() As Complex, B() As Complex, X() As Complex, Optional Info As Long, Optional Base As Long=-1) Incomplete LU decomposition (ILU) preconditioner (Complex matrices) (CSC)
Sub ZCscIlu(N As Long, Val() As Complex, Colptr() As Long, Rowind() As Long, P As Long, Val2() As Complex, Colptr2() As Long, Rowind2() As Long, D() As Complex, Optional Info As Long, Optional Base As Long=-1, Optional Base2 As Long=0) Initialize incomplete LU decomposition with level (ILU(p)) preconditioner (Complex matrices) (CSC)
- Example Results
X =
(0.59,-0.28)
(-0.2,-3.99999999999999E-02)
(0.24,-0.49)
Iter = 1, Res = 2.58729587846748E-16, Info = 0
|