NAG CL Interface
f11dcc (real_​gen_​solve_​ilu)

Settings help

CL Name Style:


1 Purpose

f11dcc 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, with incomplete LU preconditioning.

2 Specification

#include <nag.h>
void  f11dcc (Nag_SparseNsym_Method method, Integer n, Integer nnz, const double a[], Integer la, const Integer irow[], const Integer icol[], const Integer ipivp[], const Integer ipivq[], const Integer istr[], const Integer idiag[], 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: f11dcc, nag_sparse_real_gen_solve_ilu or nag_sparse_nsym_fac_sol.

3 Description

f11dcc solves a real sparse nonsymmetric linear system of equations:
Ax = b ,  
using a preconditioned RGMRES (see Saad and Schultz (1986)), CGS (see Sonneveld (1989)), or Bi-CGSTAB () method (see Van der Vorst (1989), Sleijpen and Fokkema (1993)).
f11dcc uses the incomplete LU factorization determined by f11dac as the preconditioning matrix. A call to f11dcc must always be preceded by a call to f11dac. Alternative preconditioners for the same storage scheme are available by calling f11dec.
The matrix A , and the preconditioning matrix M , are represented in coordinate storage (CS) format (see the F11 Chapter Introduction) in the arrays a, irow and icol, as returned from f11dac. The array a holds the nonzero entries in these matrices, 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
Salvini S A and Shaw G J (1996) An evaluation of new NAG Library solvers for large sparse unsymmetric linear systems NAG Technical Report TR2/96
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

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
Then the bi-conjugate gradient stabilised () method is used.
Constraint: method=Nag_SparseNsym_RGMRES, Nag_SparseNsym_CGS or Nag_SparseNsym_BiCGSTAB.
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 f11dac.
Constraint: n1 .
3: nnz Integer Input
On entry: the number of nonzero-elements in the matrix A . This must be the same value as was supplied in the preceding call to f11dac.
Constraint: 1 nnz n 2 .
4: a[la] const double Input
On entry: the values returned in the array a by a previous call to f11dac.
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 f11dac.
Constraint: la 2 × nnz .
6: irow[la] const Integer Input
7: icol[la] const Integer Input
8: ipivp[n] const Integer Input
9: ipivq[n] const Integer Input
10: istr[n+1] const Integer Input
11: idiag[n] const Integer Input
On entry: the values returned in the arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to f11dac.
12: b[n] const double Input
On entry: the right-hand side vector b .
13: 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) .
14: 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 .
15: maxitn Integer Input
On entry: the maximum number of iterations allowed.
Constraint: maxitn1 .
16: 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 .
17: rnorm double * Output
On exit: the final value of the residual norm r k , where k is the output value of itn.
18: itn Integer * Output
On exit: the number of iterations carried out.
19: 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.
20: 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.
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_ALG_FAIL
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
NE_ALLOC_FAIL
Dynamic memory allocation failed.
NE_BAD_PARAM
On entry, argument method 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_INVALID_CS
On entry, the CS representation of A is invalid. Check that the call to f11dcc has been preceded by a valid call to f11dac, and that the arrays a, irow and icol have not been corrupted between the two calls.
NE_INVALID_CS_PRECOND
On entry, the CS representation of the preconditioning matrix M is invalid. Check that the call to f11dcc has been preceded by a valid call to f11dac, and that the arrays a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between the two calls.
NE_NOT_REQ_ACC
The required accuracy has not been obtained in maxitn iterations.
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

f11dcc is not threaded in any implementation.

9 Further Comments

The time taken by f11dcc for each iteration is roughly proportional to the value of nnzc returned from the preceding call to f11dac.
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 coefficient matrix, A ¯ = M −1 A .
Some illustrations of the application of f11dcc to linear systems arising from the discretization of two-dimensional elliptic partial differential equations, and to random-valued randomly structured linear systems, can be found in Salvini and Shaw (1996).

10 Example

This example program solves a sparse linear system of equations using the CGS method, with incomplete LU preconditioning.

10.1 Program Text

Program Text (f11dcce.c)

10.2 Program Data

Program Data (f11dcce.d)

10.3 Program Results

Program Results (f11dcce.r)