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_complex_herm_solve_ilu (f11jq)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparse_complex_herm_solve_ilu (f11jq) solves a complex sparse Hermitian system of linear equations, represented in symmetric coordinate storage format, using a conjugate gradient or Lanczos method, with incomplete Cholesky preconditioning.

Syntax

[x, rnorm, itn, ifail] = f11jq(method, nz, a, irow, icol, ipiv, istr, b, tol, maxitn, x, 'n', n, 'la', la)
[x, rnorm, itn, ifail] = nag_sparse_complex_herm_solve_ilu(method, nz, a, irow, icol, ipiv, istr, b, tol, maxitn, x, 'n', n, 'la', la)

Description

nag_sparse_complex_herm_solve_ilu (f11jq) solves a complex sparse Hermitian linear system of equations
Ax=b,  
using a preconditioned conjugate gradient method (see Meijerink and Van der Vorst (1977)), or a preconditioned Lanczos method based on the algorithm SYMMLQ (see Paige and Saunders (1975)). The conjugate gradient method is more efficient if A 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).
nag_sparse_complex_herm_solve_ilu (f11jq) uses the incomplete Cholesky factorization determined by nag_sparse_complex_herm_precon_ilu (f11jn) as the preconditioning matrix. A call to nag_sparse_complex_herm_solve_ilu (f11jq) must always be preceded by a call to nag_sparse_complex_herm_precon_ilu (f11jn). Alternative preconditioners for the same storage scheme are available by calling nag_sparse_complex_herm_solve_jacssor (f11js).
The matrix A and the preconditioning matrix M are represented in symmetric coordinate storage (SCS) format (see Symmetric coordinate storage (SCS) format in the F11 Chapter Introduction) in the arrays a, irow and icol, as returned from nag_sparse_complex_herm_precon_ilu (f11jn). The array a holds the nonzero entries in the lower triangular parts of these matrices, 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
Meijerink J and Van der Vorst H (1977) An iterative solution method for linear systems of which the coefficient matrix is a symmetric M-matrix Math. Comput. 31 148–162
Paige C C and Saunders M A (1975) Solution of sparse indefinite systems of linear equations SIAM J. Numer. Anal. 12 617–629

Parameters

Compulsory Input Parameters

1:     method – string
Specifies the iterative method to be used.
method='CG'
Conjugate gradient method.
method='SYMMLQ'
Lanczos method (SYMMLQ).
Constraint: method='CG' or 'SYMMLQ'.
2:     nz int64int32nag_int scalar
The number of nonzero elements in the lower triangular part of the matrix A. This must be the same value as was supplied in the preceding call to nag_sparse_complex_herm_precon_ilu (f11jn).
Constraint: 1nzn×n+1/2.
3:     ala – complex array
The values returned in the array a by a previous call to nag_sparse_complex_herm_precon_ilu (f11jn).
4:     irowla int64int32nag_int array
5:     icolla int64int32nag_int array
6:     ipivn int64int32nag_int array
7:     istrn+1 int64int32nag_int array
The values returned in arrays irow, icol, ipiv and istr by a previous call to nag_sparse_complex_herm_precon_ilu (f11jn).
8:     bn – complex array
The right-hand side vector b.
9:     tol – double scalar
The required tolerance. Let xk denote the approximate solution at iteration k, and rk the corresponding residual. The algorithm is considered to have converged at iteration k if
rkτ×b+Axk.  
If tol0.0, τ=maxε,10ε,nε is used, where ε is the machine precision. Otherwise τ=maxtol,10ε,nε is used.
Constraint: tol<1.0.
10:   maxitn int64int32nag_int scalar
The maximum number of iterations allowed.
Constraint: maxitn1.
11:   xn – complex array
An initial approximation to the solution vector x.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the arrays ipiv, b, x. (An error is raised if these dimensions are not equal.)
n, the order of the matrix A. This must be the same value as was supplied in the preceding call to nag_sparse_complex_herm_precon_ilu (f11jn).
Constraint: n1.
2:     la int64int32nag_int scalar
Default: the dimension of the arrays a, irow, icol. (An error is raised if these dimensions are not equal.)
The dimension of the arrays a, irow and icol. this must be the same value as was supplied in the preceding call to nag_sparse_complex_herm_precon_ilu (f11jn).
Constraint: la2×nz.

Output Parameters

1:     xn – complex array
An improved approximation to the solution vector x.
2:     rnorm – double scalar
The final value of the residual norm rk, where k is the output value of itn.
3:     itn int64int32nag_int scalar
The number of iterations carried out.
4:     ifail int64int32nag_int scalar
ifail=0 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.

   ifail=1
On entry,method'CG' or 'SYMMLQ',
orn<1,
ornz<1,
ornz>n×n+1/2,
orla too small,
ortol1.0,
ormaxitn<1,
orlwork too small.
   ifail=2
On entry, the SCS representation of A is invalid. Further details are given in the error message. Check that the call to nag_sparse_complex_herm_solve_ilu (f11jq) has been preceded by a valid call to nag_sparse_complex_herm_precon_ilu (f11jn), and that the arrays a, irow, and icol have not been corrupted between the two calls.
   ifail=3
On entry, the SCS representation of M is invalid. Further details are given in the error message. Check that the call to nag_sparse_complex_herm_solve_ilu (f11jq) has been preceded by a valid call to nag_sparse_complex_herm_precon_ilu (f11jn), and that the arrays a, irow, icol, ipiv and istr have not been corrupted between the two calls.
W  ifail=4
The required accuracy could not be obtained. However, a reasonable accuracy has been obtained and further iterations could not improve the result.
   ifail=5
Required accuracy not obtained in maxitn iterations.
   ifail=6
The preconditioner appears not to be positive definite.
   ifail=7
The matrix of the coefficients appears not to be positive definite (conjugate gradient method only).
   ifail=8
A serious error has occurred in an internal call to an auxiliary function. Check all function calls and array sizes. Seek expert help.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

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

Further Comments

The time taken by nag_sparse_complex_herm_solve_ilu (f11jq) for each iteration is roughly proportional to the value of nnzc returned from the preceding call to nag_sparse_complex_herm_precon_ilu (f11jn). 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 easily be determined a priori, as it can depend dramatically on the conditioning and spectrum of the preconditioned matrix of the coefficients A-=M-1A.

Example

This example solves a complex sparse Hermitian positive definite system of equations using the conjugate gradient method, with incomplete Cholesky preconditioning.
function f11jq_example


fprintf('f11jq example results\n\n');

% Solve sparse Hermitian system Ax = b using CG method with
% Incomplete Cholesky preconditioning (IC)

% Define A and b 
n  = int64(9);
nz = int64(23);
a    = zeros(3*nz,1);
irow = zeros(3*nz, 1, 'int64');
icol = zeros(3*nz, 1, 'int64');
a(1:nz) = [ 6 + 0.i; -1 + 1.i;  6 + 0.i;  0 + 1.i;
            5 + 0.i;  5 + 0.i;  2 - 2.i;  4 + 0.i;
            1 + 1.i;  2 + 0.i;  6 + 0.i; -4 + 3.i;
            0 + 1.i; -1 + 0.i;  6 + 0.i; -1 - 1.i;
            0 - 1.i;  9 + 0.i;  1 + 3.i;  1 + 2.i;
           -1 + 0.i;  1 + 4.i;  9 + 0.i];
b       = [ 8 + 54i;-10 - 92i; 25 + 27i; 26 - 28i;
           54 + 12i; 26 - 22i; 47 + 65i; 71 - 57i;
           60 + 70i];
irow(1:nz) = int64([1;2;2;3;3;4;5;5;6;6;6;7;7;7;7;8;8;8;9;9;9;9;9]);
icol(1:nz) = int64([1;1;2;2;3;4;1;5;3;4;6;2;5;6;7;4;6;8;1;5;6;8;9]);

% Calculate incomplete Cholesky factorization
lfill  = int64(0);
dtol   = 0;
mic    = 'N';
dscale = 0;
ipiv   = zeros(n, 1, 'int64');

[a, irow, icol, ipiv, istr, nnzc, npivm, ifail] = ...
  f11jn( ...
         nz, a, irow, icol, lfill, dtol, mic, dscale, ipiv);

% Solve Ax = b
method = 'CG';
tol    = 1e-06;
maxitn = int64(100);
x      = complex(zeros(n, 1));

[x, rnorm, itn, ifail] = ...
  f11jq( ...
         method, nz, a, irow, icol, ipiv, istr, b, tol, maxitn, x);

fprintf('Converged in %d iterations\n', itn);
fprintf('Final redidual norm = %16.3d\n\n', rnorm);
disp('Solution');
disp(x);


f11jq example results

Converged in 5 iterations
Final redidual norm =        3.197e-14

Solution
   1.0000 + 9.0000i
   2.0000 - 8.0000i
   3.0000 + 7.0000i
   4.0000 - 6.0000i
   5.0000 + 5.0000i
   6.0000 - 4.0000i
   7.0000 + 3.0000i
   8.0000 - 2.0000i
   9.0000 + 1.0000i


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–2015