NAG C Library Function Document
nag_sparse_nherm_sol (f11dsc)
1
Purpose
nag_sparse_nherm_sol (f11dsc) 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, without preconditioning, with Jacobi, or with SSOR preconditioning.
2
Specification
#include <nag.h> |
#include <nagf11.h> |
void |
nag_sparse_nherm_sol (Nag_SparseNsym_Method method,
Nag_SparseNsym_PrecType precon,
Integer n,
Integer nnz,
const Complex a[],
const Integer irow[],
const Integer icol[],
double omega,
const Complex b[],
Integer m,
double tol,
Integer maxitn,
Complex x[],
double *rnorm,
Integer *itn,
NagError *fail) |
|
3
Description
nag_sparse_nherm_sol (f11dsc) solves a complex sparse non-Hermitian 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.
nag_sparse_nherm_sol (f11dsc) 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_nherm_fac_sol (f11dqc).
The matrix
is represented in coordinate storage (CS) format (see
Section 2.1.1 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_nherm_sol (f11dsc) is a Black Box function which calls
nag_sparse_nherm_basic_setup (f11brc),
nag_sparse_nherm_basic_solver (f11bsc) and
nag_sparse_nherm_basic_diagnostic (f11btc). 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.
4
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
5
Arguments
- 1:
– Nag_SparseNsym_MethodInput
-
On entry: 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:
– Nag_SparseNsym_PrecTypeInput
-
On entry: specifies the type of preconditioning to be used.
- No preconditioning.
- Jacobi.
- Symmetric successive-over-relaxation (SSOR).
Constraint:
, or .
- 3:
– IntegerInput
-
On entry: , the order of the matrix .
Constraint:
.
- 4:
– IntegerInput
-
On entry: the number of nonzero elements in the matrix .
Constraint:
.
- 5:
– const ComplexInput
-
On entry: 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_nherm_sort (f11znc) may be used to order the elements in this way.
- 6:
– const IntegerInput
- 7:
– const IntegerInput
-
On entry: 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_nherm_sort (f11znc)):
- and , for ;
- either or both and , for .
- 8:
– doubleInput
-
On entry: if
,
omega is the relaxation parameter
to be used in the SSOR method. Otherwise
omega need not be initialized and is not referenced.
Constraint:
.
- 9:
– const ComplexInput
-
On entry: the right-hand side vector .
- 10:
– 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 , .
- 11:
– 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:
.
- 12:
– IntegerInput
-
On entry: the maximum number of iterations allowed.
Constraint:
.
- 13:
– ComplexInput/Output
-
On entry: an initial approximation to the solution vector .
On exit: an improved approximation to the solution vector .
- 14:
– double *Output
-
On exit: the final value of the residual norm
, where
is the output value of
itn.
- 15:
– Integer *Output
-
On exit: the number of iterations carried out.
- 16:
– NagError *Input/Output
-
The NAG error argument (see
Section 3.7 in How to Use the NAG Library and its Documentation).
6
Error Indicators and Warnings
- 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.
- A nonzero element has been supplied which does not lie in the matrix , is out of order, or has duplicate row and column indices. Consider calling nag_sparse_nherm_sort (f11znc) to reorder and sum or remove duplicates.
- Jacobi and SSOR preconditioners are not appropriate for this problem.
- NE_ACCURACY
-
The required accuracy could not be obtained. However, a reasonable accuracy may have been achieved.
- 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.
See
Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
- NE_BAD_PARAM
-
On entry, argument had an illegal value.
- NE_CONVERGENCE
-
The solution has not converged after iterations.
- NE_ENUM_INT_2
-
On entry, and .
Constraint: .
- NE_INT
-
On entry, .
Constraint:
On entry, .
Constraint: .
On entry, .
Constraint: .
- NE_INT_2
-
On entry, and .
Constraint: .
- 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 contact
NAG for assistance.
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact
NAG for assistance.
See
Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
- NE_INVALID_CS
-
On entry, , and .
Constraint: and .
On entry, , and .
Constraint: and .
- NE_NO_LICENCE
-
Your licence key may have expired or may not have been installed correctly.
See
Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.
- NE_NOT_STRICTLY_INCREASING
-
On entry, is out of order: .
On entry, the location () is a duplicate: .
- NE_REAL
-
On entry, .
Constraint:
On entry, .
Constraint: .
- NE_ZERO_DIAG_ELEM
-
The matrix has a zero diagonal entry in row .
The matrix has no diagonal entry in row .
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_nherm_sol (f11dsc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_sparse_nherm_sol (f11dsc) makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the
x06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this function. Please also consult the
Users' Note for your implementation for any additional implementation-specific information.
The time taken by
nag_sparse_nherm_sol (f11dsc) for each iteration is roughly proportional to
nnz.
The number of iterations required to achieve a prescribed accuracy cannot easily be determined a priori, as it can depend dramatically on the conditioning and spectrum of the preconditioned coefficient matrix , for some preconditioning matrix .
10
Example
This example solves a complex sparse non-Hermitian system of equations using the CGS method, with no preconditioning.
10.1
Program Text
Program Text (f11dsce.c)
10.2
Program Data
Program Data (f11dsce.d)
10.3
Program Results
Program Results (f11dsce.r)