PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_sparse_real_gen_solve_ilu (f11dc)
Purpose
nag_sparse_real_gen_solve_ilu (f11dc) 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 incomplete preconditioning.
Syntax
[
x,
rnorm,
itn,
ifail] = f11dc(
method,
nz,
a,
irow,
icol,
ipivp,
ipivq,
istr,
idiag,
b,
m,
tol,
maxitn,
x, 'n',
n, 'la',
la)
[
x,
rnorm,
itn,
ifail] = nag_sparse_real_gen_solve_ilu(
method,
nz,
a,
irow,
icol,
ipivp,
ipivq,
istr,
idiag,
b,
m,
tol,
maxitn,
x, 'n',
n, 'la',
la)
Description
nag_sparse_real_gen_solve_ilu (f11dc) solves a real sparse nonsymmetric linear system of equations:
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_real_gen_solve_ilu (f11dc) uses the incomplete
factorization determined by
nag_sparse_real_gen_precon_ilu (f11da) as the preconditioning matrix. A call to
nag_sparse_real_gen_solve_ilu (f11dc) must always be preceded by a call to
nag_sparse_real_gen_precon_ilu (f11da). Alternative preconditioners for the same storage scheme are available by calling
nag_sparse_real_gen_solve_jacssor (f11de).
The matrix
, and the preconditioning matrix
, 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_real_gen_precon_ilu (f11da). The array
a holds the nonzero entries in these matrices, while
irow and
icol hold the corresponding row and column indices.
nag_sparse_real_gen_solve_ilu (f11dc) is a Black Box function which calls
nag_sparse_real_gen_basic_setup (f11bd),
nag_sparse_real_gen_basic_solver (f11be) and
nag_sparse_real_gen_basic_diag (f11bf). 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
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
Parameters
Compulsory Input Parameters
- 1:
– string
-
Specifies the iterative method to be used.
- Restarted generalized minimum residual method.
- Conjugate gradient squared method.
- Bi-conjugate gradient stabilized () method.
- Transpose-free quasi-minimal residual method.
Constraint:
, , or .
- 2:
– int64int32nag_int scalar
-
The number of nonzero elements in the matrix
. This
must be the same value as was supplied in the preceding call to
nag_sparse_real_gen_precon_ilu (f11da).
Constraint:
.
- 3:
– double array
-
The values returned in the array
a by a previous call to
nag_sparse_real_gen_precon_ilu (f11da).
- 4:
– int64int32nag_int array
- 5:
– int64int32nag_int array
- 6:
– int64int32nag_int array
- 7:
– int64int32nag_int array
- 8:
– int64int32nag_int array
- 9:
– int64int32nag_int array
-
The values returned in arrays
irow,
icol,
ipivp,
ipivq,
istr and
idiag by a previous call to
nag_sparse_real_gen_precon_ilu (f11da).
ipivp and
ipivq are restored on exit.
- 10:
– double array
-
The right-hand side vector .
- 11:
– int64int32nag_int scalar
-
If
,
m is the dimension of the restart subspace.
If
,
m is the order
of the polynomial Bi-CGSTAB method; otherwise,
m is not referenced.
Constraints:
- if , ;
- if , .
- 12:
– double scalar
-
The required tolerance. Let
denote the approximate solution at iteration
, and
the corresponding residual. The algorithm is considered to have converged at iteration
if
If
,
is used, where
is the
machine precision. Otherwise
is used.
Constraint:
.
- 13:
– int64int32nag_int scalar
-
The maximum number of iterations allowed.
Constraint:
.
- 14:
– double array
-
An initial approximation to the solution vector .
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the arrays
ipivp,
ipivq,
idiag,
b,
x. (An error is raised if these dimensions are not equal.)
, the order of the matrix
. This
must be the same value as was supplied in the preceding call to
nag_sparse_real_gen_precon_ilu (f11da).
Constraint:
.
- 2:
– 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_real_gen_precon_ilu (f11da).
Constraint:
.
Output Parameters
- 1:
– double array
-
An improved approximation to the solution vector .
- 2:
– double scalar
-
The final value of the residual norm
, where
is the output value of
itn.
- 3:
– int64int32nag_int scalar
-
The number of iterations carried out.
- 4:
– int64int32nag_int scalar
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.
-
-
On entry, | , or 'TFQMR', |
or | , |
or | , |
or | , |
or | , |
or | and or , |
or | , with , |
or | , with , |
or | , |
or | , |
or | lwork too small. |
-
-
On entry, the CS representation of
is invalid. Further details are given in the error message. Check that the call to
nag_sparse_real_gen_solve_ilu (f11dc) has been preceded by a valid call to
nag_sparse_real_gen_precon_ilu (f11da), and that the arrays
a,
irow, and
icol have not been corrupted between the two calls.
-
-
On entry, the CS representation of the preconditioning matrix
is invalid. Further details are given in the error message. Check that the call to
nag_sparse_real_gen_solve_ilu (f11dc) has been preceded by a valid call to
nag_sparse_real_gen_precon_ilu (f11da) and that the arrays
a,
irow,
icol,
ipivp,
ipivq,
istr and
idiag have not been corrupted between the two calls.
- W
-
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.
-
-
Required accuracy not obtained in
maxitn iterations.
-
-
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
- (nag_sparse_real_gen_basic_setup (f11bd), nag_sparse_real_gen_basic_solver (f11be) or nag_sparse_real_gen_basic_diag (f11bf))
-
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.
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
On successful termination, the final residual
, where
, satisfies the termination criterion
The value of the final residual norm is returned in
rnorm.
Further Comments
The time taken by
nag_sparse_real_gen_solve_ilu (f11dc) for each iteration is roughly proportional to the value of
nnzc returned from the preceding call to
nag_sparse_real_gen_precon_ilu (f11da).
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 .
Some illustrations of the application of
nag_sparse_real_gen_solve_ilu (f11dc) 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).
Example
This example solves a sparse linear system of equations using the CGS method, with incomplete preconditioning.
Open in the MATLAB editor:
f11dc_example
function f11dc_example
fprintf('f11dc example results\n\n');
nz = int64(24);
a = zeros(3*nz, 1);
irow = zeros(3*nz, 1, 'int64');
icol = irow;
a(1:nz) = [ 2; -1; 1; 4; -3; 2; -7; 2; 3; -4; 5; 5;
-1; 8; -3; -6; 5; 2; -5; -1; 6; -1; 2; 3];
irow(1:nz) = [ 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) = [ 1; 4; 8; 1; 2; 5; 3; 6; 1; 3; 4; 7;
2; 5; 7; 1; 3; 6; 3; 5; 7; 2; 6; 8];
m = int64(4);
method = 'CGS';
lfill = int64(0);
dtol = 0;
pstrat = 'C';
milu = 'N';
tol = 1e-12;
maxitn = int64(100);
ipivp = zeros(8, 1, 'int64');
ipivq = zeros(8, 1, 'int64');
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
f11da( ...
nz, a, irow, icol, lfill, dtol, milu, ipivp, ipivq);
b = [6; 8; -9; 46; 17; 21; 22; 34];
x = zeros(8, 1);
[x, rnorm, itn, ifail] = f11dc( ...
method, nz, a, irow, icol, ipivp, ipivq, ...
istr, idiag, b, m, tol, maxitn, x);
fprintf('Converged in %4d iterations\n', itn);
if rnorm>=tol
fprintf('Final residual norm = %15.2e\n\n', rnorm);
else
fprintf('Final residual norm < tolerance (%10.3e)\n\n',tol);
end
disp('solution');
disp(x);
f11dc example results
Converged in 4 iterations
Final residual norm < tolerance ( 1.000e-12)
solution
1.0000
2.0000
3.0000
4.0000
5.0000
6.0000
7.0000
8.0000
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015