F11JCF (PDF version)
F11 Chapter Contents
F11 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F11JCF

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

F11JCF solves a real sparse symmetric system of linear equations, represented in symmetric coordinate storage format, using a conjugate gradient or Lanczos method, with incomplete Cholesky preconditioning.

2  Specification

SUBROUTINE F11JCF ( METHOD, N, NNZ, A, LA, IROW, ICOL, IPIV, ISTR, B, TOL, MAXITN, X, RNORM, ITN, WORK, LWORK, IFAIL)
INTEGER  N, NNZ, LA, IROW(LA), ICOL(LA), IPIV(N), ISTR(N+1), MAXITN, ITN, LWORK, IFAIL
REAL (KIND=nag_wp)  A(LA), B(N), TOL, X(N), RNORM, WORK(LWORK)
CHARACTER(*)  METHOD

3  Description

F11JCF solves a real sparse symmetric linear system of equations
Ax=b,
using a preconditioned conjugate gradient method (see Meijerink and Van der Vorst (1977)), or a preconditioned Lanczos method based on the algorithm SYMMLQ (see Paige and Saunders (1975)). The conjugate gradient method is more efficient if A is positive definite, but may fail to converge for indefinite matrices. In this case the Lanczos method should be used instead. For further details see Barrett et al. (1994).
F11JCF uses the incomplete Cholesky factorization determined by F11JAF as the preconditioning matrix. A call to F11JCF must always be preceded by a call to F11JAF. Alternative preconditioners for the same storage scheme are available by calling F11JEF.
The matrix A, and the preconditioning matrix M, are represented in symmetric coordinate storage (SCS) format (see Section 2.1.2 in the F11 Chapter Introduction) in the arrays A, IROW and ICOL, as returned from F11JAF. The array A holds the nonzero entries in the lower triangular parts of these matrices, while IROW and ICOL hold the corresponding row and column indices.

4  References

Barrett R, Berry M, Chan T F, Demmel J, Donato J, Dongarra J, Eijkhout V, Pozo R, Romine C and Van der Vorst H (1994) Templates for the Solution of Linear Systems: Building Blocks for Iterative Methods SIAM, Philadelphia
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
Paige C C and Saunders M A (1975) Solution of sparse indefinite systems of linear equations SIAM J. Numer. Anal. 12 617–629
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

5  Parameters

1:     METHOD – CHARACTER(*)Input
On entry: specifies the iterative method to be used.
METHOD='CG'
Conjugate gradient method.
METHOD='SYMMLQ'
Lanczos method (SYMMLQ).
Constraint: METHOD='CG' or 'SYMMLQ'.
2:     N – INTEGERInput
On entry: n, the order of the matrix A. This must be the same value as was supplied in the preceding call to F11JAF.
Constraint: N1.
3:     NNZ – INTEGERInput
On entry: the number of nonzero elements in the lower triangular part of the matrix A. This must be the same value as was supplied in the preceding call to F11JAF.
Constraint: 1NNZN×N+1/2.
4:     A(LA) – REAL (KIND=nag_wp) arrayInput
On entry: the values returned in the array A by a previous call to F11DAF.
5:     LA – INTEGERInput
On entry: the dimension of the arrays A, IROW and ICOL as declared in the (sub)program from which F11JCF is called. This must be the same value as was supplied in the preceding call to F11JAF.
Constraint: LA2×NNZ.
6:     IROW(LA) – INTEGER arrayInput
7:     ICOL(LA) – INTEGER arrayInput
8:     IPIV(N) – INTEGER arrayInput
9:     ISTR(N+1) – INTEGER arrayInput
On entry: the values returned in arrays IROW, ICOL, IPIV and ISTR by a previous call to F11JAF.
10:   B(N) – REAL (KIND=nag_wp) arrayInput
On entry: the right-hand side vector b.
11:   TOL – REAL (KIND=nag_wp)Input
On entry: the required tolerance. Let xk denote the approximate solution at iteration k, and rk the corresponding residual. The algorithm is considered to have converged at iteration k if
rkτ×b+Axk.
If TOL0.0, τ=maxε,nε is used, where ε is the machine precision. Otherwise τ=maxTOL,10ε,nε is used.
Constraint: TOL<1.0.
12:   MAXITN – INTEGERInput
On entry: the maximum number of iterations allowed.
Constraint: MAXITN1.
13:   X(N) – REAL (KIND=nag_wp) arrayInput/Output
On entry: an initial approximation to the solution vector x.
On exit: an improved approximation to the solution vector x.
14:   RNORM – REAL (KIND=nag_wp)Output
On exit: the final value of the residual norm rk, where k is the output value of ITN.
15:   ITN – INTEGEROutput
On exit: the number of iterations carried out.
16:   WORK(LWORK) – REAL (KIND=nag_wp) arrayWorkspace
17:   LWORK – INTEGERInput
On entry: the dimension of the array WORK as declared in the (sub)program from which F11JCF is called.
Constraints:
  • if METHOD='CG', LWORK6×N+120;
  • if METHOD='SYMMLQ', LWORK7×N+120.
18:   IFAIL – INTEGERInput/Output
On entry: IFAIL must be set to 0, -1​ or ​1. If you are unfamiliar with this parameter you should refer to Section 3.3 in the Essential Introduction 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 parameter, 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,METHOD'CG' or 'SYMMLQ',
orN<1,
orNNZ<1,
orNNZ>N×N+1/2,
orLA too small,
orTOL1.0,
orMAXITN<1,
orLWORK too small.
IFAIL=2
On entry, the SCS representation of A is invalid. Further details are given in the error message. Check that the call to F11JCF has been preceded by a valid call to F11JAF, and that the arrays A, IROW, and ICOL have not been corrupted between the two calls.
IFAIL=3
On entry, the SCS representation of the preconditioning matrix M is invalid. Further details are given in the error message. Check that the call to F11JCF has been preceded by a valid call to F11JAF and that the arrays A, IROW, ICOL, IPIV and ISTR have not been corrupted between the two calls.
IFAIL=4
The required accuracy could not be obtained. However, a reasonable accuracy has been obtained and further iterations could not improve the result.
IFAIL=5
Required accuracy not obtained in MAXITN iterations.
IFAIL=6
The preconditioner appears not to be positive definite.
IFAIL=7
The matrix of the coefficients appears not to be positive definite (conjugate gradient method only).
IFAIL=8 (F11GDF, F11GEF or F11GFF)
A serious error has occurred in an internal call to one of the specified routines. Check all subroutine calls and array sizes. Seek expert help.

7  Accuracy

On successful termination, the final residual rk=b-Axk, where k=ITN, satisfies the termination criterion
rkτ×b+Axk.
The value of the final residual norm is returned in RNORM.

8  Further Comments

The time taken by F11JCF for each iteration is roughly proportional to the value of NNZC returned from the preceding call to F11JAF. One iteration with the Lanczos method (SYMMLQ) requires a slightly larger number of operations than one iteration with the conjugate gradient method.
The number of iterations required to achieve a prescribed accuracy cannot be easily determined a priori, as it can depend dramatically on the conditioning and spectrum of the preconditioned matrix of the coefficients A-=M-1A.
Some illustrations of the application of F11JCF 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  Example

This example solves a symmetric positive definite system of equations using the conjugate gradient method, with incomplete Cholesky preconditioning.

9.1  Program Text

Program Text (f11jcfe.f90)

9.2  Program Data

Program Data (f11jcfe.d)

9.3  Program Results

Program Results (f11jcfe.r)


F11JCF (PDF version)
F11 Chapter Contents
F11 Chapter Introduction
NAG Library Manual

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