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_gen_precon_ilu_solve (f11db)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparse_real_gen_precon_ilu_solve (f11db) solves a system of linear equations involving the incomplete LU preconditioning matrix generated by nag_sparse_real_gen_precon_ilu (f11da).

Syntax

[x, ifail] = f11db(trans, a, irow, icol, ipivp, ipivq, istr, idiag, check, y, 'n', n, 'la', la)
[x, ifail] = nag_sparse_real_gen_precon_ilu_solve(trans, a, irow, icol, ipivp, ipivq, istr, idiag, check, y, 'n', n, 'la', la)

Description

nag_sparse_real_gen_precon_ilu_solve (f11db) solves a system of linear equations
Mx=y, or MTx=y,  
according to the value of the argument trans, where the matrix M=PLDUQ, corresponds to an incomplete LU decomposition of a sparse matrix stored in coordinate storage (CS) format (see Coordinate storage (CS) format in the F11 Chapter Introduction), as generated by nag_sparse_real_gen_precon_ilu (f11da).
In the above decomposition L is a lower triangular sparse matrix with unit diagonal elements, D is a diagonal matrix, U is an upper triangular sparse matrix with unit diagonal elements and, P and Q are permutation matrices. L, D and U are supplied to nag_sparse_real_gen_precon_ilu_solve (f11db) through the matrix
C=L+D-1+U-2I  
which is an n by n sparse matrix, stored in CS format, as returned by nag_sparse_real_gen_precon_ilu (f11da). The permutation matrices P and Q are returned from nag_sparse_real_gen_precon_ilu (f11da) via the arrays ipivp and ipivq.
It is envisaged that a common use of nag_sparse_real_gen_precon_ilu_solve (f11db) will be to carry out the preconditioning step required in the application of nag_sparse_real_gen_basic_solver (f11be) to sparse linear systems. nag_sparse_real_gen_precon_ilu_solve (f11db) is used for this purpose by the Black Box function nag_sparse_real_gen_solve_ilu (f11dc).
nag_sparse_real_gen_precon_ilu_solve (f11db) may also be used in combination with nag_sparse_real_gen_precon_ilu (f11da) to solve a sparse system of linear equations directly (see Direct Solution of Sparse Linear Systems in nag_sparse_real_gen_precon_ilu (f11da)). This use of nag_sparse_real_gen_precon_ilu_solve (f11db) is demonstrated in Example.

References

None.

Parameters

Compulsory Input Parameters

1:     trans – string (length ≥ 1)
Specifies whether or not the matrix M is transposed.
trans='N'
Mx=y is solved.
trans='T'
MTx=y is solved.
Constraint: trans='N' or 'T'.
2:     ala – double array
The values returned in the array a by a previous call to nag_sparse_real_gen_precon_ilu (f11da).
3:     irowla int64int32nag_int array
4:     icolla int64int32nag_int array
5:     ipivpn int64int32nag_int array
6:     ipivqn int64int32nag_int array
7:     istrn+1 int64int32nag_int array
8:     idiagn int64int32nag_int array
The values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_real_gen_precon_ilu (f11da).
9:     check – string (length ≥ 1)
Specifies whether or not the CS representation of the matrix M should be checked.
check='C'
Checks are carried on the values of n, irow, icol, ipivp, ipivq, istr and idiag.
check='N'
None of these checks are carried out.
See also Use of check.
Constraint: check='C' or 'N'.
10:   yn – double array
The right-hand side vector y.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the arrays ipivp, ipivq, istr, idiag, y. (An error is raised if these dimensions are not equal.)
n, the order of the matrix M. This must be the same value as was supplied in the preceding call to nag_sparse_real_gen_precon_ilu (f11da).
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 returned by the preceding call to nag_sparse_real_gen_precon_ilu (f11da).

Output Parameters

1:     xn – double array
The solution vector x.
2:     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:
   ifail=1
On entry,trans'N' or 'T',
orcheck'C' or 'N'.
   ifail=2
On entry,n<1.
   ifail=3
On entry, the CS representation of the preconditioning matrix M is invalid. Further details are given in the error message. Check that the call to nag_sparse_real_gen_precon_ilu_solve (f11db) has been preceded by a valid call to nag_sparse_real_gen_precon_ilu (f11da) and that the arrays a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between the two calls.
   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

If trans='N' the computed solution x is the exact solution of a perturbed system of equations M+δMx=y, where
δMcnεPLDUQ,  
cn is a modest linear function of n, and ε is the machine precision. An equivalent result holds when trans='T'.

Further Comments

Timing

The time taken for a call to nag_sparse_real_gen_precon_ilu_solve (f11db) is proportional to the value of nnzc returned from nag_sparse_real_gen_precon_ilu (f11da).

Use of check

It is expected that a common use of nag_sparse_real_gen_precon_ilu_solve (f11db) will be to carry out the preconditioning step required in the application of nag_sparse_real_gen_basic_solver (f11be) to sparse linear systems. In this situation nag_sparse_real_gen_precon_ilu_solve (f11db) is likely to be called many times with the same matrix M. In the interests of both reliability and efficiency, you are recommended to set check='C' for the first of such calls, and for all subsequent calls set check='N'.

Example

This example reads in a sparse nonsymmetric matrix A and a vector y. It then calls nag_sparse_real_gen_precon_ilu (f11da), with lfill=-1 and dtol=0.0, to compute the complete LU decomposition
A=PLDUQ.  
Finally it calls nag_sparse_real_gen_precon_ilu_solve (f11db) to solve the system
PLDUQx=y.  
function f11db_example


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

% Solve sparse system Ax = b using complete LU decomposition

% Define sparse matrix A and RHS B
n    = int64(4);
nz   = int64(11);
a    = zeros(3*nz, 1);
irow = zeros(3*nz, 1, 'int64');
icol = irow;
sparseA = [  1   1    2;
             1   1    3;
            -1   2    1;
             2   2    3;
             2   2    4;
             3   3    1;
            -2   3    4;
             1   4    1;
            -2   4    2;
             1   4    3;
             1   4    4];
a(1:nz)    = sparseA(:,1);
irow(1:nz) = int64(sparseA(:,2));
icol(1:nz) = int64(sparseA(:,3));

b          = [5; 13; -5; 4];

% ILU preconditioner
%     Input parameters
lfill = int64(-1);
dtol  = 0;
milu  = 'No modification';
ipivp = zeros(n, 1, 'int64');
ipivq = zeros(n, 1, 'int64');
% Compute complete factorization
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
    f11da(...
          nz, a, irow, icol, lfill, dtol, milu, ipivp, ipivq);

% Check whether factorization is complete
if npivm>0
  error('Factorization is not complete');
end

% Solution Setup
%      Input parameters
[x, ifail] = f11db( ...
                    'N', a, irow, icol, ipivp, ipivq, istr, idiag, 'C', b);

disp('Solution of linear system');
disp(x);


f11db example results

Solution of linear system
     1
     2
     3
     4


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