NAG Library Function Document

nag_sparse_nherm_precon_bdilu_solve (f11duc)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

nag_sparse_nherm_precon_bdilu_solve (f11duc) 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 block Jacobi or additive Schwarz preconditioning.

2
Specification

#include <nag.h>
#include <nagf11.h>
void  nag_sparse_nherm_precon_bdilu_solve (Nag_SparseNsym_Method method, Integer n, Integer nnz, const Complex a[], Integer la, const Integer irow[], const Integer icol[], Integer nb, const Integer istb[], const Integer indb[], Integer lindb, const Integer ipivp[], const Integer ipivq[], const Integer istr[], const Integer idiag[], const Complex b[], Integer m, double tol, Integer maxitn, Complex x[], double *rnorm, Integer *itn, NagError *fail)

3
Description

nag_sparse_nherm_precon_bdilu_solve (f11duc) 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.
nag_sparse_nherm_precon_bdilu_solve (f11duc) uses the incomplete (possibly overlapping) block LU factorization determined by nag_sparse_nherm_precon_bdilu (f11dtc) as the preconditioning matrix. A call to nag_sparse_nherm_precon_bdilu_solve (f11duc) must always be preceded by a call to nag_sparse_nherm_precon_bdilu (f11dtc). Alternative preconditioners for the same storage scheme are available by calling nag_sparse_nherm_fac_sol (f11dqc) or nag_sparse_nherm_sol (f11dsc).
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 nag_sparse_nherm_precon_bdilu (f11dtc). The array a holds the nonzero entries in these matrices, while irow and icol hold the corresponding row and column indices.
nag_sparse_nherm_precon_bdilu_solve (f11duc) is a Black Box function which calls nag_sparse_nherm_basic_setup (f11brc), nag_sparse_nherm_basic_solver (f11bsc) and nag_sparse_nherm_basic_diagnostic (f11btc). If you wish to use an alternative storage scheme, preconditioner, or termination criterion, or require additional diagnostic information, you should call these underlying functions 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 Nag_SparseNsym_MethodInput
On entry: specifies the iterative method to be used.
method=Nag_SparseNsym_RGMRES
Restarted generalized minimum residual method.
method=Nag_SparseNsym_CGS
Conjugate gradient squared method.
method=Nag_SparseNsym_BiCGSTAB
Bi-conjugate gradient stabilized () method.
method=Nag_SparseNsym_TFQMR
Transpose-free quasi-minimal residual method.
Constraint: method=Nag_SparseNsym_RGMRES, Nag_SparseNsym_CGS, Nag_SparseNsym_BiCGSTAB or Nag_SparseNsym_TFQMR.
2:     n IntegerInput
3:     nnz IntegerInput
4:     a[la] const ComplexInput
5:     la IntegerInput
6:     irow[la] const IntegerInput
7:     icol[la] const IntegerInput
8:     nb IntegerInput
9:     istb[nb+1] const IntegerInput
10:   indb[lindb] const IntegerInput
11:   lindb IntegerInput
12:   ipivp[lindb] const IntegerInput
13:   ipivq[lindb] const IntegerInput
14:   istr[lindb+1] const IntegerInput
15:   idiag[lindb] const IntegerInput
On entry: the values returned in arrays a, irow, icol, ipivp, ipivq and idiag by a previous call to nag_sparse_nherm_precon_bdilu (f11dtc).
The arrays istb and indb together with the scalars n, nnz, la, nb and lindb must be the same values that were supplied in the preceding call to nag_sparse_nherm_precon_bdilu (f11dtc).
16:   b[n] const ComplexInput
On entry: the right-hand side vector b.
17:   m IntegerInput
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<mminn,50;
  • if method=Nag_SparseNsym_BiCGSTAB, 0<mminn,10.
18:   tol doubleInput
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:   x[n] ComplexInput/Output
On entry: an initial approximation to the solution vector x.
On exit: an improved approximation to the solution vector x.
21:   rnorm double *Output
On exit: the final value of the residual norm rk, where k is the output value of itn.
22:   itn Integer *Output
On exit: the number of iterations carried out.
23:   fail NagError *Input/Output
The NAG error argument (see Section 3.7 in How to Use the NAG Library and its Documentation).

6
Error Indicators and Warnings

NE_ACCURACY
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.
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.
See Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_CONVERGENCE
The solution has not converged after value iterations.
NE_INT
On entry, istb[0]=value.
Constraint: istb[0]1.
On entry, maxitn=value.
Constraint: maxitn1.
On entry, n=value.
Constraint: n1.
On entry, nnz=value.
Constraint: nnz1.
NE_INT_2
On entry, la=value and nnz=value.
Constraint: la2×nnz.
On entry, nb=value and n=value.
Constraint: 1nbn.
On entry, nnz=value and n=value.
Constraint: nnzn2.
NE_INT_3
On entry, lindb=value, istb[nb]-1=value and nb=value.
Constraint: lindbistb[nb]-1.
On entry, m=value and n=value.
Constraint: if method=Nag_SparseNsym_RGMRES, 1mminn,value.
If method=Nag_SparseNsym_BiCGSTAB, 1mminn,value.
NE_INT_ARRAY
On entry, indb[value]=value and n=value.
Constraint: 1indb[m-1]n, for m=istb[0],istb[0]+1,,istb[nb]-1.
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.
See Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
NE_INVALID_CS
On entry, icol[value]=value and n=value.
Constraint: 1icol[i-1]n, for i=1,2,,nnz.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_bdilu (f11dtc) and nag_sparse_nherm_precon_bdilu_solve (f11duc).
On entry, irow[value]=value and n=value.
Constraint: 1irow[i-1]n, for i=1,2,,nnz.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_bdilu (f11dtc) and nag_sparse_nherm_precon_bdilu_solve (f11duc).
NE_INVALID_CS_PRECOND
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 nag_sparse_nherm_precon_bdilu (f11dtc) and nag_sparse_nherm_precon_bdilu_solve (f11duc).
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.
NE_NOT_STRICTLY_INCREASING
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 nag_sparse_nherm_precon_bdilu (f11dtc) and nag_sparse_nherm_precon_bdilu_solve (f11duc).
On entry, for b=value, istb[b]=value and istb[b-1]=value.
Constraint: istb[b]>istb[b-1], for b=1,2,,nb.
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 nag_sparse_nherm_precon_bdilu (f11dtc) and nag_sparse_nherm_precon_bdilu_solve (f11duc).
NE_REAL
On entry, tol=value.
Constraint: tol<1.0.

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

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

9
Further Comments

The time taken by nag_sparse_nherm_precon_bdilu_solve (f11duc) for each iteration is roughly proportional to the value of nnzc returned from the preceding call to nag_sparse_nherm_precon_bdilu (f11dtc).
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 program reads in a sparse matrix A and a vector b. It calls nag_sparse_nherm_precon_bdilu (f11dtc), 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 nag_sparse_nherm_precon_bdilu_solve (f11duc) which uses the RGMRES method to solve Ax=b.

10.1
Program Text

Program Text (f11duce.c)

10.2
Program Data

Program Data (f11duce.d)

10.3
Program Results

Program Results (f11duce.r)

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