|
NAG[f11dec] NAG[nag_sparse_nsym_sol] - Solver with no Jacobi/SSOR preconditioning (nonsymmetric)
|
|
Calling Sequence
f11dec(method, precon, a, irow, icol, omega, b, m, tol, maxitn, x, rnorm, itn, comm, 'n'=n, 'nnz'=nnz, 'fail'=fail)
nag_sparse_nsym_sol(. . .)
Parameters
|
method - String;
|
|
|
On entry: specifies the iterative method to be used.
|
|
The restarted generalized minimum residual method is used.
|
|
The conjugate gradient squared method is used.
|
|
The bi-conjugate gradient stabilised method is used.
|
|
Constraint: "Nag_SparseNsym_RGMRES", "Nag_SparseNsym_CGS" or "Nag_SparseNsym_BiCGSTAB". .
|
|
|
precon - String;
|
|
|
On entry: specifies the type of preconditioning to be used. The possible choices are :
|
|
if , no preconditioning;
|
|
if , symmetric successive-over-relaxation;
|
|
if , Jacobi.
|
|
Constraint: "Nag_SparseNsym_NoPrec", "Nag_SparseNsym_SSORPrec" or "Nag_SparseNsym_JacPrec". .
|
|
|
a - Vector(1..nnz, datatype=float[8]);
|
|
|
On entry: the non-zero 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 f11zac (nag_sparse_nsym_sort) may be used to order the elements in this way.
|
|
|
irow - Vector(1..nnz, datatype=integer[kernelopts('wordsize')/8]);
icol - Vector(1..nnz, datatype=integer[kernelopts('wordsize')/8]);
|
|
|
On entry: the row and column indices of the non-zero elements supplied in a.
|
|
|
omega - float;
|
|
|
On entry: if , omega is the relaxation argument to be used in the SSOR method. Otherwise omega need not be initialized and is not referenced.
|
|
Constraint: . .
|
|
|
b - Vector(1..n, datatype=float[8]);
|
|
|
On entry: the right-hand side vector .
|
|
|
m - integer;
|
|
|
On entry: 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.
|
|
if , ;
|
|
if , .
|
|
|
tol - float;
|
|
|
Constraint: . .
|
|
|
maxitn - integer;
|
|
|
On entry: the maximum number of iterations allowed.
|
|
Constraint: . .
|
|
|
x - Vector(1..n, datatype=float[8]);
|
|
|
On entry: an initial approximation to the solution vector .
|
|
On exit: an improved approximation to the solution vector .
|
|
|
rnorm - assignable;
|
|
|
Note: On exit the variable rnorm will have a value of type float.
|
|
On exit: the final value of the residual norm , where is the output value of itn.
|
|
|
itn - assignable;
|
|
|
Note: On exit the variable itn will have a value of type integer.
|
|
On exit: the number of iterations carried out.
|
|
|
comm - Vector;
|
|
|
A Maple Vector, which should be generated using NAG[Nag_Sparse_Comm], corresponding to the Nag_Sparse_Comm structure.
|
|
|
'n'=n - integer; (optional)
|
|
|
Default value: the first dimension of the arrays b, x.
|
|
On entry: the order of the matrix .
|
|
Constraint: . .
|
|
|
'nnz'=nnz - integer; (optional)
|
|
|
Default value: the first dimension of the arrays a, irow, icol.
|
|
On entry: the number of non-zero elements in the matrix .
|
|
Constraint: . .
|
|
|
'fail'=fail - table; (optional)
|
|
|
The NAG error argument, see the documentation for NagError.
|
|
|
|
Description
|
|
|
Purpose
|
|
nag_sparse_nsym_sol (f11dec) 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), or stabilized bi-conjugate gradient (Bi-CGSTAB) method, without preconditioning, with Jacobi, or with SSOR preconditioning.
|
|
Description
|
|
nag_sparse_nsym_sol (f11dec) solves a real sparse nonsymmetric system of linear equations:
using an RGMRES (see Saad and Schultz (1986)), CGS (see Sonneveld (1989)), or Bi-CGSTAB method (see Van der Vorst (1989), Sleijpen and Fokkema (1993)).
The function allows the following choices for the preconditioner:
|
Jacobi preconditioning (see Young (1971));
|
|
symmetric successive-over-relaxation (SSOR) preconditioning (see Young (1971)).
|
For incomplete (ILU) preconditioning see f11dcc (nag_sparse_nsym_fac_sol).
The matrix is represented in coordinate storage (CS) format (see the f11 Chapter Introduction) in the arrays a, irow and icol. The array a holds the non-zero entries in the matrix, while irow and icol hold the corresponding row and column indices.
|
|
Error Indicators and Warnings
|
|
"NE_ACC_LIMIT"
The required accuracy could not be obtained. However, a reasonable accuracy has been obtained and further iterations cannot 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.
"NE_ALLOC_FAIL"
Dynamic memory allocation failed.
"NE_BAD_PARAM"
On entry, argument method had an illegal value.
"NE_INT_2"
On entry, , . Constraint:
"NE_INT_ARG_LT"
On entry, n must not be less than 1: .
"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 consult NAG for assistance.
"NE_NONSYMM_MATRIX_DUP"
A non-zero matrix element has been supplied which does not lie within the matrix , is out of order or has duplicate row and column indices, i.e., one or more of the following constraints has been violated:
|
, or
|
Call f11zac (nag_sparse_nsym_sort) to reorder and sum or remove duplicates.
"NE_NOT_REQ_ACC"
The required accuracy has not been obtained in maxitn iterations.
"NE_REAL"
On entry, . Constraint: when .
"NE_REAL_ARG_GE"
On entry, tol must not be greater than or equal to 1: .
"NE_ZERO_DIAGONAL_ELEM"
On entry, the matrix a has a zero diagonal element. Jacobi and SSOR preconditioners are not appropriate for this problem.
|
|
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_nsym_sol (f11dec) for each iteration is roughly proportional to nnz.
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 matrix of the coefficients .
|
|
|
Examples
|
|
>
|
method := "Nag_SparseNsym_RGMRES":
precon := "Nag_SparseNsym_SSORPrec":
n := 5:
nnz := 16:
omega := 1.05:
m := 1:
tol := 1e-10:
maxitn := 1000:
comm := NAG:-Nag_Sparse_Comm():
a := Vector([2, 1, -1, -3, -2, 1, 1, 5, 3, 1, -2, -3, -1, 4, -2, -6], datatype=float[8]):
irow := Vector([1, 1, 1, 2, 2, 2, 3, 3, 3, 3, 4, 4, 4, 5, 5, 5], datatype=integer[kernelopts('wordsize')/8]):
icol := Vector([1, 2, 4, 2, 3, 5, 1, 3, 4, 5, 1, 4, 5, 2, 3, 5], datatype=integer[kernelopts('wordsize')/8]):
b := Vector([0, -7, 33, -19, -28], datatype=float[8]):
x := Vector([0, 0, 0, 0, 0], datatype=float[8]):
NAG:-f11dec(method, precon, a, irow, icol, omega, b, m, tol, maxitn, x, rnorm, itn, comm, 'n' = n, 'nnz' = nnz):
|
|
|
See Also
|
|
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
f11 Chapter Introduction.
NAG Toolbox Overview.
NAG Web Site.
|
|