NAG C Library Function Document

nag_sparse_nherm_sol (f11dsc)

1
Purpose

nag_sparse_nherm_sol (f11dsc) 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, without preconditioning, with Jacobi, or with SSOR preconditioning.

2
Specification

#include <nag.h>
#include <nagf11.h>
void  nag_sparse_nherm_sol (Nag_SparseNsym_Method method, Nag_SparseNsym_PrecType precon, Integer n, Integer nnz, const Complex a[], const Integer irow[], const Integer icol[], double omega, const Complex b[], Integer m, double tol, Integer maxitn, Complex x[], double *rnorm, Integer *itn, NagError *fail)

3
Description

nag_sparse_nherm_sol (f11dsc) solves a complex sparse non-Hermitian system of linear equations:
Ax=b,  
using an 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_sol (f11dsc) allows the following choices for the preconditioner:
For incomplete LU (ILU) preconditioning see nag_sparse_nherm_fac_sol (f11dqc).
The matrix A is represented in coordinate storage (CS) format (see Section 2.1.1 in the f11 Chapter Introduction) in the arrays a, irow and icol. The array a holds the nonzero entries in the matrix, while irow and icol hold the corresponding row and column indices.
nag_sparse_nherm_sol (f11dsc) 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
Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York

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:     precon Nag_SparseNsym_PrecTypeInput
On entry: specifies the type of preconditioning to be used.
precon=Nag_SparseNsym_NoPrec
No preconditioning.
precon=Nag_SparseNsym_JacPrec
Jacobi.
precon=Nag_SparseNsym_SSORPrec
Symmetric successive-over-relaxation (SSOR).
Constraint: precon=Nag_SparseNsym_NoPrec, Nag_SparseNsym_JacPrec or Nag_SparseNsym_SSORPrec.
3:     n IntegerInput
On entry: n, the order of the matrix A.
Constraint: n1.
4:     nnz IntegerInput
On entry: the number of nonzero elements in the matrix A.
Constraint: 1nnzn2.
5:     a[nnz] const ComplexInput
On entry: the nonzero elements of the matrix A, ordered by increasing row index, and by increasing column index within each row. Multiple entries for the same row and column indices are not permitted. The function nag_sparse_nherm_sort (f11znc) may be used to order the elements in this way.
6:     irow[nnz] const IntegerInput
7:     icol[nnz] const IntegerInput
On entry: the row and column indices of the nonzero elements supplied in a.
Constraints:
irow and icol must satisfy the following constraints (which may be imposed by a call to nag_sparse_nherm_sort (f11znc)):
  • 1irow[i]n and 1icol[i]n, for i=0,1,,nnz-1;
  • either irow[i-1]<irow[i] or both irow[i-1]=irow[i] and icol[i-1]<icol[i], for i=1,2,,nnz-1.
8:     omega doubleInput
On entry: if precon=Nag_SparseNsym_SSORPrec, omega is the relaxation parameter ω to be used in the SSOR method. Otherwise omega need not be initialized and is not referenced.
Constraint: 0.0<omega<2.0.
9:     b[n] const ComplexInput
On entry: the right-hand side vector b.
10:   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.
11:   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ε,10ε,nε is used, where ε is the machine precision. Otherwise τ=maxtol,10ε,nε is used.
Constraint: tol<1.0.
12:   maxitn IntegerInput
On entry: the maximum number of iterations allowed.
Constraint: maxitn1.
13:   x[n] ComplexInput/Output
On entry: an initial approximation to the solution vector x.
On exit: an improved approximation to the solution vector x.
14:   rnorm double *Output
On exit: the final value of the residual norm rk, where k is the output value of itn.
15:   itn Integer *Output
On exit: the number of iterations carried out.
16:   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

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.
A nonzero element has been supplied which does not lie in the matrix A, is out of order, or has duplicate row and column indices. Consider calling nag_sparse_nherm_sort (f11znc) to reorder and sum or remove duplicates.
Jacobi and SSOR preconditioners are not appropriate for this problem.
NE_ACCURACY
The required accuracy could not be obtained. However, a reasonable accuracy may have been achieved.
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_ENUM_INT_2
On entry, m=value and n=value.
Constraint: 0<mminn,value.
NE_INT
On entry, maxitn=value.
Constraint: maxitn1 
On entry, n=value.
Constraint: n1.
On entry, nnz=value.
Constraint: nnz1.
NE_INT_2
On entry, nnz=value and n=value.
Constraint: 1nnzn2.
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.
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, i=value, icol[i-1]=value and n=value.
Constraint: icol[i-1]1 and icol[i-1]n.
On entry, i=value, irow[i-1]=value and n=value.
Constraint: irow[i-1]1 and irow[i-1]n.
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, a[i-1] is out of order: i=value.
On entry, the location (irow[I-1],icol[I-1]) is a duplicate: I=value.
NE_REAL
On entry, omega=value.
Constraint: 0.0<omega<2.0 
On entry, tol=value.
Constraint: tol<1.0.
NE_ZERO_DIAG_ELEM
The matrix A has a zero diagonal entry in row value.
The matrix A has no diagonal entry in row value.

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_sol (f11dsc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_sparse_nherm_sol (f11dsc) 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_sol (f11dsc) for each iteration is roughly proportional to nnz.
The number of iterations required to achieve a prescribed accuracy cannot easily be determined a priori, as it can depend dramatically on the conditioning and spectrum of the preconditioned coefficient matrix A-=M-1A, for some preconditioning matrix M.

10
Example

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

10.1
Program Text

Program Text (f11dsce.c)

10.2
Program Data

Program Data (f11dsce.d)

10.3
Program Results

Program Results (f11dsce.r)