hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_sparse_real_symm_solve_jacssor (f11je)

Purpose

nag_sparse_real_symm_solve_jacssor (f11je) solves a real sparse symmetric system of linear equations, represented in symmetric coordinate storage format, using a conjugate gradient or Lanczos method, without preconditioning, with Jacobi or with SSOR preconditioning.

Syntax

[x, rnorm, itn, ifail] = f11je(method, precon, a, irow, icol, omega, b, tol, maxitn, x, 'n', n, 'nnz', nnz)
[x, rnorm, itn, ifail] = nag_sparse_real_symm_solve_jacssor(method, precon, a, irow, icol, omega, b, tol, maxitn, x, 'n', n, 'nnz', nnz)

Description

nag_sparse_real_symm_solve_jacssor (f11je) solves a real sparse symmetric linear system of equations
Ax = b,
Ax=b,
using a preconditioned conjugate gradient method (see Barrett et al. (1994)), or a preconditioned Lanczos method based on the algorithm SYMMLQ (see Paige and Saunders (1975)). The conjugate gradient method is more efficient if AA is positive definite, but may fail to converge for indefinite matrices. In this case the Lanczos method should be used instead. For further details see Barrett et al. (1994).
The function allows the following choices for the preconditioner:
For incomplete Cholesky (IC) preconditioning see nag_sparse_real_symm_solve_ichol (f11jc).
The matrix AA is represented in symmetric coordinate storage (SCS) format (see Section [Symmetric coordinate storage (SCS) format] in the F11 Chapter Introduction) in the arrays a, irow and icol. The array a holds the nonzero entries in the lower triangular part of the matrix, while irow and icol hold the corresponding row and column indices.

References

Barrett R, Berry M, Chan T F, Demmel J, Donato J, Dongarra J, Eijkhout V, Pozo R, Romine C and Van der Vorst H (1994) Templates for the Solution of Linear Systems: Building Blocks for Iterative Methods SIAM, Philadelphia
Paige C C and Saunders M A (1975) Solution of sparse indefinite systems of linear equations SIAM J. Numer. Anal. 12 617–629
Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York

Parameters

Compulsory Input Parameters

1:     method – string
Specifies the iterative method to be used.
method = 'CG'method='CG'
Conjugate gradient method.
method = 'SYMMLQ'method='SYMMLQ'
Lanczos method (SYMMLQ).
Constraint: method = 'CG'method='CG' or 'SYMMLQ''SYMMLQ'.
2:     precon – string (length ≥ 1)
Specifies the type of preconditioning to be used.
precon = 'N'precon='N'
No preconditioning.
precon = 'J'precon='J'
Jacobi.
precon = 'S'precon='S'
Symmetric successive-over-relaxation (SSOR).
Constraint: precon = 'N'precon='N', 'J''J' or 'S''S'.
3:     a(nnz) – double array
nnz, the dimension of the array, must satisfy the constraint 1nnzn × (n + 1) / 21nnzn×(n+1)/2.
The nonzero elements of the lower triangular part of the matrix AA, 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_symm_sort (f11zb) may be used to order the elements in this way.
4:     irow(nnz) – int64int32nag_int array
5:     icol(nnz) – int64int32nag_int array
nnz, the dimension of the array, must satisfy the constraint 1nnzn × (n + 1) / 21nnzn×(n+1)/2.
The row and column indices of the nonzero elements supplied in array a.
Constraints:
irow and icol must satisfy these constraints (which may be imposed by a call to nag_sparse_real_symm_sort (f11zb)):
  • 1irow(i)n1irowin and 1icol(i)irow(i)1icoliirowi, for i = 1,2,,nnzi=1,2,,nnz;
  • irow(i1) < irow(i)irowi-1<irowi or irow(i1) = irow(i)irowi-1=irowi and icol(i1) < icol(i)icoli-1<icoli, for i = 2,3,,nnzi=2,3,,nnz.
6:     omega – double scalar
If precon = 'S'precon='S', omega is the relaxation parameter ωω to be used in the SSOR method. Otherwise omega need not be initialized.
Constraint: 0.0 < omega < 2.00.0<omega<2.0.
7:     b(n) – double array
n, the dimension of the array, must satisfy the constraint n1n1.
The right-hand side vector bb.
8:     tol – double scalar
The required tolerance. Let xkxk denote the approximate solution at iteration kk, and rkrk the corresponding residual. The algorithm is considered to have converged at iteration kk if
rkτ × (b + Axk).
rkτ×(b+Axk).
If tol0.0tol0.0, τ = max (sqrt(ε),sqrt(n)ε)τ=max(ε,nε) is used, where εε is the machine precision. Otherwise τ = max (tol,10ε,sqrt(n)ε)τ=max(tol,10ε,nε) is used.
Constraint: tol < 1.0tol<1.0.
9:     maxitn – int64int32nag_int scalar
The maximum number of iterations allowed.
Constraint: maxitn1maxitn1.
10:   x(n) – double array
n, the dimension of the array, must satisfy the constraint n1n1.
An initial approximation to the solution vector xx.

Optional Input Parameters

1:     n – int64int32nag_int scalar
Default: The dimension of the arrays b, x. (An error is raised if these dimensions are not equal.)
nn, the order of the matrix AA.
Constraint: n1n1.
2:     nnz – 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 lower triangular part of the matrix AA.
Constraint: 1nnzn × (n + 1) / 21nnzn×(n+1)/2.

Input Parameters Omitted from the MATLAB Interface

work lwork iwork

Output Parameters

1:     x(n) – double array
An improved approximation to the solution vector xx.
2:     rnorm – double scalar
The final value of the residual norm rkrk, where kk is the output value of itn.
3:     itn – int64int32nag_int scalar
The number of iterations carried out.
4:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:
  ifail = 1ifail=1
On entry,method'CG'method'CG' or 'SYMMLQ''SYMMLQ',
orprecon'N'precon'N', 'J''J' or 'S''S',
orn < 1n<1,
ornnz < 1nnz<1,
ornnz > n × (n + 1) / 2nnz>n×(n+1)/2,
oromega lies outside the interval (0.0,2.0)(0.0,2.0),
ortol1.0tol1.0,
ormaxitn < 1maxitn<1,
orlwork too small.
  ifail = 2ifail=2
On entry, the arrays irow and icol fail to satisfy the following constraints:
  • 1irow(i)n1irowin and 1icol(i)irow(i)1icoliirowi, for i = 1,2,,nnzi=1,2,,nnz;
  • irow(i1) < irow(i)irowi-1<irowi, or irow(i1) = irow(i)irowi-1=irowi and icol(i1) < icol(i)icoli-1<icoli, for i = 2,3,,nnzi=2,3,,nnz.
Therefore a nonzero element has been supplied which does not lie in the lower triangular part of AA, is out of order, or has duplicate row and column indices. Call nag_sparse_real_symm_sort (f11zb) to reorder and sum or remove duplicates.
  ifail = 3ifail=3
On entry, the matrix AA has a zero diagonal element. Jacobi and SSOR preconditioners are not appropriate for this problem.
  ifail = 4ifail=4
The required accuracy could not be obtained. However, a reasonable accuracy has been obtained and further iterations could not improve the result.
  ifail = 5ifail=5
Required accuracy not obtained in maxitn iterations.
  ifail = 6ifail=6
The preconditioner appears not to be positive definite.
  ifail = 7ifail=7
The matrix of the coefficients appears not to be positive definite (conjugate gradient method only).
  ifail = 8ifail=8 (nag_sparse_real_symm_basic_setup (f11gd), nag_sparse_real_symm_basic_solver (f11ge) or nag_sparse_real_symm_basic_diag (f11gf))
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.

Accuracy

On successful termination, the final residual rk = bAxkrk=b-Axk, where k = itnk=itn, satisfies the termination criterion
rkτ × (b + Axk).
rkτ×(b+Axk).
The value of the final residual norm is returned in rnorm.

Further Comments

The time taken by nag_sparse_real_symm_solve_jacssor (f11je) for each iteration is roughly proportional to nnz. One iteration with the Lanczos method (SYMMLQ) requires a slightly larger number of operations than one iteration with the conjugate gradient method.
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 A = M1AA-=M-1A.

Example

function nag_sparse_real_symm_solve_jacssor_example
method = 'CG';
precon = 'S';
a = [4;
     1;
     5;
     2;
     2;
     3;
     -1;
     1;
     4;
     1;
     -2;
     3;
     2;
     -1;
     -2;
     5];
irow = [int64(1);2;2;3;4;4;5;5;5;6;6;6;7;7;7;7];
icol = [int64(1);1;2;3;2;4;1;4;5;2;5;6;1;2;3;7];
omega = 1.1;
b = [15;
     18;
     -8;
     21;
     11;
     10;
     29];
tol = 1e-06;
maxitn = int64(100);
x = [0;
     0;
     0;
     0;
     0;
     0;
     0];
[xOut, rnorm, itn, ifail] = ...
    nag_sparse_real_symm_solve_jacssor(method, precon, a, irow, icol, omega, b, tol, maxitn, x)
 

xOut =

    1.0000
    2.0000
    3.0000
    4.0000
    5.0000
    6.0000
    7.0000


rnorm =

   5.0257e-06


itn =

                    6


ifail =

                    0


function f11je_example
method = 'CG';
precon = 'S';
a = [4;
     1;
     5;
     2;
     2;
     3;
     -1;
     1;
     4;
     1;
     -2;
     3;
     2;
     -1;
     -2;
     5];
irow = [int64(1);2;2;3;4;4;5;5;5;6;6;6;7;7;7;7];
icol = [int64(1);1;2;3;2;4;1;4;5;2;5;6;1;2;3;7];
omega = 1.1;
b = [15;
     18;
     -8;
     21;
     11;
     10;
     29];
tol = 1e-06;
maxitn = int64(100);
x = [0;
     0;
     0;
     0;
     0;
     0;
     0];
[xOut, rnorm, itn, ifail] = ...
    f11je(method, precon, a, irow, icol, omega, b, tol, maxitn, x)
 

xOut =

    1.0000
    2.0000
    3.0000
    4.0000
    5.0000
    6.0000
    7.0000


rnorm =

   5.0257e-06


itn =

                    6


ifail =

                    0



PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013