PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_sparse_real_gen_solve_jacssor (f11de)
Purpose
nag_sparse_real_gen_solve_jacssor (f11de) 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, without preconditioning, with Jacobi, or with SSOR preconditioning.
Syntax
[
x,
rnorm,
itn,
ifail] = f11de(
method,
precon,
a,
irow,
icol,
omega,
b,
m,
tol,
maxitn,
x, 'n',
n, 'nz',
nz)
[
x,
rnorm,
itn,
ifail] = nag_sparse_real_gen_solve_jacssor(
method,
precon,
a,
irow,
icol,
omega,
b,
m,
tol,
maxitn,
x, 'n',
n, 'nz',
nz)
Description
nag_sparse_real_gen_solve_jacssor (f11de) solves a real sparse nonsymmetric system of linear equations
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.
The function allows the following choices for the preconditioner:
- no preconditioning;
- Jacobi preconditioning (see Young (1971));
- symmetric successive-over-relaxation (SSOR) preconditioning (see Young (1971)).
For incomplete
(ILU) preconditioning see
nag_sparse_real_gen_solve_ilu (f11dc).
The matrix
is represented in coordinate storage (CS) format (see
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_real_gen_solve_jacssor (f11de) 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
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:
– string
-
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:
– string (length ≥ 1)
-
Specifies the type of preconditioning to be used.
- No preconditioning.
- Jacobi.
- Symmetric successive-over-relaxation.
Constraint:
, or .
- 3:
– double array
-
The nonzero elements of the matrix
, 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_real_gen_sort (f11za) may be used to order the elements in this way.
- 4:
– int64int32nag_int array
- 5:
– int64int32nag_int array
-
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_real_gen_sort (f11za)):
- and , for ;
- or and , for .
- 6:
– double scalar
-
If
,
omega is the relaxation parameter
to be used in the SSOR method. Otherwise
omega need not be initialized and is not referenced.
Constraint:
.
- 7:
– double array
-
The right-hand side vector .
- 8:
– 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 , .
- 9:
– 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:
.
- 10:
– int64int32nag_int scalar
-
The maximum number of iterations allowed.
Constraint:
.
- 11:
– double array
-
An initial approximation to the solution vector .
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the arrays
b,
x. (An error is raised if these dimensions are not equal.)
, the order of the matrix .
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 number of nonzero elements in the matrix .
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 | , |
or | and omega lies outside the interval , |
or | , |
or | , with , |
or | , with , |
or | , |
or | , |
or | lwork too small. |
-
-
On entry, the arrays
irow and
icol fail to satisfy the following constraints:
- and , for ;
- , or and , for .
Therefore a nonzero element has been supplied which does not lie within the matrix
, is out of order, or has duplicate row and column indices. Call
nag_sparse_real_gen_sort (f11za) to reorder and sum or remove duplicates.
-
-
On entry, the matrix has a zero diagonal element. Jacobi and SSOR preconditioners are not appropriate for this problem.
- 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_jacssor (f11de) for each iteration is roughly proportional to
nz.
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 .
Example
This example solves a sparse nonsymmetric system of equations using the RGMRES method, with SSOR preconditioning.
Open in the MATLAB editor:
f11de_example
function f11de_example
fprintf('f11de example results\n\n');
a = [ 2; 1;-1;-3;-2; 1; 1; 5; 3; 1;-2;-3;-1; 4;-2;-6];
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];
b = [0; -7; 33; -19; -28];
x = zeros(5, 1);
method = 'RGMRES';
precon = 'S';
m = int64(1);
omega = 1.05;
tol = 1e-10;
maxitn = int64(1000);
[x, rnorm, itn, ifail] = ...
f11de( ...
method, precon, a, irow, icol, omega, b, m, tol, maxitn, x);
fprintf('Converged in %d iterations\n', itn);
fprintf('Final residual norm = %16.1e\n\n', rnorm);
disp('Solution');
disp(x);
f11de example results
Converged in 13 iterations
Final residual norm = 5.1e-09
Solution
1.0000
2.0000
3.0000
4.0000
5.0000
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015