hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_sparse_complex_gen_solve_jacssor (f11ds)

Purpose

nag_sparse_complex_gen_solve_jacssor (f11ds) 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.

Syntax

[x, rnorm, itn, ifail] = f11ds(method, precon, a, irow, icol, omega, b, m, tol, maxitn, x, 'n', n, 'nnz', nnz)
[x, rnorm, itn, ifail] = nag_sparse_complex_gen_solve_jacssor(method, precon, a, irow, icol, omega, b, m, tol, maxitn, x, 'n', n, 'nnz', nnz)

Description

nag_sparse_complex_gen_solve_jacssor (f11ds) solves a complex sparse non-Hermitian system of linear equations:
Ax = b,
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_complex_gen_solve_jacssor (f11ds) allows the following choices for the preconditioner:
For incomplete LULU (ILU) preconditioning see nag_sparse_complex_gen_solve_ilu (f11dq).
The matrix AA is represented in coordinate storage (CS) format (see Section [Coordinate storage (CS) format] 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_complex_gen_solve_jacssor (f11ds) is a Black Box function which calls nag_sparse_complex_gen_basic_setup (f11br), nag_sparse_complex_gen_basic_solver (f11bs) and nag_sparse_complex_gen_basic_diag (f11bt). 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.

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

Parameters

Compulsory Input Parameters

1:     method – string
Specifies the iterative method to be used.
method = 'RGMRES'method='RGMRES'
Restarted generalized minimum residual method.
method = 'CGS'method='CGS'
Conjugate gradient squared method.
method = 'BICGSTAB'method='BICGSTAB'
Bi-conjugate gradient stabilized () method.
method = 'TFQMR'method='TFQMR'
Transpose-free quasi-minimal residual method.
Constraint: method = 'RGMRES'method='RGMRES', 'CGS''CGS', 'BICGSTAB''BICGSTAB' or 'TFQMR''TFQMR'.
2:     precon – string (length ≥ 1)
Specifies the type of preconditioning to be used.
precon = 'N'precon='N'
No preconditioning.
precon = 'J'precon='J'
Jacobi.
precon = 'S'precon='S'
Symmetric successive-over-relaxation (SSOR).
Constraint: precon = 'N'precon='N', 'J''J' or 'S''S'.
3:     a(nnz) – complex array
nnz, the dimension of the array, must satisfy the constraint 1nnzn21nnzn2.
The nonzero elements of the matrix AA, 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_complex_gen_sort (f11zn) may be used to order the elements in this way.
4:     irow(nnz) – int64int32nag_int array
5:     icol(nnz) – int64int32nag_int array
nnz, the dimension of the array, must satisfy the constraint 1nnzn21nnzn2.
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_complex_gen_sort (f11zn)):
  • 1irow(i)n1irowin and 1icol(i)n1icolin, for i = 1,2,,nnzi=1,2,,nnz;
  • either irow(i1) < irow(i)irowi-1<irowi or both irow(i1) = irow(i)irowi-1=irowi and icol(i1) < icol(i)icoli-1<icoli, for i = 2,3,,nnzi=2,3,,nnz.
6:     omega – double scalar
If precon = 'S'precon='S', 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.00.0<omega<2.0.
7:     b(n) – complex array
n, the dimension of the array, must satisfy the constraint n1n1.
The right-hand side vector bb.
8:     m – int64int32nag_int scalar
If method = 'RGMRES'method='RGMRES', m is the dimension of the restart subspace.
If method = 'BICGSTAB'method='BICGSTAB', m is the order  of the polynomial Bi-CGSTAB method.
Otherwise, m is not referenced.
Constraints:
  • if method = 'RGMRES'method='RGMRES', 0 < mmin (n,50)0<mmin(n,50);
  • if method = 'BICGSTAB'method='BICGSTAB', 0 < mmin (n,10)0<mmin(n,10).
9:     tol – double scalar
The required tolerance. Let xkxk denote the approximate solution at iteration kk, and rkrk the corresponding residual. The algorithm is considered to have converged at iteration kk if
rkτ × (b + Axk).
rkτ×(b+Axk).
If tol0.0tol0.0, τ = max (sqrt(ε),sqrt(n)ε)τ=max(ε,nε) is used, where εε is the machine precision. Otherwise τ = max (tol,10ε,sqrt(n)ε)τ=max(tol,10ε,nε) is used.
Constraint: tol < 1.0tol<1.0.
10:   maxitn – int64int32nag_int scalar
The maximum number of iterations allowed.
Constraint: maxitn1maxitn1.
11:   x(n) – complex array
n, the dimension of the array, must satisfy the constraint n1n1.
An initial approximation to the solution vector xx.

Optional Input Parameters

1:     n – int64int32nag_int scalar
Default: The dimension of the arrays b, x. (An error is raised if these dimensions are not equal.)
nn, the order of the matrix AA.
Constraint: n1n1.
2:     nnz – int64int32nag_int scalar
Default: The dimension of the arrays a, irow, icol. (An error is raised if these dimensions are not equal.)
The number of nonzero elements in the matrix AA.
Constraint: 1nnzn21nnzn2.

Input Parameters Omitted from the MATLAB Interface

work lwork iwork

Output Parameters

1:     x(n) – complex array
An improved approximation to the solution vector xx.
2:     rnorm – double scalar
The final value of the residual norm rkrk, where kk is the output value of itn.
3:     itn – int64int32nag_int scalar
The number of iterations carried out.
4:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:

Cases prefixed with W are classified as warnings and do not generate an error of type NAG:error_n. See nag_issue_warnings.

  ifail = 1ifail=1
On entry,method'RGMRES'method'RGMRES', 'CGS''CGS', 'BICGSTAB''BICGSTAB' or 'TFQMR''TFQMR',
orprecon'N'precon'N', 'J''J' or 'S''S',
orn < 1n<1,
ornnz < 1nnz<1,
ornnz > n2nnz>n2,
orprecon = 'S'precon='S' and omega lies outside the interval (0.0,2.0)(0.0,2.0),
orm < 1m<1,
orm > min (n,50)m>min(n,50), when method = 'RGMRES'method='RGMRES',
orm > min (n,10)m>min(n,10), when method = 'BICGSTAB'method='BICGSTAB',
ortol1.0tol1.0,
ormaxitn < 1maxitn<1,
orlwork is too small.
  ifail = 2ifail=2
On entry, the arrays irow and icol fail to satisfy the following constraints:
  • 1irow(i)n1irowin and 1icol(i)n1icolin, for i = 1,2,,nnzi=1,2,,nnz;
  • irow(i1) < irow(i)irowi-1<irowi, or irow(i1) = irow(i)irowi-1=irowi and icol(i1) < icol(i)icoli-1<icoli, for i = 2,3,,nnzi=2,3,,nnz.
Therefore a nonzero element has been supplied which does not lie within the matrix AA, is out of order, or has duplicate row and column indices. Call nag_sparse_complex_gen_sort (f11zn) to reorder and sum or remove duplicates.
  ifail = 3ifail=3
On entry, the matrix AA has a zero diagonal element. Jacobi and SSOR preconditioners are therefore not appropriate for this problem.
W ifail = 4ifail=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 = 5ifail=5
Required accuracy not obtained in maxitn iterations.
  ifail = 6ifail=6
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
  ifail = 7ifail=7 (nag_sparse_complex_gen_basic_setup (f11br), nag_sparse_complex_gen_basic_solver (f11bs) or nag_sparse_complex_gen_basic_diag (f11bt))
A serious error has occurred in an internal call to one of the specified functions. Check all function calls and array sizes. Seek expert help.

Accuracy

On successful termination, the final residual rk = bAxkrk=b-Axk, where k = itnk=itn, satisfies the termination criterion
rkτ × (b + Axk).
rkτ×(b+Axk).
The value of the final residual norm is returned in rnorm.

Further Comments

The time taken by nag_sparse_complex_gen_solve_jacssor (f11ds) 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 = M1AA-=M-1A, for some preconditioning matrix MM.

Example

function nag_sparse_complex_gen_solve_jacssor_example
method = 'CGS';
precon = 'N';
a = [ 2 + 3i;
      1 - 1i;
      -1 + 0i;
      0 + 2i;
      -2 + 1i;
      1 + 0i;
      0 - 1i;
      5 + 4i;
      3 - 1i;
      1 + 0i;
      -2 + 2i;
      -3 + 1i;
      0 + 3i;
      4 - 2i;
      -2 + 0i;
      -6 + 1i];
irow = [int64(1); 1; 1; 2; 2; 2; 3; 3; 3; 3; 4; 4; 4; 5; 5; 5];
icol = [int64(1); 2; 4; 2; 3; 5; 1; 3; 4; 5; 1; 4; 5; 2; 3; 5];
omega = 1.05;
b = [ -3 + 3i;
      -11 + 5i;
      23 + 48i;
      -41 + 2i;
      -28 - 31i];
m = int64(1);
tol = 1e-10;
maxitn = int64(1000);
x = complex(zeros(5, 1));
[x, rnorm, itn, ifail] = nag_sparse_complex_gen_solve_jacssor(method, ...
                    precon, a, irow, icol, omega, b, m, tol, maxitn, x);

fprintf('\nConverged in %d iterations\n', itn);
fprintf('Final residual norm = %16.4e\n', rnorm);
fprintf('\n          X\n');
disp(x);
 

Converged in 5 iterations
Final residual norm =       1.0905e-10

          X
   1.0000 + 2.0000i
   2.0000 + 3.0000i
   3.0000 + 4.0000i
   4.0000 + 5.0000i
   5.0000 + 6.0000i


function f11ds_example
method = 'CGS';
precon = 'N';
a = [ 2 + 3i;
      1 - 1i;
      -1 + 0i;
      0 + 2i;
      -2 + 1i;
      1 + 0i;
      0 - 1i;
      5 + 4i;
      3 - 1i;
      1 + 0i;
      -2 + 2i;
      -3 + 1i;
      0 + 3i;
      4 - 2i;
      -2 + 0i;
      -6 + 1i];
irow = [int64(1); 1; 1; 2; 2; 2; 3; 3; 3; 3; 4; 4; 4; 5; 5; 5];
icol = [int64(1); 2; 4; 2; 3; 5; 1; 3; 4; 5; 1; 4; 5; 2; 3; 5];
omega = 1.05;
b = [ -3 + 3i;
      -11 + 5i;
      23 + 48i;
      -41 + 2i;
      -28 - 31i];
m = int64(1);
tol = 1e-10;
maxitn = int64(1000);
x = complex(zeros(5, 1));
[x, rnorm, itn, ifail] = ...
    f11ds(method, precon, a, irow, icol, omega, b, m, tol, maxitn, x);

fprintf('\nConverged in %d iterations\n', itn);
fprintf('Final residual norm = %16.4e\n', rnorm);
fprintf('\n          X\n');
disp(x);
 

Converged in 5 iterations
Final residual norm =       1.4052e-10

          X
   1.0000 + 2.0000i
   2.0000 + 3.0000i
   3.0000 + 4.0000i
   4.0000 + 5.0000i
   5.0000 + 6.0000i



PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013