XLPack 7.0
XLPack Numerical Library (C API) Reference Manual
Loading...
Searching...
No Matches

◆ z_hsr_ic0()

void z_hsr_ic0 ( char  uplo,
int  n,
const doublecomplex  val[],
const int  rowptr[],
const int  colind[],
int  base,
doublecomplex  val2[],
int  idiag[],
doublecomplex  work[],
int *  info 
)

Incomplete Cholesky decomposition without fill-in (IC0) (Hermitian positive definite matrix) (CSR)

Purpose
This routine computes the incomplete Cholesky decomposition, without fill-in, of the Hermitian positive definite coefficient matrix of the sparse linear equations.
A = L * D * L^H + R or A = U^H * D * U + R
It is called that fill-in occurs when a position that is zero in coefficient matrix A but it becomes non-zero after Cholesky decomposition is performed. To avoid this, such positions in L or U in the above equation are forcibly set to zeros. R shows the differences from the complete Cholesky decomposition arising from this modification. Assuming that R is not large, it is expected that this incomplete decomposition can be used to obtain a good approximate solution of linear equations. This decomposition is used as the preconditioner matrix M.
M = L * D * L^H or M = U^H * D * U
This routine outputs L or U and D in val2[], and indices of diagonal elements in idiag[]. val2[] and idiag[] will be used by z_hsr_ic_solve().
Parameters
[in]uploSpecifies whether the upper or lower triangular part of input matrix A is stored.
= "U": Upper triangular part (to be factorized to U^H*D*U).
= "L": Lower triangular part (to be factorized to L*D*L^H).
[in]nDimension 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. (only the diagonal and lower triangular elements are referenced)
[in]rowptr[]Array rowptr[lrowptr] (lrowptr >= n + 1)
Row pointers of matrix A.
[in]colind[]Array colind[lcolind] (lcolind >= nnz)
Column indices of matrix A.
[in]baseIndexing of rowptr[] and colind[].
= 0: Zero-based (C style) indexing: Starting index is 0.
= 1: One-based (Fortran style) indexing: Starting index is 1.
[out]val2[]Array val2[lval2] (lval2 >= nnz)
Values of non-zero elements of the lower triangular matrix L and the diagonal matrix D. (nnz is same with the number of non-zero elements of A. The values are stored in the same location of diagonal and lower triangular elements of A.)
[out]idiag[]Array idiag[lidiag] (lidiag >= n)
Indices of diagonal elements.
[out]work[]Array work[lwork] (lwork >= n)
Work array. (Not references when uplo = 'U'.)
[out]info= 0: Successful exit.
= i < 0: The (-i)-th argument is invalid.
= i > 0: Matrix is singular (i-th diagonal element is zero).