NAG Library Function Document
nag_sparse_nsym_fac_sol (f11dcc)
1 Purpose
nag_sparse_nsym_fac_sol (f11dcc) 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, with incomplete preconditioning.
2 Specification
#include <nag.h> |
#include <nagf11.h> |
void |
nag_sparse_nsym_fac_sol (Nag_SparseNsym_Method method,
Integer n,
Integer nnz,
const double a[],
Integer la,
const Integer irow[],
const Integer icol[],
const Integer ipivp[],
const Integer ipivq[],
const Integer istr[],
const Integer idiag[],
const double b[],
Integer m,
double tol,
Integer maxitn,
double x[],
double *rnorm,
Integer *itn,
Nag_Sparse_Comm *comm,
NagError *fail) |
|
3 Description
nag_sparse_nsym_fac_sol (f11dcc) solves a real sparse nonsymmetric linear system of equations:
using a preconditioned RGMRES (see
Saad and Schultz (1986)), CGS (see
Sonneveld (1989)), or Bi-CGSTAB
method (see
Van der Vorst (1989),
Sleijpen and Fokkema (1993)).
nag_sparse_nsym_fac_sol (f11dcc) uses the incomplete
factorization determined by
nag_sparse_nsym_fac (f11dac) as the preconditioning matrix. A call to nag_sparse_nsym_fac_sol (f11dcc) must always be preceded by a call to
nag_sparse_nsym_fac (f11dac). Alternative preconditioners for the same storage scheme are available by calling
nag_sparse_nsym_sol (f11dec).
The matrix
, and the preconditioning matrix
, are represented in coordinate storage (CS) format (see the
f11 Chapter Introduction) in the arrays
a,
irow and
icol, as returned from
nag_sparse_nsym_fac (f11dac). The array
a holds the nonzero entries in these matrices, while
irow and
icol hold the corresponding row and column indices.
4 References
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
5 Arguments
- 1:
– Nag_SparseNsym_MethodInput
-
On entry: specifies the iterative method to be used.
- The restarted generalized minimum residual method is used.
- The conjugate gradient squared method is used.
- Then the bi-conjugate gradient stabilised method is used.
Constraint:
, or .
- 2:
– IntegerInput
-
On entry: the order of the matrix
. This
must be the same value as was supplied in the preceding call to
nag_sparse_nsym_fac (f11dac).
Constraint:
.
- 3:
– IntegerInput
-
On entry: the number of nonzero-elements in the matrix
. This
must be the same value as was supplied in the preceding call to
nag_sparse_nsym_fac (f11dac).
Constraint:
.
- 4:
– const doubleInput
-
On entry: the values returned in the array
a by a previous call to
nag_sparse_nsym_fac (f11dac).
- 5:
– IntegerInput
-
On entry: the
second
dimension of the arrays
a,
irow and
icol.This must be the same value as returned by a previous call to
nag_sparse_nsym_fac (f11dac).
Constraint:
.
- 6:
– const IntegerInput
- 7:
– const IntegerInput
- 8:
– const IntegerInput
- 9:
– const IntegerInput
- 10:
– const IntegerInput
- 11:
– const IntegerInput
-
On entry: the values returned in the arrays
irow,
icol,
ipivp,
ipivq,
istr and
idiag by a previous call to
nag_sparse_nsym_fac (f11dac).
- 12:
– const doubleInput
-
On entry: the right-hand side vector .
- 13:
– IntegerInput
-
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.
Constraints:
- if , ;
- if , .
- 14:
– doubleInput
-
On entry: 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:
.
- 15:
– IntegerInput
-
On entry: the maximum number of iterations allowed.
Constraint:
.
- 16:
– doubleInput/Output
-
On entry: an initial approximation to the solution vector .
On exit: an improved approximation to the solution vector .
- 17:
– double *Output
-
On exit: the final value of the residual norm
, where
is the output value of
itn.
- 18:
– Integer *Output
-
On exit: the number of iterations carried out.
- 19:
– Nag_Sparse_Comm *Input/Output
-
On entry/exit: a pointer to a structure of type Nag_Sparse_Comm whose members are used by the iterative solver.
- 20:
– NagError *Input/Output
-
The NAG error argument (see
Section 2.7 in How to Use the NAG Library and its Documentation).
6 Error Indicators and Warnings
- NE_2_INT_ARG_LT
-
On entry, while . These arguments must satisfy .
- 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_ALG_FAIL
-
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
- NE_ALLOC_FAIL
-
Dynamic memory allocation failed.
- NE_BAD_PARAM
-
On entry, argument
method had an illegal value.
- NE_INT_2
-
On entry, , .
Constraint: when .
On entry, , .
Constraint: when .
On entry, , .
Constraint: .
- NE_INT_ARG_LT
-
On entry, .
Constraint: .
On entry, .
Constraint: .
- NE_INVALID_CS
-
On entry, the CS representation of
is invalid. Check that the call to nag_sparse_nsym_fac_sol (f11dcc) has been preceded by a valid call to
nag_sparse_nsym_fac (f11dac), and that the arrays
a,
irow and
icol have not been corrupted between the two calls.
- NE_INVALID_CS_PRECOND
-
On entry, the CS representation of the preconditioning matrix M is invalid. Check that the call to nag_sparse_nsym_fac_sol (f11dcc) has been preceded by a valid call to
nag_sparse_nsym_fac (f11dac), and that the arrays
a,
irow,
icol,
ipivp,
ipivq,
istr and
idiag have not been corrupted between the two calls.
- NE_NOT_REQ_ACC
-
The required accuracy has not been obtained in
maxitn iterations.
- NE_REAL_ARG_GE
-
On entry,
tol must not be greater than or equal to 1.0:
.
7 Accuracy
On successful termination, the final residual
, where
, satisfies the termination criterion
The value of the final residual norm is returned in
rnorm.
8 Parallelism and Performance
nag_sparse_nsym_fac_sol (f11dcc) is not threaded in any implementation.
The time taken by nag_sparse_nsym_fac_sol (f11dcc) for each iteration is roughly proportional to the value of
nnzc returned from the preceding call to
nag_sparse_nsym_fac (f11dac).
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_nsym_fac_sol (f11dcc) 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).
10 Example
This example program solves a sparse linear system of equations using the CGS method, with incomplete preconditioning.
10.1 Program Text
Program Text (f11dcce.c)
10.2 Program Data
Program Data (f11dcce.d)
10.3 Program Results
Program Results (f11dcce.r)