NAG FL Interface
f11dqf (complex_​gen_​solve_​ilu)

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

Fortran Interface
Subroutine f11dqf ( method, n, nnz, a, la, irow, icol, ipivp, ipivq, istr, idiag, b, m, tol, maxitn, x, rnorm, itn, work, lwork, ifail)
Integer, Intent (In) :: n, nnz, la, irow(la), icol(la), istr(n+1), idiag(n), m, maxitn, lwork
Integer, Intent (Inout) :: ipivp(n), ipivq(n), ifail
Integer, Intent (Out) :: itn
Real (Kind=nag_wp), Intent (In) :: tol
Real (Kind=nag_wp), Intent (Out) :: rnorm
Complex (Kind=nag_wp), Intent (In) :: a(la), b(n)
Complex (Kind=nag_wp), Intent (Inout) :: x(n)
Complex (Kind=nag_wp), Intent (Out) :: work(lwork)
Character (*), Intent (In) :: method
C Header Interface
#include <nag.h>
void  f11dqf_ (const char *method, const Integer *n, const Integer *nnz, const Complex a[], const Integer *la, const Integer irow[], const Integer icol[], Integer ipivp[], Integer ipivq[], const Integer istr[], const Integer idiag[], const Complex b[], const Integer *m, const double *tol, const Integer *maxitn, Complex x[], double *rnorm, Integer *itn, Complex work[], const Integer *lwork, Integer *ifail, const Charlen length_method)
The routine may be called by the names f11dqf or nagf_sparse_complex_gen_solve_ilu.

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 Arguments

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 Integer Input
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 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 f11dnf.
Constraint: 1nnzn2.
4: ala Complex (Kind=nag_wp) array Input
On entry: the values returned in the array a by a previous call to f11dnf.
5: la Integer Input
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 array Input
7: icolla Integer array Input
8: ipivpn Integer array Input
9: ipivqn Integer array Input
10: istrn+1 Integer array Input
11: idiagn Integer array Input
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) array Input
On entry: the right-hand side vector b.
13: m Integer Input
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 Integer Input
On entry: the maximum number of iterations allowed.
Constraint: maxitn1.
16: xn Complex (Kind=nag_wp) array Input/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 Integer Output
On exit: the number of iterations carried out.
19: worklwork Complex (Kind=nag_wp) array Workspace
20: lwork Integer Input
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 Integer Input/Output
On entry: ifail must be set to 0, -1 or 1 to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of 0 causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of -1 means that an error message is printed while a value of 1 means that it is not.
If halting is not appropriate, the value -1 or 1 is recommended. If message printing is undesirable, then the value 1 is recommended. Otherwise, the value 0 is recommended. 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, la=value and nnz=value.
Constraint: la2×nnz.
On entry, lwork is too small: lwork=value. Minimum required value of lwork=value.
On entry, m=value and n=value.
Constraint: m1 and mminn,value.
On entry, maxitn=value.
Constraint: maxitn1.
On entry, method=value.
Constraint: method='RGMRES', 'CGS' or 'BICGSTAB'.
On entry, n=value.
Constraint: n1.
On entry, nnz=value.
Constraint: nnz1.
On entry, nnz=value and n=value.
Constraint: nnzn2.
On entry, tol=value.
Constraint: tol<1.0.
ifail=2
On entry, ai is out of order: i=value.
On entry, i=value, icoli=value, and n=value.
Constraint: icoli1 and icolin.
On entry, i=value, irowi=value, n=value.
Constraint: irowi1 and irowin.
On entry, the location (irowi,icoli) is a duplicate: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to f11dqf and f11dnf.
ifail=3
The CS representation of the preconditioner is invalid.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to f11dnf and f11dqf.
ifail=4
The required accuracy could not be obtained. However, a reasonable accuracy may have been achieved.
ifail=5
The solution has not converged after value iterations.
ifail=6
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
ifail=7
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.
A serious error, code value, has occurred in an internal call. 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 7 in the Introduction to the NAG Library FL Interface for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library FL Interface for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface 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)