PDF version (NAG web site
, 64-bit version, 64-bit version)
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
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
according to the value of the argument
trans, where the matrix
, corresponds to an incomplete
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
is a lower triangular sparse matrix with unit diagonal elements,
is a diagonal matrix,
is an upper triangular sparse matrix with unit diagonal elements and,
and
are permutation matrices.
,
and
are supplied to
nag_sparse_real_gen_precon_ilu_solve (f11db) through the matrix
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
and
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:
– string (length ≥ 1)
-
Specifies whether or not the matrix
is transposed.
- is solved.
- is solved.
Constraint:
or .
- 2:
– double array
-
The values returned in the array
a by a previous call to
nag_sparse_real_gen_precon_ilu (f11da).
- 3:
– int64int32nag_int array
- 4:
– int64int32nag_int array
- 5:
– int64int32nag_int array
- 6:
– int64int32nag_int array
- 7:
– int64int32nag_int array
- 8:
– 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:
– string (length ≥ 1)
-
Specifies whether or not the CS representation of the matrix
should be checked.
- Checks are carried on the values of n, irow, icol, ipivp, ipivq, istr and idiag.
- None of these checks are carried out.
Constraint:
or .
- 10:
– double array
-
The right-hand side vector .
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the arrays
ipivp,
ipivq,
istr,
idiag,
y. (An error is raised if these dimensions are not equal.)
, the order of the matrix
. This
must be the same value as was supplied in the preceding call to
nag_sparse_real_gen_precon_ilu (f11da).
Constraint:
.
- 2:
– 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:
– double array
-
The solution vector .
- 2:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
Errors or warnings detected by the function:
-
-
On entry, | or , |
or | or . |
-
-
-
-
On entry, the CS representation of the preconditioning matrix
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.
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
If
the computed solution
is the exact solution of a perturbed system of equations
, where
is a modest linear function of
, and
is the
machine precision. An equivalent result holds when
.
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
. In the interests of both reliability and efficiency, you are recommended to set
for the first of such calls, and for all subsequent calls set
.
Example
This example reads in a sparse nonsymmetric matrix
and a vector
. It then calls
nag_sparse_real_gen_precon_ilu (f11da), with
and
, to compute the
complete
decomposition
Finally it calls
nag_sparse_real_gen_precon_ilu_solve (f11db) to solve the system
Open in the MATLAB editor:
f11db_example
function f11db_example
fprintf('f11db example results\n\n');
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];
lfill = int64(-1);
dtol = 0;
milu = 'No modification';
ipivp = zeros(n, 1, 'int64');
ipivq = zeros(n, 1, 'int64');
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
f11da(...
nz, a, irow, icol, lfill, dtol, milu, ipivp, ipivq);
if npivm>0
error('Factorization is not complete');
end
[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)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015