NAG Library Routine Document

f11jaf  (real_symm_precon_ichol)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

f11jaf 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 f11gef or f11jcf.

2
Specification

Fortran Interface
Subroutine f11jaf ( n, nnz, a, la, irow, icol, lfill, dtol, mic, dscale, pstrat, ipiv, istr, nnzc, npivm, iwork, liwork, ifail)
Integer, Intent (In):: n, nnz, la, lfill, liwork
Integer, Intent (Inout):: irow(la), icol(la), ipiv(n), ifail
Integer, Intent (Out):: istr(n+1), nnzc, npivm, iwork(liwork)
Real (Kind=nag_wp), Intent (In):: dtol, dscale
Real (Kind=nag_wp), Intent (Inout):: a(la)
Character (1), Intent (In):: mic, pstrat
C Header Interface
#include nagmk26.h
void  f11jaf_ ( const Integer *n, const Integer *nnz, double a[], const Integer *la, Integer irow[], Integer icol[], const Integer *lfill, const double *dtol, const char *mic, const double *dscale, const char *pstrat, Integer ipiv[], Integer istr[], Integer *nnzc, Integer *npivm, Integer iwork[], const Integer *liwork, Integer *ifail, const Charlen length_mic, const Charlen length_pstrat)

3
Description

f11jaf computes an incomplete Cholesky factorization (see Meijerink and Van der Vorst (1977)) of a real sparse symmetric n by 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 with one of the symmetric iterative solvers f11gef or f11jcf.
The decomposition is written in the form
A=M+R  
where
M=PLDLTPT  
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 – IntegerInput
On entry: n, the order of the matrix A.
Constraint: n1.
2:     nnz – IntegerInput
On entry: the number of nonzero elements in the lower triangular part of the matrix A.
Constraint: 1nnzn×n+1/2.
3:     ala – Real (Kind=nag_wp) arrayInput/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 routine f11zbf 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 – IntegerInput
On entry: the dimension of the arrays a, irow and icol as declared in the (sub)program from which f11jaf is called. These arrays must be of sufficient size to store both A (nnz elements) and C (nnzc elements).
Constraint: la2×nnz.
5:     irowla – Integer arrayInput/Output
6:     icolla – Integer arrayInput/Output
On entry: the row and column indices of the nonzero elements supplied in a.
Constraints:
irow and icol must satisfy these constraints (which may be imposed by a call to f11zbf):
  • 1irowin and 1icoliirowi, for i=1,2,,nnz;
  • irowi-1<irowi or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nnz.
On exit: the row and column indices of the nonzero elements returned in a.
7:     lfill – IntegerInput
On entry: if lfill0 its value is the maximum level of fill allowed in the decomposition (see Section 9.2). A negative value of lfill indicates that dtol will be used to control the fill instead.
8:     dtol – Real (Kind=nag_wp)Input
On entry: if lfill<0, dtol is used as a drop tolerance to control the fill-in (see Section 9.2); otherwise dtol is not referenced.
Constraint: if lfill<0, dtol0.0.
9:     mic – Character(1)Input
On entry: indicates whether or not the factorization should be modified to preserve row sums (see Section 9.3).
mic='M'
The factorization is modified.
mic='N'
The factorization is not modified.
Constraint: mic='M' or 'N'.
10:   dscale – Real (Kind=nag_wp)Input
On entry: the diagonal scaling parameter. 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.3.
11:   pstrat – Character(1)Input
On entry: specifies the pivoting strategy to be adopted.
pstrat='N'
No pivoting is carried out.
pstrat='M'
Diagonal pivoting aimed at minimizing fill-in is carried out, using the Markowitz strategy.
pstrat='U'
Diagonal pivoting is carried out according to the user-defined input value of ipiv.
Suggested value: pstrat='M'.
Constraint: pstrat='N', 'M' or 'U'.
12:   ipivn – Integer arrayInput/Output
On entry: if pstrat='U', ipivi 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='U', ipiv must contain a valid permutation of the integers on [1,n].
On exit: the pivot indices. If ipivi=j then the diagonal element in row j was used as the pivot at elimination stage i.
13:   istrn+1 – Integer arrayOutput
On exit: istri, for i=1,2,,n, is the starting address in the arrays a, irow and icol of row i of the matrix C. istrn+1 is the address of the last nonzero element in C plus one.
14:   nnzc – IntegerOutput
On exit: the number of nonzero elements in the lower triangular matrix C.
15:   npivm – IntegerOutput
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 f11jaf again with an increased value of either lfill or dscale. See also Section 9.4.
16:   iworkliwork – Integer arrayWorkspace
17:   liwork – IntegerInput
On entry: the dimension of the array iwork as declared in the (sub)program from which f11jaf is called.
Constraints:
the minimum permissible value of liwork depends on lfill as follows:
  • if lfill0, liwork2×la-3×nnz+7×n+1;
  • if lfill<0, liworkla-nnz+7×n+1.
18:   ifail – IntegerInput/Output
On entry: ifail must be set to 0, -1​ or ​1. If you are unfamiliar with this argument you should refer to Section 3.4 in How to Use the NAG Library and its Documentation for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1​ or ​1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this argument, the recommended value is 0. When the value -1​ or ​1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6
Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
ifail=1
On entry,n<1,
ornnz<1,
ornnz>n×n+1/2,
orla<2×nnz,
ordtol<0.0,
ormic'M' or 'N',
orpstrat'N', 'M' or 'U',
orliwork<2×la-3×nnz+7×n+1, and lfill0,
orliwork<la-nnz+7×n+1, and lfill<0.
ifail=2
On entry, the arrays irow and icol fail to satisfy the following constraints:
  • 1irowin and 1icoliirowi, for i=1,2,,nnz;
  • irowi-1<irowi, or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nnz.
Therefore a nonzero element has been supplied which does not lie in the lower triangular part of A, is out of order, or has duplicate row and column indices. Call f11zbf to reorder and sum or remove duplicates.
ifail=3
On entry, pstrat='U', but ipiv does not represent a valid permutation of the integers in 1,n. An input value of ipiv is either out of range or repeated.
ifail=4
la is too small, resulting in insufficient storage space for fill-in elements. The decomposition has been terminated before completion. Either increase la or reduce the amount of fill by setting pstrat='M', reducing lfill, or increasing dtol.
ifail=5 (f11zbf)
A serious error has occurred in an internal call to the specified routine. Check all subroutine calls and array sizes. Seek expert help.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.9 in How to Use the NAG Library and its Documentation for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.8 in How to Use the NAG Library and its Documentation for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 3.7 in How to Use the NAG Library and its Documentation for further information.

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 f11jaf is used in combination with f11gef or f11jcf, 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

f11jaf makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9
Further Comments

9.1
Timing

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

9.2
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=maxke,kc+1,  
where ke is the level of fill of the element being eliminated, and kc 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 aij occurring in row i and column j will not be included if
aij<dtol×aiiajj.  
For either method of control, any elements which are not included are discarded if mic='N', or subtracted from the diagonal element in the elimination row if mic='M'.

9.3
Choice of Arguments

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='N' 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='M' to see if any improvement can be obtained by using modified incomplete Cholesky.
f11jaf 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 routine.
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 (see Meijerink and Van der Vorst (1977)). In this case a good preconditioner can generally be expected by setting lfill=0, mic='M' and dscale=0.0.
For certain mesh-based problems involving M-matrices it can be shown in theory that setting mic='M', 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 (see Chan (1991)). In practise this property often holds even with dscale=0.0, although an improvement in condition can result from increasing dscale slightly (see Van der Vorst (1990)).
Some illustrations of the application of f11jaf 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.4
Direct Solution of positive definite Systems

Although it is not their primary purpose, f11jaf and f11jbf may be used together to obtain a direct solution to a symmetric positive definite linear system. To achieve this the call to f11jbf should be preceded by a complete Cholesky factorization
A=PLDLTPT=M.  
A complete factorization is obtained from a call to f11jaf 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 f11jaf and f11jbf as a direct method is illustrated in Section 10 in f11jbf.

10
Example

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

10.1
Program Text

Program Text (f11jafe.f90)

10.2
Program Data

Program Data (f11jafe.d)

10.3
Program Results

Program Results (f11jafe.r)

© The Numerical Algorithms Group Ltd, Oxford, UK. 2017