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_ilu (f11dq)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparse_complex_gen_solve_ilu (f11dq) 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.

Syntax

[x, rnorm, itn, ifail] = f11dq(method, nz, a, irow, icol, ipivp, ipivq, istr, idiag, b, m, tol, maxitn, x, 'n', n, 'la', la)
[x, rnorm, itn, ifail] = nag_sparse_complex_gen_solve_ilu(method, nz, a, irow, icol, ipivp, ipivq, istr, idiag, b, m, tol, maxitn, x, 'n', n, 'la', la)

Description

nag_sparse_complex_gen_solve_ilu (f11dq) 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_complex_gen_solve_ilu (f11dq) uses the incomplete LU factorization determined by nag_sparse_complex_gen_precon_ilu (f11dn) as the preconditioning matrix. A call to nag_sparse_complex_gen_solve_ilu (f11dq) must always be preceded by a call to nag_sparse_complex_gen_precon_ilu (f11dn). Alternative preconditioners for the same storage scheme are available by calling nag_sparse_complex_gen_solve_jacssor (f11ds).
The matrix A, and the preconditioning matrix M, are represented in coordinate storage (CS) format (see Coordinate storage (CS) format in the F11 Chapter Introduction) in the arrays a, irow and icol, as returned from nag_sparse_complex_gen_precon_ilu (f11dn). The array a holds the nonzero entries in these matrices, while irow and icol hold the corresponding row and column indices.
nag_sparse_complex_gen_solve_ilu (f11dq) 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

Parameters

Compulsory Input Parameters

1:     method – string
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:     nz int64int32nag_int scalar
The number of nonzero elements in the matrix A. This must be the same value as was supplied in the preceding call to nag_sparse_complex_gen_precon_ilu (f11dn).
Constraint: 1nzn2.
3:     ala – complex array
The values returned in the array a by a previous call to nag_sparse_complex_gen_precon_ilu (f11dn).
4:     irowla int64int32nag_int array
5:     icolla int64int32nag_int array
6:     ipivpn int64int32nag_int array
7:     ipivqn int64int32nag_int array
8:     istrn+1 int64int32nag_int array
9:     idiagn int64int32nag_int array
The values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_complex_gen_precon_ilu (f11dn).
ipivp and ipivq are restored on exit.
10:   bn – complex array
The right-hand side vector b.
11:   m int64int32nag_int scalar
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.
12:   tol – double scalar
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.
13:   maxitn int64int32nag_int scalar
The maximum number of iterations allowed.
Constraint: maxitn1.
14:   xn – complex array
An initial approximation to the solution vector x.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the arrays ipivp, ipivq, idiag, b, x. (An error is raised if these dimensions are not equal.)
n, the order of the matrix A. This must be the same value as was supplied in the preceding call to nag_sparse_complex_gen_precon_ilu (f11dn).
Constraint: n1.
2:     la int64int32nag_int scalar
Default: the dimension of the arrays a, irow, icol. (An error is raised if these dimensions are not equal.)
The dimension of the arrays a, irow and icol. this must be the same value as was supplied in the preceding call to nag_sparse_complex_gen_precon_ilu (f11dn).
Constraint: la2×nz.

Output Parameters

1:     xn – complex array
An improved approximation to the solution vector x.
2:     rnorm – double scalar
The final value of the residual norm rk, where k is the output value of itn.
3:     itn int64int32nag_int scalar
The number of iterations carried out.
4:     ifail int64int32nag_int scalar
ifail=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=1
On entry,method'RGMRES','CGS','BICGSTAB', or 'TFQMR',
orn<1,
ornz<1,
ornz>n2,
orla<2×nz,
orm<1 and method='RGMRES' or method='BICGSTAB',
orm>minn,50, with method='RGMRES',
orm>minn,10, with method='BICGSTAB',
ortol1.0,
ormaxitn<1,
orlwork too small.
   ifail=2
On entry, the CS representation of A is invalid. Further details are given in the error message. Check that the call to nag_sparse_complex_gen_solve_ilu (f11dq) has been preceded by a valid call to nag_sparse_complex_gen_precon_ilu (f11dn), and that the arrays a, irow, and icol have not been corrupted between the two calls.
   ifail=3
On entry, the CS representation of the preconditioning matrix M is invalid. Further details are given in the error message. Check that the call to nag_sparse_complex_gen_solve_ilu (f11dq) has been preceded by a valid call to nag_sparse_complex_gen_precon_ilu (f11dn) and that the arrays a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between the two calls.
W  ifail=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=5
Required accuracy not obtained in maxitn iterations.
   ifail=6
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
   ifail=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.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

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.

Further Comments

The time taken by nag_sparse_complex_gen_solve_ilu (f11dq) for each iteration is roughly proportional to the value of nnzc returned from the preceding call to nag_sparse_complex_gen_precon_ilu (f11dn).
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.

Example

This example solves a complex sparse non-Hermitian linear system of equations using the CGS method, with incomplete LU preconditioning.
function f11dq_example


fprintf('f11dq example results\n\n');

% Define sparse matrix A and RHS B
nz   = int64(24);
n    = int64(8);
a    = complex(zeros(3*nz,1));
irow = zeros(3*nz,1,'int64');
icol = zeros(3*nz,1,'int64');
a(1:nz) = [ 2 + 1i, -1 + 1i,  1 - 3i,  4 + 7i, -3 + 0i,  2 + 4i, ...
           -7 - 5i,  2 + 1i,  3 + 2i, -4 + 2i,  0 + 1i,  5 - 3i, ...
           -1 + 2i,  8 + 6i, -3 - 4i, -6 - 2i,  5 - 2i,  2 + 0i, ...
            0 - 5i, -1 + 5i,  6 + 2i, -1 + 4i,  2 + 0i,  3 + 3i];
b       = [ 7 + 11i;
            1 + 24i;
          -13 - 18i;
          -10 +  3i;
           23 + 14i;
           17 -  7i;
           15 -  3i;
           -3 + 20i];
irow(1:nz) = int64(...
             [1; 1; 1; 2; 2; 2; 3; 3; 4; 4; 4; 4;
              5; 5; 5; 6; 6; 6; 7; 7; 7; 8; 8; 8]);
icol(1:nz) = int64(...
             [1; 4; 8; 1; 2; 5; 3; 6; 1; 3; 4; 7;
              2; 5; 7; 1; 3; 6; 3; 5; 7; 2; 6; 8]);

% ILU preconditioner
%     Input parameters
lfill = int64(0);
dtol  = 0;
milu  = 'No modification';
ipivp = zeros(n, 1, 'int64');
ipivq = zeros(n, 1, 'int64');
% Compute preconditioner
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
    f11dn(...
          nz, a, irow, icol, lfill, dtol, milu, ipivp, ipivq);

% Iterative Solution Setup
%      Input parameters
method = 'CGS';
m      = int64(0);
tol    = 1e-10;
maxitn = int64(100);
x      = complex(zeros(n,1));

% Solve Ax = b
[x, rnorm, itn, ifail] = ...
  f11dq( ...
         method, nz, a, irow, icol, ipivp, ipivq, ...
         istr, idiag, b, m, tol, maxitn, x);

fprintf('Converged in %d iterations\n', itn);
fprintf('Final residual norm = %16.3e\n\n', rnorm);
disp('Solution');
disp(x);


f11dq example results

Converged in 4 iterations
Final residual norm =        1.348e-11

Solution
   1.0000 + 1.0000i
   2.0000 - 1.0000i
   3.0000 + 1.0000i
   4.0000 - 1.0000i
   3.0000 - 1.0000i
   2.0000 + 1.0000i
   1.0000 - 1.0000i
  -0.0000 + 3.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–2015