NAG CL Interface
f11jcc (real_​symm_​solve_​ichol)

Settings help

CL Name Style:


1 Purpose

f11jcc 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

#include <nag.h>
void  f11jcc (Nag_SparseSym_Method method, Integer n, Integer nnz, const double a[], Integer la, const Integer irow[], const Integer icol[], const Integer ipiv[], const Integer istr[], const double b[], double tol, Integer maxitn, double x[], double *rnorm, Integer *itn, Nag_Sparse_Comm *comm, NagError *fail)
The function may be called by the names: f11jcc, nag_sparse_real_symm_solve_ichol or nag_sparse_sym_chol_sol.

3 Description

f11jcc solves a real sparse symmetric linear system of equations:
Ax = b ,  
using a preconditioned conjugate gradient method (Meijerink and Van der Vorst (1977)), or a preconditioned Lanczos method based on the algorithm SYMMLQ (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).
f11jcc uses the incomplete Cholesky factorization determined by f11jac as the preconditioning matrix. A call to f11jcc must always be preceded by a call to f11jac. Alternative preconditioners for the same storage scheme are available by calling f11jec.
The matrix A , and the preconditioning matrix M , are represented in symmetric coordinate storage (SCS) format (see the F11 Chapter Introduction) in the arrays a, irow and icol, as returned from f11jac. 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 Arguments

1: method Nag_SparseSym_Method Input
On entry: specifies the iterative method to be used.
method=Nag_SparseSym_CG
The conjugate gradient method is used.
method=Nag_SparseSym_Lanczos
The Lanczos method, SYMMLQ is used.
Constraint: method=Nag_SparseSym_CG or Nag_SparseSym_Lanczos.
2: n Integer Input
On entry: the order of the matrix A . This must be the same value as was supplied in the preceding call to f11jac.
Constraint: n1 .
3: nnz Integer Input
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 f11jac.
Constraint: 1 nnz n × (n+1) / 2 .
4: a[la] const double Input
On entry: the values returned in array a by a previous call to f11jac.
5: la Integer Input
On entry: the second dimension of the arrays a, irow and icol.This must be the same value as returned by a previous call to f11jac.
Constraint: la 2 × nnz .
6: irow[la] const Integer Input
7: icol[la] const Integer Input
8: ipiv[n] const Integer Input
9: istr[n+1] const Integer Input
On entry: the values returned in the arrays irow, icol, ipiv and istr by a previous call to f11jac.
10: b[n] const double Input
On entry: the right-hand side vector b .
11: tol double Input
On entry: the required tolerance. Let x k denote the approximate solution at iteration k , and r k the corresponding residual. The algorithm is considered to have converged at iteration k if:
r k τ × ( b + A x k ) .  
If tol0.0 , τ = max( ε , n ε ) is used, where ε is the machine precision. Otherwise τ = max(tol,10ε, n ,ε) is used.
Constraint: tol<1.0 .
12: maxitn Integer Input
On entry: the maximum number of iterations allowed.
Constraint: maxitn1 .
13: x[n] double Input/Output
On entry: an initial approximation to the solution vector x .
On exit: an improved approximation to the solution vector x .
14: rnorm double * Output
On exit: the final value of the residual norm r k , where k is the output value of itn.
15: itn Integer * Output
On exit: the number of iterations carried out.
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_ACC_LIMIT
The required accuracy could not be obtained. However, a reasonable accuracy has been obtained and further iterations cannot improve the result.
NE_ALLOC_FAIL
Dynamic memory allocation failed.
NE_BAD_PARAM
On entry, argument method had an illegal value.
NE_COEFF_NOT_POS_DEF
The matrix of coefficients appears not to be positive definite.
NE_INT_2
On entry, nnz=value , n=value .
Constraint: 1 nnz n × (n+1) / 2 .
NE_INT_ARG_LT
On entry, maxitn=value.
Constraint: maxitn1.
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_SCS
The SCS representation of the matrix A is invalid. Check that the call to f11jcc has been preceded by a valid call to f11jac, and that the arrays a, irow and icol have not been corrupted between the two calls.
NE_INVALID_SCS_PRECOND
The SCS representation of the preconditioning matrix M is invalid. Check that the call to f11jcc has been preceded by a valid call to f11jac, and that the arrays a, irow, icol, ipiv and istr have not been corrupted between the two calls.
NE_NOT_REQ_ACC
The required accuracy has not been obtained in maxitn iterations.
NE_PRECOND_NOT_POS_DEF
The preconditioner appears not to be positive definite.
NE_REAL_ARG_GE
On entry, tol must not be greater than or equal to 1.0: tol=value .

7 Accuracy

On successful termination, the final residual r k = b - A x k , where k=itn , satisfies the termination criterion
r k τ × ( b + A x k ) .  
The value of the final residual norm is returned in rnorm.

8 Parallelism and Performance

f11jcc is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f11jcc 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 function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

The time taken by f11jcc for each iteration is roughly proportional to the value of nnzc returned from the preceding call to f11jac. 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 −1 A .
Some illustrations of the application of f11jcc 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).

10 Example

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

10.1 Program Text

Program Text (f11jcce.c)

10.2 Program Data

Program Data (f11jcce.d)

10.3 Program Results

Program Results (f11jcce.r)