NAG FL Interface
f11jnf (complex_​herm_​precon_​ichol)

1 Purpose

f11jnf 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 f11jqf.

2 Specification

Fortran Interface
Subroutine f11jnf ( n, nnz, a, la, irow, icol, lfill, dtol, mic, dscale, pstrat, ipiv, istr, nnzc, npivm, iwork, liwork, ifail)
Integer, Intent (In) :: n, nnz, la, lfill, liwork
Integer, Intent (Inout) :: irow(la), icol(la), ipiv(n), ifail
Integer, Intent (Out) :: istr(n+1), nnzc, npivm, iwork(liwork)
Real (Kind=nag_wp), Intent (In) :: dtol, dscale
Complex (Kind=nag_wp), Intent (Inout) :: a(la)
Character (1), Intent (In) :: mic, pstrat
C Header Interface
#include <nag.h>
void  f11jnf_ (const Integer *n, const Integer *nnz, Complex a[], const Integer *la, Integer irow[], Integer icol[], const Integer *lfill, const double *dtol, const char *mic, const double *dscale, const char *pstrat, Integer ipiv[], Integer istr[], Integer *nnzc, Integer *npivm, Integer iwork[], const Integer *liwork, Integer *ifail, const Charlen length_mic, const Charlen length_pstrat)
The routine may be called by the names f11jnf or nagf_sparse_complex_herm_precon_ichol.

3 Description

f11jnf 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 f11jqf.
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 Section 9.
The sparse matrix A is represented in symmetric coordinate storage (SCS) format (see Section 2.1.2 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.  

4 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

5 Arguments

1: n Integer Input
On entry: n, the order of the matrix A.
Constraint: n1.
2: nnz Integer Input
On entry: the number of nonzero elements in the lower triangular part of the matrix A.
Constraint: 1nnzn×n+1/2.
3: ala Complex (Kind=nag_wp) array Input/Output
On entry: 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 routine f11zpf may be used to order the elements in this way.
On exit: the first nnz 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.
4: la Integer Input
On entry: the dimension of the arrays a, irow and icol as declared in the (sub)program from which f11jnf is called. These arrays must be of sufficient size to store both A (nnz elements) and C (nnzc elements).
Constraint: la2×nnz.
5: irowla Integer array Input/Output
6: icolla Integer array Input/Output
On entry: 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 f11zpf):
  • 1irowin and 1icoliirowi, for i=1,2,,nnz;
  • irowi-1<irowi or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nnz.
On exit: the row and column indices of the nonzero elements returned in a.
7: lfill Integer Input
On entry: if lfill0 its value is the maximum level of fill allowed in the decomposition (see Section 9.2). A negative value of lfill indicates that dtol will be used to control the fill instead.
8: dtol Real (Kind=nag_wp) Input
On entry: if lfill<0, dtol is used as a drop tolerance to control the fill-in (see Section 9.2); otherwise dtol is not referenced.
Constraint: if lfill<0, dtol0.0.
9: mic Character(1) Input
On entry: indicates whether or not the factorization should be modified to preserve row sums (see Section 9.3).
mic='M'
The factorization is modified.
mic='N'
The factorization is not modified.
Constraint: mic='M' or 'N'.
10: dscale Real (Kind=nag_wp) Input
On entry: 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 Section 9.3.
11: pstrat Character(1) Input
On entry: 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.
Suggested value: pstrat='M'.
Constraint: pstrat='N', 'M' or 'U'.
12: ipivn Integer array Input/Output
On entry: 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.
On exit: the pivot indices. If ipivi=j, the diagonal element in row j was used as the pivot at elimination stage i.
13: istrn+1 Integer array Output
On exit: 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.
14: nnzc Integer Output
On exit: the number of nonzero elements in the lower triangular matrix C.
15: npivm Integer Output
On exit: 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 f11jnf again with an increased value of either lfill or dscale. See also Sections 9.3 and 9.4.
16: iworkliwork Integer array Workspace
17: liwork Integer Input
On entry: the dimension of the array iwork as declared in the (sub)program from which f11jnf is called.
Constraints:
the minimum permissible value of liwork depends on lfill as follows:
  • if lfill0, liwork2×la-3×nnz+7×n+1;
  • otherwise liworkla-nnz+7×n+1.
18: ifail Integer Input/Output
On entry: ifail must be set to 0, -1 or 1. If you are unfamiliar with this argument you should refer to Section 4 in the Introduction to the NAG Library FL Interface for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1 or 1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this argument, the recommended value is 0. When the value -1 or 1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6 Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
ifail=1
On entry, dtol=value.
Constraint: dtol0.0
On entry, la=value and nnz=value.
Constraint: la2×nnz.
On entry, liwork is too small: liwork=value. Minimum required value of liwork=value.
On entry, mic'M' or 'N': mic=value.
On entry, n=value.
Constraint: n1.
On entry, nnz=value.
Constraint: nnz1.
On entry, nnz=value and n=value.
Constraint: nnzn×n+1/2
On entry, pstrat'N', 'U' or 'M': pstrat=value.
ifail=2
On entry, ai is out of order: i=value.
On entry, I=value, icolI=value and irowI=value.
Constraint: icolI1 and icolIirowI.
On entry, I=value, irowI=value and n=value.
Constraint: irowI1 and irowIn.
On entry, the location (irowI,icolI) is a duplicate: I=value.
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. Consider calling f11zpf to reorder and sum or remove duplicates.
ifail=3
On entry, a user-supplied value of ipiv is repeated.
On entry, a user-supplied value of ipiv lies outside the range 1,n.
ifail=4
The number of nonzero entries in the decomposition is too large. The decomposition has been terminated before completion. Either increase la, or reduce the fill by setting pstrat='M', reducing lfill, or increasing dtol.
ifail=5
A serious error has occurred in an internal call. Check all subroutine calls and array sizes. Seek expert help.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 7 in the Introduction to the NAG Library FL Interface for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library FL Interface for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.

7 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 f11jnf is used in combination with f11jqf, the more accurate the factorization the fewer iterations will be required. However, the cost of the decomposition will also generally increase.

8 Parallelism and Performance

f11jnf makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

9.1 Timing

The time taken for a call to f11jnf is roughly proportional to nnzc2/n.

9.2 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'.

9.3 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.
f11jnf 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 routine.
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'.

9.4 Direct Solution of positive definite Systems

Although it is not their primary purpose, f11jnf and f11jpf may be used together to obtain a direct solution to a complex Hermitian positive definite linear system. To achieve this the call to f11jpf should be preceded by a complete Cholesky factorization
A=PLDLHPT=M.  
A complete factorization is obtained from a call to f11jnf 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 f11jnf and f11jpf as a direct method is illustrated in f11jpf.

10 Example

This example reads in a complex sparse Hermitian matrix A and calls f11jnf to compute an incomplete Cholesky factorization. It then outputs the nonzero elements of both A and C=L+D-1-I.
The call to f11jnf 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.

10.1 Program Text

Program Text (f11jnfe.f90)

10.2 Program Data

Program Data (f11jnfe.d)

10.3 Program Results

Program Results (f11jnfe.r)