nag_sparse_nherm_precon_ilu_solve (f11dpc) (PDF version)
f11 Chapter Contents
f11 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_sparse_nherm_precon_ilu_solve (f11dpc)

 Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_sparse_nherm_precon_ilu_solve (f11dpc) solves a system of complex linear equations involving the incomplete LU preconditioning matrix generated by nag_sparse_nherm_fac (f11dnc).

2  Specification

#include <nag.h>
#include <nagf11.h>
void  nag_sparse_nherm_precon_ilu_solve (Nag_TransType trans, Integer n, const Complex a[], Integer la, const Integer irow[], const Integer icol[], const Integer ipivp[], const Integer ipivq[], const Integer istr[], const Integer idiag[], Nag_SparseNsym_CheckData check, const Complex y[], Complex x[], NagError *fail)

3  Description

nag_sparse_nherm_precon_ilu_solve (f11dpc) solves a system of complex 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 complex sparse matrix stored in coordinate storage (CS) format (see Section 2.1.1 in the f11 Chapter Introduction), as generated by nag_sparse_nherm_fac (f11dnc).
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_nherm_precon_ilu_solve (f11dpc) 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_nherm_fac (f11dnc). The permutation matrices P and Q are returned from nag_sparse_nherm_fac (f11dnc) via the arrays ipivp and ipivq.
It is envisaged that a common use of nag_sparse_nherm_precon_ilu_solve (f11dpc) will be to carry out the preconditioning step required in the application of nag_sparse_nherm_basic_solver (f11bsc) to sparse complex linear systems. nag_sparse_nherm_precon_ilu_solve (f11dpc) is used for this purpose by the Black Box function nag_sparse_nherm_fac_sol (f11dqc).
nag_sparse_nherm_precon_ilu_solve (f11dpc) may also be used in combination with nag_sparse_nherm_fac (f11dnc) to solve a sparse system of complex linear equations directly (see Section 9.5 in nag_sparse_nherm_fac (f11dnc)). This use of nag_sparse_nherm_precon_ilu_solve (f11dpc) is illustrated in Section 10.

4  References

None.

5  Arguments

1:     trans Nag_TransTypeInput
On entry: specifies whether or not the matrix M is transposed.
trans=Nag_NoTrans
Mx=y is solved.
trans=Nag_Trans
MTx=y is solved.
Constraint: trans=Nag_NoTrans or Nag_Trans.
2:     n IntegerInput
On entry: n, the order of the matrix M. This must be the same value as was supplied in the preceding call to nag_sparse_nherm_fac (f11dnc).
Constraint: n1.
3:     a[la] const ComplexInput
On entry: the values returned in the array a by a previous call to nag_sparse_nherm_fac (f11dnc).
4:     la IntegerInput
On entry: the dimension of the arrays a, irow and icol. This must be the same value supplied in the preceding call to nag_sparse_nherm_fac (f11dnc).
5:     irow[la] const IntegerInput
6:     icol[la] const IntegerInput
7:     ipivp[n] const IntegerInput
8:     ipivq[n] const IntegerInput
9:     istr[n+1] const IntegerInput
10:   idiag[n] const IntegerInput
On entry: the values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_nherm_fac (f11dnc).
11:   check Nag_SparseNsym_CheckDataInput
On entry: specifies whether or not the CS representation of the matrix M should be checked.
check=Nag_SparseNsym_Check
Checks are carried on the values of n, irow, icol, ipivp, ipivq, istr and idiag.
check=Nag_SparseNsym_NoCheck
None of these checks are carried out.
See also Section 9.2.
Constraint: check=Nag_SparseNsym_Check or Nag_SparseNsym_NoCheck.
12:   y[n] const ComplexInput
On entry: the right-hand side vector y.
13:   x[n] ComplexOutput
On exit: the solution vector x.
14:   fail NagError *Input/Output
The NAG error argument (see Section 3.6 in the Essential Introduction).

6  Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.2.1.2 in the Essential Introduction for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_INT
On entry, n=value.
Constraint: n1.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
An unexpected error has been triggered by this function. Please contact NAG.
See Section 3.6.6 in the Essential Introduction for further information.
NE_INVALID_CS
On entry, i=value, icol[i-1]=value, and n=value.
Constraint: icol[i-1]1 and icol[i-1]n.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, i=value, irow[i-1]=value, n=value.
Constraint: irow[i-1]1 and irow[i-1]n.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
NE_INVALID_CS_PRECOND
On entry, idiag[i-1] appears to be incorrect: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, istr appears to be invalid.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, istr[i-1] is inconsistent with irow: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
NE_INVALID_ROWCOL_PIVOT
On entry, i=value, ipivp[i-1]=value, n=value.
Constraint: ipivp[i-1]1 and ipivp[i-1]n.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, i=value, ipivq[i-1]=value, n=value.
Constraint: ipivq[i-1]1 and ipivq[i-1]n.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, ipivp[i-1] is a repeated value: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, ipivq[i-1] is a repeated value: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 3.6.5 in the Essential Introduction for further information.
NE_NOT_STRICTLY_INCREASING
On entry, a[i-1] is out of order: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).
On entry, the location (irow[i-1],icol[i-1]) is a duplicate: i=value.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_nherm_precon_ilu_solve (f11dpc) and nag_sparse_nherm_fac (f11dnc).

7  Accuracy

If trans=Nag_NoTrans 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=Nag_Trans.

8  Parallelism and Performance

Not applicable.

9  Further Comments

9.1  Timing

The time taken for a call to nag_sparse_nherm_precon_ilu_solve (f11dpc) is proportional to the value of nnzc returned from nag_sparse_nherm_fac (f11dnc).

9.2  Use of check

It is expected that a common use of nag_sparse_nherm_precon_ilu_solve (f11dpc) will be to carry out the preconditioning step required in the application of nag_sparse_nherm_basic_solver (f11bsc) to sparse complex linear systems. In this situation nag_sparse_nherm_precon_ilu_solve (f11dpc) 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=Nag_SparseNsym_Check for the first of such calls, and to set check=Nag_SparseNsym_NoCheck for all subsequent calls.

10  Example

This example reads in a complex sparse non-Hermitian matrix A and a vector y. It then calls nag_sparse_nherm_fac (f11dnc), with lfill=-1 and dtol=0.0, to compute the complete LU decomposition
A=PLDUQ.  
Finally it calls nag_sparse_nherm_precon_ilu_solve (f11dpc) to solve the system
PLDUQx=y.  

10.1  Program Text

Program Text (f11dpce.c)

10.2  Program Data

Program Data (f11dpce.d)

10.3  Program Results

Program Results (f11dpce.r)


nag_sparse_nherm_precon_ilu_solve (f11dpc) (PDF version)
f11 Chapter Contents
f11 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2015