NAG CL Interface
f11jec (real_​symm_​solve_​jacssor)

Settings help

CL Name Style:


1 Purpose

f11jec solves a real sparse symmetric system of linear equations, represented in symmetric coordinate storage format, using a conjugate gradient or Lanczos method, without preconditioning, with Jacobi or with SSOR preconditioning.

2 Specification

#include <nag.h>
void  f11jec (Nag_SparseSym_Method method, Nag_SparseSym_PrecType precon, Integer n, Integer nnz, const double a[], const Integer irow[], const Integer icol[], double omega, 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: f11jec, nag_sparse_real_symm_solve_jacssor or nag_sparse_sym_sol.

3 Description

f11jec solves a real sparse symmetric linear system of equations:
Ax = b ,  
using a preconditioned conjugate gradient method (see Barrett et al. (1994)), 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).
The function allows the following choices for the preconditioner:
For incomplete Cholesky (IC) preconditioning see f11jcc.
The matrix A is represented in symmetric coordinate storage (SCS) format (see the F11 Chapter Introduction) in the arrays a, irow and icol. The array a holds the nonzero entries in the lower triangular part of the matrix, 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
Paige C C and Saunders M A (1975) Solution of sparse indefinite systems of linear equations SIAM J. Numer. Anal. 12 617–629
Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York

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: precon Nag_SparseSym_PrecType Input
On entry: specifies the type of preconditioning to be used.
precon=Nag_SparseSym_NoPrec
No preconditioning is used.
precon=Nag_SparseSym_SSORPrec
Symmetric successive-over-relaxation is used.
precon=Nag_SparseSym_JacPrec
Jacobi preconditioning is used.
Constraint: precon=Nag_SparseSym_NoPrec, Nag_SparseSym_SSORPrec or Nag_SparseSym_JacPrec.
3: n Integer Input
On entry: the order of the matrix A .
Constraint: n1 .
4: 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 .
5: a[nnz] const double Input
On entry: the nonzero elements of 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.
6: irow[nnz] const Integer Input
7: icol[nnz] const Integer Input
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] irow[i] , 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.
8: omega double Input
On entry: if precon=Nag_SparseSym_SSORPrec, omega is the relaxation argument ω to be used in the SSOR method. Otherwise omega need not be initialized.
Constraint: 0.0 omega 2.0 .
9: b[n] const double Input
On entry: the right-hand side vector b .
10: 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 .
11: maxitn Integer Input
On entry: the maximum number of iterations allowed.
Constraint: maxitn1 .
12: x[n] double Input/Output
On entry: an initial approximation of the solution vector x .
On exit: an improved approximation to the solution vector x .
13: rnorm double * Output
On exit: the final value of the residual norm r k , where k is the output value of itn.
14: itn Integer * Output
On exit: the number of iterations carried out.
15: 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.
16: 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_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.
On entry, argument precon had an illegal value.
NE_COEFF_NOT_POS_DEF
The matrix of coefficients appears not to be positive definite (conjugate gradient method only).
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_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
On entry, omega=value .
Constraint: 0.0 omega 2.0 .
NE_REAL_ARG_GE
On entry, tol must not be greater than or equal to 1.0: tol=value .
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] irow[i] , 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.
NE_ZERO_DIAGONAL_ELEM
The matrix A has a zero diagonal element. Jacobi and SSOR preconditioners are not appropriate for this problem.

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

Background information to multithreading can be found in the Multithreading documentation.
f11jec is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f11jec 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 f11jec for each iteration is roughly proportional to nnz. 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 .

10 Example

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

10.1 Program Text

Program Text (f11jece.c)

10.2 Program Data

Program Data (f11jece.d)

10.3 Program Results

Program Results (f11jece.r)