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)

Purpose

nag_sparse_real_gen_precon_ilu_solve (f11db) solves a system of linear equations involving the incomplete LULU 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,
Mx=y, or MTx=y,
according to the value of the parameter trans, where the matrix M = PLDUQM=PLDUQ, corresponds to an incomplete LULU decomposition of a sparse matrix stored in coordinate storage (CS) format (see Section [Coordinate storage (CS) format] in the F11 Chapter Introduction), as generated by nag_sparse_real_gen_precon_ilu (f11da).
In the above decomposition LL is a lower triangular sparse matrix with unit diagonal elements, DD is a diagonal matrix, UU is an upper triangular sparse matrix with unit diagonal elements and, PP and QQ are permutation matrices. LL, DD and UU are supplied to nag_sparse_real_gen_precon_ilu_solve (f11db) through the matrix
C = L + D1 + U2I
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 PP and QQ 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 Section [Direct Solution of Sparse Linear Systems] in (f11da)). This use of nag_sparse_real_gen_precon_ilu_solve (f11db) is demonstrated in Section [Example].

References

None.

Parameters

Compulsory Input Parameters

1:     trans – string (length ≥ 1)
Specifies whether or not the matrix MM is transposed.
trans = 'N'trans='N'
Mx = yMx=y is solved.
trans = 'T'trans='T'
MTx = yMTx=y is solved.
Constraint: trans = 'N'trans='N' or 'T''T'.
2:     a(la) – double array
The values returned in the array a by a previous call to nag_sparse_real_gen_precon_ilu (f11da).
3:     irow(la) – int64int32nag_int array
4:     icol(la) – int64int32nag_int array
5:     ipivp(n) – int64int32nag_int array
6:     ipivq(n) – int64int32nag_int array
7:     istr(n + 1n+1) – int64int32nag_int array
8:     idiag(n) – 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 MM should be checked.
check = 'C'check='C'
Checks are carried on the values of n, irow, icol, ipivp, ipivq, istr and idiag.
check = 'N'check='N'
None of these checks are carried out.
Constraint: check = 'C'check='C' or 'N''N'.
10:   y(n) – double array
n, the dimension of the array, must satisfy the constraint n1n1.
The right-hand side vector yy.

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.)
nn, the order of the matrix MM. This must be the same value as was supplied in the preceding call to nag_sparse_real_gen_precon_ilu (f11da).
Constraint: n1n1.
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 as declared in the (sub)program from which nag_sparse_real_gen_precon_ilu_solve (f11db) is called. This must be the same value returned by the preceding call to nag_sparse_real_gen_precon_ilu (f11da).

Input Parameters Omitted from the MATLAB Interface

None.

Output Parameters

1:     x(n) – double array
The solution vector xx.
2:     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,trans'N'trans'N' or 'T''T',
orcheck'C'check'C' or 'N''N'.
  ifail = 2ifail=2
On entry,n < 1n<1.
  ifail = 3ifail=3
On entry, the CS representation of the preconditioning matrix MM 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.

Accuracy

If trans = 'N'trans='N' the computed solution xx is the exact solution of a perturbed system of equations (M + δM)x = y(M+δM)x=y, where
|δM|c(n)εP|L||D||U|Q,
|δM|c(n)εP|L||D||U|Q,
c(n)c(n) is a modest linear function of nn, and εε is the machine precision. An equivalent result holds when trans = 'T'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 MM. In the interests of both reliability and efficiency, you are recommended to set check = 'C'check='C' for the first of such calls, and for all subsequent calls set check = 'N'check='N'.

Example

function nag_sparse_real_gen_precon_ilu_solve_example
nz = int64(11);
a = zeros(2*nz, 1);
a(1:nz) = [1; 1; -1; 2; 2; 3; -2; 1; -2; 1; 1];
irow = zeros(2*nz, 1, 'int64');
irow(1:nz) = [1; 1; 2; 2; 2; 3; 3; 4; 4; 4; 4];
icol = zeros(2*nz, 1, 'int64');
icol(1:nz) = [2; 3; 1; 3; 4; 1; 4; 1; 2; 3; 4];
lfill = int64(-1);
dtol = 0;
pstrat = 'C';
milu = 'N';
ipivp = zeros(4, 1, 'int64');
ipivq = zeros(4, 1, 'int64');
check = 'C';
trans = 'N';
y = [5; 13; -5; 4];
% Calculate LU factorization
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
    nag_sparse_real_gen_precon_ilu(nz, a, irow, icol, lfill, dtol, milu, ...
                                   ipivp, ipivq);

if npivm > 0
  fprintf('\nFactorization is not complete.\n');
else
  % Solve P L D U x = y
  [x, ifail] = nag_sparse_real_gen_precon_ilu_solve(trans, a, irow, icol, ...
                                          ipivp, ipivq, istr, idiag, check, y)
end
 

x =

     1
     2
     3
     4


ifail =

                    0


function f11db_example
nz = int64(11);
a = zeros(2*nz, 1);
a(1:nz) = [1; 1; -1; 2; 2; 3; -2; 1; -2; 1; 1];
irow = zeros(2*nz, 1, 'int64');
irow(1:nz) = [1; 1; 2; 2; 2; 3; 3; 4; 4; 4; 4];
icol = zeros(2*nz, 1, 'int64');
icol(1:nz) = [2; 3; 1; 3; 4; 1; 4; 1; 2; 3; 4];
lfill = int64(-1);
dtol = 0;
pstrat = 'C';
milu = 'N';
ipivp = zeros(4, 1, 'int64');
ipivq = zeros(4, 1, 'int64');
check = 'C';
trans = 'N';
y = [5; 13; -5; 4];
% Calculate LU factorization
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
    f11da(nz, a, irow, icol, lfill, dtol, milu, ipivp, ipivq);

if npivm > 0
  fprintf('\nFactorization is not complete.\n');
else
  % Solve P L D U x = y
  [x, ifail] = f11db(trans, a, irow, icol, ipivp, ipivq, istr, idiag, check, y)
end
 

x =

     1
     2
     3
     4


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