NAG CL Interface
f11jac (real_​symm_​precon_​ichol)

Settings help

CL Name Style:


1 Purpose

f11jac computes an incomplete Cholesky factorization of a real sparse symmetric matrix, represented in symmetric coordinate storage format. This factorization may be used as a preconditioner in combination with f11jcc.

2 Specification

#include <nag.h>
void  f11jac (Integer n, Integer nnz, double *a[], Integer *la, Integer *irow[], Integer *icol[], Integer lfill, double dtol, Nag_SparseSym_Fact mic, double dscale, Nag_SparseSym_Piv pstrat, Integer ipiv[], Integer istr[], Integer *nnzc, Integer *npivm, Nag_Sparse_Comm *comm, NagError *fail)
The function may be called by the names: f11jac, nag_sparse_real_symm_precon_ichol or nag_sparse_sym_chol_fac.

3 Description

This function computes an incomplete Cholesky factorization (see Meijerink and Van der Vorst (1977)) of a real sparse symmetric n × n matrix A . It is designed specifically for positive definite matrices, but may also work for some mildly indefinite cases. The factorization is intended primarily for use as a preconditioner for the symmetric iterative solver f11jcc.
The decomposition is written in the form
A = M + R  
where
M = PLDLT PT  
and P is a permutation matrix, L is lower triangular with unit diagonal elements, D is diagonal and R is a remainder matrix.
The amount of fill-in occurring in the factorization can vary from zero to complete fill, and can be controlled by specifying either the maximum level of fill lfill, or the drop tolerance dtol. The factorization may be modified in order to preserve row sums, and the diagonal elements may be perturbed to ensure that the preconditioner is positive definite. Diagonal pivoting may optionally be employed, either with a user-defined ordering, or using the Markowitz strategy (see Markowitz (1957)) which aims to minimize fill-in. For further details see Section 9.
The sparse matrix A is represented in symmetric coordinate storage (SCS) format (see Section 2.1.2 in the F11 Chapter Introduction). The array a stores all the nonzero elements of the lower triangular part of A , while arrays irow and icol store the corresponding row and column indices respectively. Multiple nonzero elements may not be specified for the same row and column index.
The preconditioning matrix M is returned in terms of the SCS representation of the lower triangular matrix
C = L + D −1 - I .  

4 References

Chan T F (1991) Fourier analysis of relaxed incomplete factorization preconditioners SIAM J. Sci. Statist. Comput. 12(2) 668–680
Markowitz H M (1957) The elimination form of the inverse and its application to linear programming Management Sci. 3 255–269
Meijerink J and Van der Vorst H (1977) An iterative solution method for linear systems of which the coefficient matrix is a symmetric M-matrix Math. Comput. 31 148–162
Salvini S A and Shaw G J (1995) An evaluation of new NAG Library solvers for large sparse symmetric linear systems NAG Technical Report TR1/95
Van der Vorst H A (1990) The convergence behaviour of preconditioned CG and CG-S in the presence of rounding errors Lecture Notes in Mathematics (eds O Axelsson and L Y Kolotilina) 1457 Springer–Verlag

5 Arguments

1: n Integer Input
On entry: the order of the matrix A .
Constraint: n1 .
2: nnz Integer Input
On entry: the number of nonzero elements in the lower triangular part of the matrix A .
Constraint: 1 < nnz n × (n+1) / 2 .
3: a[la] double * Input/Output
On entry: the nonzero elements in the lower triangular part of the matrix A , ordered by increasing row index, and by increasing column index within each row. Multiple entries for the same row and column indices are not permitted. The function f11zbc may be used to order the elements in this way.
On exit: the first nnz elements of a contain the nonzero elements of A and the next nnzc elements contain the elements of the lower triangular matrix C . Matrix elements are ordered by increasing row index, and by increasing column index within each row.
4: la Integer * Input/Output
On entry: the dimension of the arrays a, irow and icol.
These arrays must be of sufficient size to store both A (nnz elements) and C (nnzc elements); for this reason the length of the arrays may be changed internally by calls to realloc. It is, therefore, imperative that these arrays are allocated using malloc and not declared as automatic arrays.
On exit: if internal allocation has taken place then la is set to nnz+nnzc , otherwise it remains unchanged.
Constraint: la 2 × nnz .
5: irow[la] Integer * Input/Output
6: icol[la] Integer * Input/Output
On entry: the row and column indices of the nonzero elements supplied in a.
Constraints:
  • irow and icol must satisfy the following constraints (which may be imposed by a call to f11zbc):;
  • 1 irow[i] n and 1 icol[i] n , for i=0,1,,nnz-1;
  • irow[i-1] < irow[i] or irow[i-1] = irow[i] and icol[i-1] < icol[i] , for i=1,2,,nnz-1.
On exit: the row and column indices of the nonzero elements returned in a.
7: lfill Integer Input
On entry: if lfill0 its value is the maximum level of fill allowed in the decomposition (see Section 9.1). A negative value of lfill indicates that dtol will be used to control the fill instead.
8: dtol double Input
On entry: if lfill<0 then dtol is used as a drop tolerance to control the fill-in (see Section 9.1). Otherwise dtol is not referenced.
Constraint: if lfill<0 , dtol0.0 .
9: mic Nag_SparseSym_Fact Input
On entry: indicates whether or not the factorization should be modified to preserve row sums (see Section 9.2).
mic=Nag_SparseSym_ModFact
The factorization is modified (MIC).
mic=Nag_SparseSym_UnModFact
The factorization is not modified.
Constraint: mic=Nag_SparseSym_ModFact or Nag_SparseSym_UnModFact.
10: dscale double Input
On entry: the diagonal scaling argument. All diagonal elements are multiplied by the factor (1+dscale) at the start of the factorization. This can be used to ensure that the preconditioner is positive definite. See Section 9.2.
11: pstrat Nag_SparseSym_Piv Input
On entry: specifies the pivoting strategy to be adopted as follows:
  • if pstrat=Nag_SparseSym_NoPiv then no pivoting is carried out;
  • if pstrat=Nag_SparseSym_MarkPiv then diagonal pivoting aimed at minimizing fill-in is carried out, using the Markowitz strategy;
  • if pstrat=Nag_SparseSym_UserPiv then diagonal pivoting is carried out according to the user-defined input value of ipiv.
Suggested value: pstrat=Nag_SparseSym_MarkPiv.
Constraint: pstrat=Nag_SparseSym_NoPiv, Nag_SparseSym_MarkPiv or Nag_SparseSym_UserPiv.
12: ipiv[n] Integer Input/Output
On entry: if pstrat=Nag_SparseSym_UserPiv, then ipiv[i-1] must specify the row index of the diagonal element used as a pivot at elimination stage i . Otherwise ipiv need not be initialized.
Constraint: if pstrat=Nag_SparseSym_UserPiv, then ipiv must contain a valid permutation of the integers on [1,n] .
On exit: the pivot indices. If ipiv[i-1] = j then the diagonal element in row j was used as the pivot at elimination stage i .
13: istr[n+1] Integer Output
On exit: istr[i]-1 , for i=0,1,,n - 1, is the starting address in the arrays a, irow and icol of row i of the matrix C . istr[n]-1 is the address of the last nonzero element in C plus one.
14: nnzc Integer * Output
On exit: the number of nonzero elements in the lower triangular matrix C .
15: npivm Integer * Output
On exit: the number of pivots which were modified during the factorization to ensure that M was positive definite. The quality of the preconditioner will generally depend on the returned value of npivm. If npivm is large the preconditioner may not be satisfactory. In this case it may be advantageous to call f11jac again with an increased value of either lfill or dscale.
16: comm Nag_Sparse_Comm * Input/Output
On entry/exit: a pointer to a structure of type Nag_Sparse_Comm whose members are used by the iterative solver.
17: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_2_INT_ARG_LT
On entry, la=value while nnz=value . These arguments must satisfy la 2 × nnz .
NE_ALLOC_FAIL
Dynamic memory allocation failed.
NE_BAD_PARAM
On entry, argument mic had an illegal value.
On entry, argument pstrat had an illegal value.
NE_INT_2
On entry, nnz=value , n=value .
Constraint: 1 nnz n × (n+1) / 2 .
NE_INT_ARG_LT
On entry, n=value.
Constraint: n1.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
NE_INVALID_ROW_PIVOT
On entry, pstrat=Nag_SparseSym_UserPiv and the array ipiv does not represent a valid permutation of integers in [1,n] . An input value of ipiv is either out of range or repeated.
NE_REAL_INT_ARG_CONS
On entry, dtol=value and lfill=value . These arguments must satisfy dtol0.0 if lfill<0 .
NE_SYMM_MATRIX_DUP
A nonzero element has been supplied which does not lie in the lower triangular part of the matrix A , is out of order, or has duplicate row and column indices, i.e., one or more of the following constraints has been violated:
  • 1 irow[i] n and 1 icol[i] n , for i=0,1,,nnz - 1.
  • irow[i-1] < irow[i] , or
  • irow[i-1] = irow[i] and icol[i-1] < icol[i] , for i=1,2,,nnz - 1.

Call f11zbc to reorder and sum or remove duplicates.

7 Accuracy

The accuracy of the factorization will be determined by the size of the elements that are dropped and the size of any modifications made to the diagonal elements. If these sizes are small then the computed factors will correspond to a matrix close to A . The factorization can generally be made more accurate by increasing lfill, or by reducing dtol with lfill<0 . If f11jac is used in combination with f11jcc, the more accurate the factorization the fewer iterations will be required. However, the cost of the decomposition will also generally increase.

8 Parallelism and Performance

f11jac is not threaded in any implementation.

9 Further Comments

The time taken for a call to f11jac is roughly proportional to nnzc 2 / n .

9.1 Control of Fill-in

If lfill0 the amount of fill-in occurring in the incomplete factorization is controlled by limiting the maximum level of fill-in to lfill. The original nonzero elements of A are defined to be of level 0. The fill level of a new nonzero location occurring during the factorization is defined as:
k = max( k e , k c ) + 1 ,  
where k e is the level of fill of the element being eliminated, and k c is the level of fill of the element causing the fill-in.
If lfill<0 the fill-in is controlled by means of the drop tolerance dtol. A potential fill-in element a ij occurring in row i and column j will not be included if:
| a ij | < dtol × | a ii a jj | .  
For either method of control, any elements which are not included are discarded if mic=Nag_SparseSym_UnModFact, or subtracted from the diagonal element in the elimination row if mic=Nag_SparseSym_ModFact.

9.2 Choice of Parameters

There is unfortunately no choice of the various algorithmic arguments which is optimal for all types of symmetric matrix, and some experimentation will generally be required for each new type of matrix encountered.
If the matrix A is not known to have any particular special properties the following strategy is recommended. Start with lfill=0 , mic=Nag_SparseSym_UnModFact and dscale=0.0 . If the value returned for npivm is significantly larger than zero, i.e., a large number of pivot modifications were required to ensure that M was positive definite, the preconditioner is not likely to be satisfactory. In this case increase either lfill or dscale until npivm falls to a value close to zero. Once suitable values of lfill and dscale have been found try setting mic=Nag_SparseSym_ModFact to see if any improvement can be obtained by using modified incomplete Cholesky.
f11jac is primarily designed for positive definite matrices, but may work for some mildly indefinite problems. If npivm cannot be satisfactorily reduced by increasing lfill or dscale then A is probably too indefinite for this function.
If A has non-positive off-diagonal elements, is nonsingular, and has only non-negative elements in its inverse, it is called an ‘M-matrix’. It can be shown that no pivot modifications are required in the incomplete Cholesky factorization of an M-matrix (Meijerink and Van der Vorst (1977)). In this case a good preconditioner can generally be expected by setting lfill=0 , mic=Nag_SparseSym_ModFact and dscale=0.0 .
For certain mesh-based problems involving M-matrices it can be shown in theory that setting mic=Nag_SparseSym_ModFact, and choosing dscale appropriately can reduce the order of magnitude of the condition number of the preconditioned matrix as a function of the mesh steplength (Chan (1991)). In practise this property often holds even with dscale=0.0 , although an improvement in condition can result from increasing dscale slightly (Van der Vorst (1990)).
Some illustrations of the application of f11jac to linear systems arising from the discretization of two-dimensional elliptic partial differential equations, and to random-valued randomly structured symmetric positive definite linear systems, can be found in Salvini and Shaw (1995).

9.3 Direct Solution of Positive Definite Systems

Although it is not their primary purpose, f11jac and f11jbc may be used together to obtain a direct solution to a symmetric positive definite linear system. To achieve this the call to f11jbc should be preceded by a complete Cholesky factorization
A = PLD LT PT = M .  
A complete factorization is obtained from a call to f11jac with lfill<0 and dtol=0.0, provided npivm=0 on exit. A nonzero value of npivm indicates that A is not positive definite, or is ill-conditioned. A factorization with nonzero npivm may serve as a preconditioner, but will not result in a direct solution. It is, therefore, essential to check the output value of npivm if a direct solution is required.
The use of f11jac and f11jbc as a direct method is illustrated in Section 10 in f11jbc.

10 Example

This example program reads in a symmetric sparse matrix A and calls f11jac to compute an incomplete Cholesky factorization. It then outputs the nonzero elements of both A and C = L + D −1 - I . The call to f11jac has lfill=0 , mic=Nag_SparseSym_UnModFact, dscale=0.0 and pstrat=Nag_SparseSym_MarkPiv, giving an unmodified zero-fill factorization of an unperturbed matrix, with Markowitz diagonal pivoting.

10.1 Program Text

Program Text (f11jace.c)

10.2 Program Data

Program Data (f11jace.d)

10.3 Program Results

Program Results (f11jace.r)