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

NAG Library Routine Document

F11DQF

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

F11DQF solves a complex sparse non-Hermitian system of linear equations, represented in coordinate storage format, using a restarted generalized minimal residual (RGMRES), conjugate gradient squared (CGS), stabilized bi-conjugate gradient (Bi-CGSTAB), or transpose-free quasi-minimal residual (TFQMR) method, with incomplete LU preconditioning.

2  Specification

SUBROUTINE F11DQF ( METHOD, N, NNZ, A, LA, IROW, ICOL, IPIVP, IPIVQ, ISTR, IDIAG, B, M, TOL, MAXITN, X, RNORM, ITN, WORK, LWORK, IFAIL)
INTEGER  N, NNZ, LA, IROW(LA), ICOL(LA), IPIVP(N), IPIVQ(N), ISTR(N+1), IDIAG(N), M, MAXITN, ITN, LWORK, IFAIL
REAL (KIND=nag_wp)  TOL, RNORM
COMPLEX (KIND=nag_wp)  A(LA), B(N), X(N), WORK(LWORK)
CHARACTER(*)  METHOD

3  Description

F11DQF solves a complex sparse non-Hermitian linear system of equations
Ax=b,  
using a preconditioned RGMRES (see Saad and Schultz (1986)), CGS (see Sonneveld (1989)), Bi-CGSTAB() (see Van der Vorst (1989) and Sleijpen and Fokkema (1993)), or TFQMR (see Freund and Nachtigal (1991) and Freund (1993)) method.
F11DQF uses the incomplete LU factorization determined by F11DNF as the preconditioning matrix. A call to F11DQF must always be preceded by a call to F11DNF. Alternative preconditioners for the same storage scheme are available by calling F11DSF.
The matrix A, and the preconditioning matrix M, are represented in coordinate storage (CS) format (see Section 2.1.1 in the F11 Chapter Introduction) in the arrays A, IROW and ICOL, as returned from F11DNF. The array A holds the nonzero entries in these matrices, while IROW and ICOL hold the corresponding row and column indices.
F11DQF is a Black Box routine which calls F11BRF, F11BSF and F11BTF. If you wish to use an alternative storage scheme, preconditioner, or termination criterion, or require additional diagnostic information, you should call these underlying routines directly.

4  References

Freund R W (1993) A transpose-free quasi-minimal residual algorithm for non-Hermitian linear systems SIAM J. Sci. Comput. 14 470–482
Freund R W and Nachtigal N (1991) QMR: a Quasi-Minimal Residual Method for Non-Hermitian Linear Systems Numer. Math. 60 315–339
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

5  Parameters

1:     METHOD – CHARACTER(*)Input
On entry: specifies the iterative method to be used.
METHOD='RGMRES'
Restarted generalized minimum residual method.
METHOD='CGS'
Conjugate gradient squared method.
METHOD='BICGSTAB'
Bi-conjugate gradient stabilized () method.
METHOD='TFQMR'
Transpose-free quasi-minimal residual method.
Constraint: METHOD='RGMRES', 'CGS', 'BICGSTAB' or 'TFQMR'.
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 F11DNF.
Constraint: N1.
3:     NNZ – INTEGERInput
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 F11DNF.
Constraint: 1NNZN2.
4:     ALA – COMPLEX (KIND=nag_wp) arrayInput
On entry: the values returned in the array A by a previous call to F11DNF.
5:     LA – INTEGERInput
On entry: the dimension of the arrays A, IROW and ICOL as declared in the (sub)program from which F11DQF is called. This must be the same value as was supplied in the preceding call to F11DNF.
Constraint: LA2×NNZ.
6:     IROWLA – INTEGER arrayInput
7:     ICOLLA – INTEGER arrayInput
8:     IPIVPN – INTEGER arrayInput
9:     IPIVQN – INTEGER arrayInput
10:   ISTRN+1 – INTEGER arrayInput
11:   IDIAGN – INTEGER arrayInput
On entry: the values returned in arrays IROW, ICOL, IPIVP, IPIVQ, ISTR and IDIAG by a previous call to F11DNF.
IPIVP and IPIVQ are restored on exit.
12:   BN – COMPLEX (KIND=nag_wp) arrayInput
On entry: the right-hand side vector b.
13:   M – INTEGERInput
On entry: if METHOD='RGMRES', M is the dimension of the restart subspace.
If METHOD='BICGSTAB', M is the order  of the polynomial Bi-CGSTAB method.
Otherwise, M is not referenced.
Constraints:
  • if METHOD='RGMRES', 0<MminN,50;
  • if METHOD='BICGSTAB', 0<MminN,10.
14:   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ε,10ε,nε is used, where ε is the machine precision. Otherwise τ=maxTOL,10ε,nε is used.
Constraint: TOL<1.0.
15:   MAXITN – INTEGERInput
On entry: the maximum number of iterations allowed.
Constraint: MAXITN1.
16:   XN – COMPLEX (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.
17:   RNORM – REAL (KIND=nag_wp)Output
On exit: the final value of the residual norm rk, where k is the output value of ITN.
18:   ITN – INTEGEROutput
On exit: the number of iterations carried out.
19:   WORKLWORK – COMPLEX (KIND=nag_wp) arrayWorkspace
20:   LWORK – INTEGERInput
On entry: the dimension of the array WORK as declared in the (sub)program from which F11DQF is called.
Constraints:
  • if METHOD='RGMRES', LWORK4×N+M×M+N+5+121;
  • if METHOD='CGS', LWORK8×N+120;
  • if METHOD='BICGSTAB', LWORK2×N×M+3+M×M+2+120;
  • if METHOD='TFQMR', LWORK11×N+120.
21:   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'RGMRES','CGS','BICGSTAB', or 'TFQMR',
orN<1,
orNNZ<1,
orNNZ>N2,
orLA<2×NNZ,
orM<1 and METHOD='RGMRES' or METHOD='BICGSTAB',
orM>minN,50, with METHOD='RGMRES',
orM>minN,10, with METHOD='BICGSTAB',
orTOL1.0,
orMAXITN<1,
orLWORK too small.
IFAIL=2
On entry, the CS representation of A is invalid. Further details are given in the error message. Check that the call to F11DQF has been preceded by a valid call to F11DNF, and that the arrays A, IROW, and ICOL have not been corrupted between the two calls.
IFAIL=3
On entry, the CS representation of the preconditioning matrix M is invalid. Further details are given in the error message. Check that the call to F11DQF has been preceded by a valid call to F11DNF and that the arrays A, IROW, ICOL, IPIVP, IPIVQ, ISTR and IDIAG have not been corrupted between the two calls.
IFAIL=4
The required accuracy could not be obtained. However, a reasonable accuracy may have been obtained, and further iterations could not 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.
IFAIL=5
Required accuracy not obtained in MAXITN iterations.
IFAIL=6
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
IFAIL=7 (F11BRF, F11BSF or F11BTF)
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.
IFAIL=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.8 in the Essential Introduction for further information.
IFAIL=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.7 in the Essential Introduction for further information.
IFAIL=-999
Dynamic memory allocation failed.
See Section 3.6 in the Essential Introduction for further information.

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  Parallelism and Performance

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

9  Further Comments

The time taken by F11DQF for each iteration is roughly proportional to the value of NNZC returned from the preceding call to F11DNF.
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-1A.

10  Example

This example solves a complex sparse non-Hermitian linear system of equations using the CGS method, with incomplete LU preconditioning.

10.1  Program Text

Program Text (f11dqfe.f90)

10.2  Program Data

Program Data (f11dqfe.d)

10.3  Program Results

Program Results (f11dqfe.r)


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

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