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_precon_ilu (f11jn)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparse_complex_herm_precon_ilu (f11jn) computes an incomplete Cholesky factorization of a complex sparse Hermitian matrix, represented in symmetric coordinate storage format. This factorization may be used as a preconditioner in combination with nag_sparse_complex_herm_solve_ilu (f11jq).

Syntax

[a, irow, icol, ipiv, istr, nnzc, npivm, ifail] = f11jn(nz, a, irow, icol, lfill, dtol, mic, dscale, ipiv, 'n', n, 'la', la, 'pstrat', pstrat)
[a, irow, icol, ipiv, istr, nnzc, npivm, ifail] = nag_sparse_complex_herm_precon_ilu(nz, a, irow, icol, lfill, dtol, mic, dscale, ipiv, 'n', n, 'la', la, 'pstrat', pstrat)

Description

nag_sparse_complex_herm_precon_ilu (f11jn) computes an incomplete Cholesky factorization (see Meijerink and Van der Vorst (1977)) of a complex sparse Hermitian n by n matrix A. It is designed specifically for positive definite matrices, but may also work for some mildly indefinite cases. The factorization is intended primarily for use as a preconditioner with the complex Hermitian iterative solver nag_sparse_complex_herm_solve_ilu (f11jq).
The decomposition is written in the form
A=M+R  
where
M=PLDLHPT  
and P is a permutation matrix, L is lower triangular complex with unit diagonal elements, D is real diagonal and R is a remainder matrix.
The amount of fill-in occurring in the factorization can vary from zero to complete fill, and can be controlled by specifying either the maximum level of fill lfill, or the drop tolerance dtol. The factorization may be modified in order to preserve row sums, and the diagonal elements may be perturbed to ensure that the preconditioner is positive definite. Diagonal pivoting may optionally be employed, either with a user-defined ordering, or using the Markowitz strategy (see Markowitz (1957)), which aims to minimize fill-in. For further details see Further Comments.
The sparse matrix A is represented in symmetric coordinate storage (SCS) format (see Symmetric coordinate storage (SCS) format in the F11 Chapter Introduction). The array a stores all the nonzero elements of the lower triangular part of A, while arrays irow and icol store the corresponding row and column indices respectively. Multiple nonzero elements may not be specified for the same row and column index.
The preconditioning matrix M is returned in terms of the SCS representation of the lower triangular matrix
C=L+D-1-I.  

References

Chan T F (1991) Fourier analysis of relaxed incomplete factorization preconditioners SIAM J. Sci. Statist. Comput. 12(2) 668–680
Markowitz H M (1957) The elimination form of the inverse and its application to linear programming Management Sci. 3 255–269
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
Salvini S A and Shaw G J (1995) An evaluation of new NAG Library solvers for large sparse symmetric linear systems NAG Technical Report TR1/95
Van der Vorst H A (1990) The convergence behaviour of preconditioned CG and CG-S in the presence of rounding errors Lecture Notes in Mathematics (eds O Axelsson and L Y Kolotilina) 1457 Springer–Verlag

Parameters

Compulsory Input Parameters

1:     nz int64int32nag_int scalar
The number of nonzero elements in the lower triangular part of the matrix A.
Constraint: 1nzn×n+1/2.
2:     ala – complex array
The nonzero elements in the lower triangular part of the matrix A, ordered by increasing row index, and by increasing column index within each row. Multiple entries for the same row and column indices are not permitted. The function nag_sparse_complex_herm_sort (f11zp) may be used to order the elements in this way.
3:     irowla int64int32nag_int array
4:     icolla int64int32nag_int array
The row and column indices of the nonzero elements supplied in a.
Constraints:
irow and icol must satisfy these constraints (which may be imposed by a call to nag_sparse_complex_herm_sort (f11zp)):
  • 1irowin and 1icoliirowi, for i=1,2,,nz;
  • irowi-1<irowi or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nz.
5:     lfill int64int32nag_int scalar
If lfill0 its value is the maximum level of fill allowed in the decomposition (see Control of Fill-in). A negative value of lfill indicates that dtol will be used to control the fill instead.
6:     dtol – double scalar
If lfill<0, dtol is used as a drop tolerance to control the fill-in (see Control of Fill-in); otherwise dtol is not referenced.
Constraint: if lfill<0, dtol0.0.
7:     mic – string (length ≥ 1)
Indicates whether or not the factorization should be modified to preserve row sums (see Choice of s).
mic='M'
The factorization is modified.
mic='N'
The factorization is not modified.
Constraint: mic='M' or 'N'.
8:     dscale – double scalar
The diagonal scaling parameter. All diagonal elements are multiplied by the factor (1.0+dscale) at the start of the factorization. This can be used to ensure that the preconditioner is positive definite. See also Choice of s.
9:     ipivn int64int32nag_int array
If pstrat='U', ipivi must specify the row index of the diagonal element to be used as a pivot at elimination stage i. Otherwise ipiv need not be initialized.
Constraint: if pstrat='U', ipiv must contain a valid permutation of the integers on 1,n.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the array ipiv.
n, the order of the matrix A.
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. these arrays must be of sufficient size to store both A (nz elements) and C (nnzc elements).
Constraint: la2×nz.
3:     pstrat – string (length ≥ 1)
Default: 'M'
Specifies the pivoting strategy to be adopted.
pstrat='N'
No pivoting is carried out.
pstrat='M'
Diagonal pivoting aimed at minimizing fill-in is carried out, using the Markowitz strategy (see Markowitz (1957)).
pstrat='U'
Diagonal pivoting is carried out according to the user-defined input array ipiv.
Constraint: pstrat='N', 'M' or 'U'.

Output Parameters

1:     ala – complex array
The first nz elements of a contain the nonzero elements of A and the next nnzc elements contain the elements of the lower triangular matrix C. Matrix elements are ordered by increasing row index, and by increasing column index within each row.
2:     irowla int64int32nag_int array
3:     icolla int64int32nag_int array
The row and column indices of the nonzero elements returned in a.
4:     ipivn int64int32nag_int array
The pivot indices. If ipivi=j, the diagonal element in row j was used as the pivot at elimination stage i.
5:     istrn+1 int64int32nag_int array
istri, for i=1,2,,n, is the starting address in the arrays a, irow and icol of row i of the matrix C. istrn+1 is the address of the last nonzero element in C plus one.
6:     nnzc int64int32nag_int scalar
The number of nonzero elements in the lower triangular matrix C.
7:     npivm int64int32nag_int scalar
The number of pivots which were modified during the factorization to ensure that M was positive definite. The quality of the preconditioner will generally depend on the returned value of npivm. If npivm is large the preconditioner may not be satisfactory. In this case it may be advantageous to call nag_sparse_complex_herm_precon_ilu (f11jn) again with an increased value of either lfill or dscale. See also Choice of s and Direct Solution of Systems.
8:     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,n<1,
ornz<1,
ornz>n×n+1/2,
orla<2×nz,
ordtol<0.0,
ormic'M' or 'N',
orpstrat'N', 'M' or 'U',
orliwork is too small.
   ifail=2
On entry, the arrays irow and icol fail to satisfy the following constraints:
  • 1irowin and 1icoliirowi, for i=1,2,,nz;
  • irowi-1<irowi, or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nz.
Therefore a nonzero element has been supplied which does not lie in the lower triangular part of A, is out of order, or has duplicate row and column indices. Call nag_sparse_complex_herm_sort (f11zp) to reorder and sum or remove duplicates.
   ifail=3
On entry, pstrat='U', but ipiv does not represent a valid permutation of the integers in 1,n. An input value of ipiv is either out of range or repeated.
   ifail=4
la is too small, resulting in insufficient storage space for fill-in elements. The decomposition has been terminated before completion. Either increase la or reduce the amount of fill by setting pstrat='M', reducing lfill, or increasing dtol.
   ifail=5 (nag_sparse_complex_herm_sort (f11zp))
A serious error has occurred in an internal call to the specified 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

The accuracy of the factorization will be determined by the size of the elements that are dropped and the size of any modifications made to the diagonal elements. If these sizes are small then the computed factors will correspond to a matrix close to A. The factorization can generally be made more accurate by increasing lfill, or by reducing dtol with lfill<0.
If nag_sparse_complex_herm_precon_ilu (f11jn) is used in combination with nag_sparse_complex_herm_solve_ilu (f11jq), the more accurate the factorization the fewer iterations will be required. However, the cost of the decomposition will also generally increase.

Further Comments

Timing

The time taken for a call to nag_sparse_complex_herm_precon_ilu (f11jn) is roughly proportional to nnzc2/n.

Control of Fill-in

If lfill0, the amount of fill-in occurring in the incomplete factorization is controlled by limiting the maximum ‘level’ of fill-in to lfill. The original nonzero elements of A are defined to be of level 0. The fill level of a new nonzero location occurring during the factorization is defined as:
k=maxke,kc+1,  
where ke is the level of fill of the element being eliminated, and kc is the level of fill of the element causing the fill-in.
If lfill<0, the fill-in is controlled by means of the ‘drop tolerance’ dtol. A potential fill-in element aij occurring in row i and column j will not be included if
aij<dtol×aiiajj.  
For either method of control, any elements which are not included are discarded if mic='N', or subtracted from the diagonal element in the elimination row if mic='M'.

Choice of Arguments

There is unfortunately no choice of the various algorithmic arguments which is optimal for all types of complex Hermitian matrix, and some experimentation will generally be required for each new type of matrix encountered.
If the matrix A is not known to have any particular special properties, the following strategy is recommended. Start with lfill=0, mic='N' and dscale=0.0. If the value returned for npivm is significantly larger than zero, i.e., a large number of pivot modifications were required to ensure that M was positive definite, the preconditioner is not likely to be satisfactory. In this case increase either lfill or dscale until npivm falls to a value close to zero. Once suitable values of lfill and dscale have been found try setting mic='M' to see if any improvement can be obtained by using modified incomplete Cholesky.
nag_sparse_complex_herm_precon_ilu (f11jn) is primarily designed for positive definite matrices, but may work for some mildly indefinite problems. If npivm cannot be satisfactorily reduced by increasing lfill or dscale then A is probably too indefinite for this function.
For certain classes of matrices (typically those arising from the discretization of elliptic or parabolic partial differential equations), the convergence rate of the preconditioned iterative solver can sometimes be significantly improved by using an incomplete factorization which preserves the row-sums of the original matrix. In these cases try setting mic='M'.

Direct Solution of positive definite Systems

Although it is not their primary purpose, nag_sparse_complex_herm_precon_ilu (f11jn) and nag_sparse_complex_herm_precon_ilu_solve (f11jp) may be used together to obtain a direct solution to a complex Hermitian positive definite linear system. To achieve this the call to nag_sparse_complex_herm_precon_ilu_solve (f11jp) should be preceded by a complete Cholesky factorization
A=PLDLHPT=M.  
A complete factorization is obtained from a call to nag_sparse_complex_herm_precon_ilu (f11jn) with lfill<0 and dtol=0.0, provided npivm=0 on exit. A nonzero value of npivm indicates that a is not positive definite, or is ill-conditioned. A factorization with nonzero npivm may serve as a preconditioner, but will not result in a direct solution. It is therefore essential to check the output value of npivm if a direct solution is required.
The use of nag_sparse_complex_herm_precon_ilu (f11jn) and nag_sparse_complex_herm_precon_ilu_solve (f11jp) as a direct method is illustrated in nag_sparse_complex_herm_precon_ilu_solve (f11jp).

Example

This example reads in a complex sparse Hermitian matrix A and calls nag_sparse_complex_herm_precon_ilu (f11jn) to compute an incomplete Cholesky factorization. It then outputs the nonzero elements of both A and C=L+D-1-I.
The call to nag_sparse_complex_herm_precon_ilu (f11jn) has lfill=0, mic='N', dscale=0.0 and pstrat='M', giving an unmodified zero-fill factorization of an unperturbed matrix, with Markowitz diagonal pivoting.
function f11jn_example


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

% The sparse matrix A
n    = 7;
nz   = int64(16);
a    = zeros(3*nz, 1);
irow = zeros(3*nz, 1, 'int64');
icol = zeros(3*nz, 1, 'int64');

a(1:nz) = [ 6 + 0i;  1 - 2i;  9 + 0i; 4 + 0i; 
            2 + 2i;  5 + 0i;  0 - 1i; 1 + 0i;
            4 + 0i;  1 + 3i;  0 - 2i; 3 + 0i;
            2 + 1i; -1 + 0i; -3 - 1i; 5 + 0i];

irow(1:nz) = int64([1; 2; 2; 3; 4; 4; 5; 5; 5; 6; 6; 6; 7; 7; 7; 7]);
icol(1:nz) = int64([1; 1; 2; 3; 2; 4; 1; 4; 5; 2; 5; 6; 1; 2; 3; 7]);

% Zero-fill, unmodified incmoplee Cholesky factorization of A
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);

% Display details of factorization
disp('Details of Incomplete Cholesky Factorization')
range = nz+1:nz+nnzc;
S = sparse(double(irow(range)),double(icol(range)),a(range));
disp(S);

disp('Pivots, ipiv:');
fprintf('%5d',ipiv);
fprintf('\n');


f11jn example results

Details of Incomplete Cholesky Factorization
   (1,1)      0.2500 + 0.0000i
   (6,1)     -0.7500 - 0.2500i
   (2,2)      0.2000 + 0.0000i
   (3,2)      0.2000 + 0.0000i
   (7,2)      0.4000 - 0.4000i
   (3,3)      0.2632 + 0.0000i
   (4,3)      0.0000 - 0.5263i
   (5,3)      0.0000 + 0.2632i
   (4,4)      0.5135 + 0.0000i
   (7,4)      0.5135 - 1.5405i
   (5,5)      0.1743 + 0.0000i
   (6,5)      0.3486 + 0.1743i
   (7,5)      0.1743 - 0.3486i
   (6,6)      0.6141 + 0.0000i
   (7,6)     -0.6141 + 0.5352i
   (7,7)      3.1974 + 0.0000i

Pivots, ipiv:
    3    4    5    6    1    7    2

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