NAG Library Routine Document

f11dgf  (real_gen_solve_bdilu)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

f11dgf 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), stabilized bi-conjugate gradient (Bi-CGSTAB), or transpose-free quasi-minimal residual (TFQMR) method, with block Jacobi or additive Schwarz preconditioning.

2
Specification

Fortran Interface
Subroutine f11dgf ( method, n, nnz, a, la, irow, icol, nb, istb, indb, lindb, ipivp, ipivq, istr, idiag, b, m, tol, maxitn, x, rnorm, itn, work, lwork, ifail)
Integer, Intent (In):: n, nnz, la, irow(la), icol(la), nb, istb(nb+1), indb(lindb), lindb, istr(lindb+1), idiag(lindb), m, maxitn, lwork
Integer, Intent (Inout):: ipivp(lindb), ipivq(lindb), ifail
Integer, Intent (Out):: itn
Real (Kind=nag_wp), Intent (In):: a(la), b(n), tol
Real (Kind=nag_wp), Intent (Inout):: x(n)
Real (Kind=nag_wp), Intent (Out):: rnorm, work(lwork)
Character (*), Intent (In):: method
C Header Interface
#include nagmk26.h
void  f11dgf_ ( const char *method, const Integer *n, const Integer *nnz, const double a[], const Integer *la, const Integer irow[], const Integer icol[], const Integer *nb, const Integer istb[], const Integer indb[], const Integer *lindb, Integer ipivp[], Integer ipivq[], const Integer istr[], const Integer idiag[], const double b[], const Integer *m, const double *tol, const Integer *maxitn, double x[], double *rnorm, Integer *itn, double work[], const Integer *lwork, Integer *ifail, const Charlen length_method)

3
Description

f11dgf solves a real sparse nonsymmetric 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.
f11dgf uses the incomplete (possibly overlapping) block LU factorization determined by f11dff as the preconditioning matrix. A call to f11dgf must always be preceded by a call to f11dff. Alternative preconditioners for the same storage scheme are available by calling f11dcf or f11def.
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 f11dff. The array a holds the nonzero entries in these matrices, while irow and icol hold the corresponding row and column indices.
f11dgf is a Black Box routine which calls f11bdf, f11bef and f11bff. 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
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 – 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
3:     nnz – IntegerInput
4:     ala – Real (Kind=nag_wp) arrayInput
5:     la – IntegerInput
6:     irowla – Integer arrayInput
7:     icolla – Integer arrayInput
8:     nb – IntegerInput
9:     istbnb+1 – Integer arrayInput
10:   indblindb – Integer arrayInput
11:   lindb – IntegerInput
12:   ipivplindb – Integer arrayInput
13:   ipivqlindb – Integer arrayInput
14:   istrlindb+1 – Integer arrayInput
15:   idiaglindb – Integer arrayInput
On entry: the values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to f11dff.
The arrays istb, indb and a together with the the scalars n, nnz, la, nb and lindb must be the same values that were supplied in the preceding call to f11dff.
16:   bn – Real (Kind=nag_wp) arrayInput
On entry: the right-hand side vector b.
17:   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.
18:   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ε,nε is used, where ε is the machine precision. Otherwise τ=maxtol,10ε,nε is used.
Constraint: tol<1.0.
19:   maxitn – IntegerInput
On entry: the maximum number of iterations allowed.
Constraint: maxitn1.
20:   xn – Real (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.
21:   rnorm – Real (Kind=nag_wp)Output
On exit: the final value of the residual norm rk, where k is the output value of itn.
22:   itn – IntegerOutput
On exit: the number of iterations carried out.
23:   worklwork – Real (Kind=nag_wp) arrayWorkspace
24:   lwork – IntegerInput
On entry: the dimension of the array work as declared in the (sub)program from which f11dgf is called.
Constraints:
  • if method='RGMRES', lwork6×n+m×m+n+5+101;
  • if method='CGS', lwork10×n+100;
  • if method='BICGSTAB', lwork2×n×m+8×n+m×m+2+100;
  • if method='TFQMR', lwork13×n+100.
25:   ifail – IntegerInput/Output
On entry: ifail must be set to 0, -1​ or ​1. If you are unfamiliar with this argument you should refer to Section 3.4 in How to Use the NAG Library and its Documentation 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 argument, 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, for b=value, istbb+1=value and istbb=value.
Constraint: istbb+1>istbb, for b=1,2,,nb.
On entry, indbvalue=value and n=value.
Constraint: 1indbmn, for m=1,2,,istbnb+1-1 
On entry, istb1=value.
Constraint: istb11.
On entry, la=value and nnz=value.
Constraint: la2×nnz.
On entry, lindb=value, istbnb+1-1=value and nb=value.
Constraint: lindbistbnb+1-1.
On entry, lwork=value.
Constraint: lworkvalue.
On entry, m=value and n=value.
Constraint: if method='RGMRES', 1mminn,value.
If method='BICGSTAB', 1mminn,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, nb=value and n=value.
Constraint: 1nbn.
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, element value of a was out of order.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to f11dff and f11dgf.
On entry, icolvalue=value and n=value.
Constraint: 1icolin, for i=1,2,,nnz.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to f11dff and f11dgf.
On entry, irowvalue=value and n=value.
Constraint: 1irowin, for i=1,2,,nnz.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to f11dff and f11dgf.
On entry, location value of irow,icol was a duplicate.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to f11dff and f11dgf.
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 f11dff and f11dgf.
ifail=4
The required accuracy could not be obtained. However a reasonable accuracy may have been achieved. 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
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=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.9 in How to Use the NAG Library and its Documentation for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.8 in How to Use the NAG Library and its Documentation for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 3.7 in How to Use the NAG Library and its Documentation 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

f11dgf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f11dgf 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 f11dgf for each iteration is roughly proportional to the value of nnzc returned from the preceding call to f11dff.
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.
Some illustrations of the application of f11dgf 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 reads in a sparse matrix A and a vector b. It calls f11dff, with the array lfill=0 and the array dtol=0.0, to compute an overlapping incomplete LU factorization. This is then used as an additive Schwarz preconditioner on a call to f11dgf which uses the Bi-CGSTAB method to solve Ax=b.

10.1
Program Text

Program Text (f11dgfe.f90)

10.2
Program Data

Program Data (f11dgfe.d)

10.3
Program Results

Program Results (f11dgfe.r)

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