NAG CL Interface
f11dec (real_​gen_​solve_​jacssor)

Settings help

CL Name Style:


1 Purpose

f11dec solves a real sparse nonsymmetric system of linear equations, represented in coordinate storage format, using a restarted generalized minimal residual (RGMRES), conjugate gradient squared (CGS), or stabilized bi-conjugate gradient (Bi-CGSTAB) method, without preconditioning, with Jacobi, or with SSOR preconditioning.

2 Specification

#include <nag.h>
void  f11dec (Nag_SparseNsym_Method method, Nag_SparseNsym_PrecType precon, Integer n, Integer nnz, const double a[], const Integer irow[], const Integer icol[], double omega, const double b[], Integer m, double tol, Integer maxitn, double x[], double *rnorm, Integer *itn, Nag_Sparse_Comm *comm, NagError *fail)
The function may be called by the names: f11dec, nag_sparse_real_gen_solve_jacssor or nag_sparse_nsym_sol.

3 Description

f11dec solves a real sparse nonsymmetric system of linear equations:
Ax = b ,  
using an RGMRES (see Saad and Schultz (1986)), CGS (see Sonneveld (1989)), or Bi-CGSTAB () method (see Van der Vorst (1989), Sleijpen and Fokkema (1993)).
The function allows the following choices for the preconditioner:
For incomplete LU (ILU) preconditioning see f11dcc.
The matrix A is represented in coordinate storage (CS) format (see the F11 Chapter Introduction) in the arrays a, irow and icol. The array a holds the nonzero entries in the matrix, while irow and icol hold the corresponding row and column indices.

4 References

Saad Y and Schultz M (1986) GMRES: a generalized minimal residual algorithm for solving nonsymmetric linear systems SIAM J. Sci. Statist. Comput. 7 856–869
Sleijpen G L G and Fokkema D R (1993) BiCGSTAB() for linear equations involving matrices with complex spectrum ETNA 1 11–32
Sonneveld P (1989) CGS, a fast Lanczos-type solver for nonsymmetric linear systems SIAM J. Sci. Statist. Comput. 10 36–52
Van der Vorst H (1989) Bi-CGSTAB, a fast and smoothly converging variant of Bi-CG for the solution of nonsymmetric linear systems SIAM J. Sci. Statist. Comput. 13 631–644
Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York

5 Arguments

1: method Nag_SparseNsym_Method Input
On entry: specifies the iterative method to be used.
method=Nag_SparseNsym_RGMRES
The restarted generalized minimum residual method is used.
method=Nag_SparseNsym_CGS
The conjugate gradient squared method is used.
method=Nag_SparseNsym_BiCGSTAB
The bi-conjugate gradient stabilised () method is used.
Constraint: method=Nag_SparseNsym_RGMRES, Nag_SparseNsym_CGS or Nag_SparseNsym_BiCGSTAB.
2: precon Nag_SparseNsym_PrecType Input
On entry: specifies the type of preconditioning to be used.
precon=Nag_SparseNsym_NoPrec
No preconditioning.
precon=Nag_SparseNsym_SSORPrec
Symmetric successive-over-relaxation.
precon=Nag_SparseNsym_JacPrec
Jacobi.
Constraint: precon=Nag_SparseNsym_NoPrec, Nag_SparseNsym_SSORPrec or Nag_SparseNsym_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 matrix A .
Constraint: 1 nnz n 2 .
5: a[nnz] const double Input
On entry: the nonzero elements 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 f11zac 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 f11zac):;
  • 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.
8: omega double Input
On entry: if precon=Nag_SparseNsym_SSORPrec, omega is the relaxation argument ω to be used in the SSOR method. Otherwise omega need not be initialized and is not referenced.
Constraint: 0.0 < omega < 2.0 .
9: b[n] const double Input
On entry: the right-hand side vector b .
10: m Integer Input
On entry: if method=Nag_SparseNsym_RGMRES, m is the dimension of the restart subspace.
If method=Nag_SparseNsym_BiCGSTAB, m is the order of the polynomial Bi-CGSTAB method; otherwise m is not referenced.
Constraints:
  • if method=Nag_SparseNsym_RGMRES, 0 < m min(n,50) ;
  • if method=Nag_SparseNsym_BiCGSTAB, 0 < m min(n,10) .
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_ACC_LIMIT
The required accuracy could not be obtained. However, a reasonable accuracy has been obtained and further iterations cannot improve the result.
You should check the output value of rnorm for acceptability. This error code usually implies that your problem has been fully and satisfactorily solved to within or close to the accuracy available on your system. Further iterations are unlikely to improve on this situation.
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_INT_2
On entry, m=value , min(n,10) = value.
Constraint: 0 < m min(n,10) when method=Nag_SparseNsym_BiCGSTAB.
On entry, m=value , min(n,50) = value.
Constraint: 0 < m min(n,50) when method=Nag_SparseNsym_RGMRES.
On entry, nnz=value , n=value .
Constraint: 1 nnz n 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_NONSYMM_MATRIX_DUP
A nonzero matrix element has been supplied which does not lie within 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 f11zac to reorder and sum or remove duplicates.
NE_NOT_REQ_ACC
The required accuracy has not been obtained in maxitn iterations.
NE_REAL
On entry, omega=value .
Constraint: 0.0 < omega < 2.0 when precon=Nag_SparseNsym_SSORPrec.
NE_REAL_ARG_GE
On entry, tol must not be greater than or equal to 1: tol=value .
NE_ZERO_DIAGONAL_ELEM
On entry, 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

f11dec is not threaded in any implementation.

9 Further Comments

The time taken by f11dec for each iteration is roughly proportional to nnz.
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 sparse nonsymmetric system of equations using the RGMRES method, with SSOR preconditioning.

10.1 Program Text

Program Text (f11dece.c)

10.2 Program Data

Program Data (f11dece.d)

10.3 Program Results

Program Results (f11dece.r)