|
NAG[f11jac] NAG[nag_sparse_sym_chol_fac] - Incomplete Cholesky factorization (symmetric)
|
|
Calling Sequence
f11jac(nnz, a, la, irow, icol, lfill, dtol, mic, dscale, ipiv, istr, nnzc, npivm, comm, 'n'=n, 'pstrat'=pstrat, 'fail'=fail)
nag_sparse_sym_chol_fac(. . .)
Parameters
|
nnz - integer;
|
|
|
On entry: the number of non-zero elements in the lower triangular part of the matrix .
|
|
Constraint: . .
|
|
|
a - Vector(1..la, datatype=float[8]);
|
|
|
On entry: the non-zero elements in the lower triangular part of the matrix , 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 f11zbc (nag_sparse_sym_sort) may be used to order the elements in this way.
|
|
On exit: the first nnz elements of a contain the non-zero elements of and the next nnzc elements contain the elements of the lower triangular matrix . Matrix elements are ordered by increasing row index, and by increasing column index within each row.
|
|
|
la - assignable;
|
|
|
Note: On exit the variable la will have a value of type integer.
|
|
Default value: the dimension of the arrays a, irow, icol.
|
|
On entry: the first dimension of the arrays a, irow and icol as declared in the function from which nag_sparse_sym_chol_fac (f11jac) is called.
|
|
On exit: if internal allocation has taken place then la is set to , otherwise it remains unchanged.
|
|
Constraint: . .
|
|
|
irow - Vector(1..la, datatype=integer[kernelopts('wordsize')/8]);
icol - Vector(1..la, datatype=integer[kernelopts('wordsize')/8]);
|
|
|
On entry: the row and column indices of the non-zero elements supplied in a.
|
|
On exit: the row and column indices of the non-zero elements returned in a.
|
|
|
lfill - integer;
|
|
|
On entry: if its value is the maximum level of fill allowed in the decomposition (see Section [Control of Fill-in ]). A negative value of lfill indicates that dtol will be used to control the fill instead.
|
|
|
dtol - float;
|
|
|
On entry: if then dtol is used as a drop tolerance to control the fill-in (see Section [Control of Fill-in ]). Otherwise dtol is not referenced.
|
|
Constraint: if , . .
|
|
|
mic - String;
|
|
|
On entry: indicates whether or not the factorization should be modified to preserve row sums (see Section [Choice of Arguments ]).
|
|
The factorization is modified (MIC).
|
|
The factorization is not modified.
|
|
Constraint: "Nag_SparseSym_ModFact" or "Nag_SparseSym_UnModFact". .
|
|
|
dscale - float;
|
|
|
On entry: the diagonal scaling argument. All diagonal elements are multiplied by the factor at the start of the factorization. This can be used to ensure that the preconditioner is positive-definite. See Section [Choice of Arguments ].
|
|
|
ipiv - Vector(1..n, datatype=integer[kernelopts('wordsize')/8]);
|
|
|
Constraint: if , then ipiv must contain a valid permutation of the integers on .
|
|
|
istr - Vector(1.. , datatype=integer[kernelopts('wordsize')/8]);
|
|
|
|
nnzc - assignable;
|
|
|
Note: On exit the variable nnzc will have a value of type integer.
|
|
On exit: the number of non-zero elements in the lower triangular matrix .
|
|
|
npivm - assignable;
|
|
|
Note: On exit the variable npivm will have a value of type integer.
|
|
On exit: the number of pivots which were modified during the factorization to ensure that 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_sym_chol_fac (f11jac) again with an increased value of either lfill or dscale.
|
|
|
comm - Vector;
|
|
|
A Maple Vector, which should be generated using NAG[Nag_Sparse_Comm], corresponding to the Nag_Sparse_Comm structure.
|
|
|
'n'=n - integer; (optional)
|
|
|
Default value: the first dimension of the array ipiv.
|
|
On entry: the order of the matrix .
|
|
Constraint: . .
|
|
|
'pstrat'=pstrat - String; (optional)
|
|
|
On entry: specifies the pivoting strategy to be adopted as follows:
|
|
if then no pivoting is carried out;
|
|
if then diagonal pivoting aimed at minimizing fill-in is carried out, using the Markowitz strategy;
|
|
if then diagonal pivoting is carried out according to the user-defined input value of ipiv.
|
|
Suggested value: . (default: )
|
|
Constraint: "Nag_SparseSym_NoPiv", "Nag_SparseSym_MarkPiv" or "Nag_SparseSym_UserPiv". .
|
|
|
'fail'=fail - table; (optional)
|
|
|
The NAG error argument, see the documentation for NagError.
|
|
|
|
Description
|
|
|
Purpose
|
|
nag_sparse_sym_chol_fac (f11jac) computes an incomplete Cholesky factorization of a real sparse symmetric matrix, represented in symmetric coordinate storage format. This factorization may be used as a preconditioner in combination with f11jcc (nag_sparse_sym_chol_sol).
|
|
Description
|
|
This function computes an incomplete Cholesky factorization (see Meijerink and Van der Vorst (1977)) of a real sparse symmetric by matrix . 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 for the symmetric iterative solver f11jcc (nag_sparse_sym_chol_sol).
The decomposition is written in the form
where
and is a permutation matrix, is lower triangular with unit diagonal elements, is diagonal and 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 [Further Comments].
The sparse matrix is represented in symmetric coordinate storage (SCS) format (see Section [Symmetric coordinate storage (SCS) format] in the f11 Chapter Introduction). The array a stores all the non-zero elements of the lower triangular part of , while arrays irow and icol store the corresponding row and column indices respectively. Multiple non-zero elements may not be specified for the same row and column index.
The preconditioning matrix is returned in terms of the SCS representation of the lower triangular matrix
|
|
Error Indicators and Warnings
|
|
"NE_2_INT_ARG_LT"
On entry, while . These arguments must satisfy .
"NE_ALLOC_FAIL"
Dynamic memory allocation failed.
"NE_BAD_PARAM"
On entry, argument mic had an illegal value.
"NE_INT_2"
On entry, , . Constraint: .
"NE_INT_ARG_LT"
On entry, n must not be less than 1: .
"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 consult NAG for assistance.
"NE_INVALID_ROW_PIVOT"
On entry, and the array ipiv does not represent a valid permutation of integers in . An input value of ipiv is either out of range or repeated.
"NE_REAL_INT_ARG_CONS"
On entry, and . These arguments must satisfy if .,
"NE_SYMM_MATRIX_DUP"
A non-zero element has been supplied which does not lie in the lower triangular part of the matrix , is out of order, or has duplicate row and column indices, i.e., one or more of the following constraints has been violated:
|
, or
|
Call f11zbc (nag_sparse_sym_sort) to reorder and sum or remove duplicates.
|
|
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 . The factorization can generally be made more accurate by increasing lfill, or by reducing dtol with . If nag_sparse_sym_chol_fac (f11jac) is used in combination with f11jcc (nag_sparse_sym_chol_sol), the more accurate the factorization the fewer iterations will be required. However, the cost of the decomposition will also generally increase.
|
|
Further Comments
|
|
The time taken for a call to nag_sparse_sym_chol_fac (f11jac) is roughly proportional to .
|
Choice of Arguments
|
|
There is unfortunately no choice of the various algorithmic arguments which is optimal for all types of symmetric matrix, and some experimentation will generally be required for each new type of matrix encountered.
If the matrix is not known to have any particular special properties the following strategy is recommended. Start with , and . If the value returned for npivm is significantly larger than zero, i.e., a large number of pivot modifications were required to ensure that 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 to see if any improvement can be obtained by using modified incomplete Cholesky.
nag_sparse_sym_chol_fac (f11jac) 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 is probably too indefinite for this function.
If has non-positive off-diagonal elements, is non-singular, and has only non-negative elements in its inverse, it is called an "M-matrix". It can be shown that no pivot modifications are required in the incomplete Cholesky factorization of an M-matrix (Meijerink and Van der Vorst (1977)). In this case a good preconditioner can generally be expected by setting , and .
For certain mesh-based problems involving M-matrices it can be shown in theory that setting , and choosing dscale appropriately can reduce the order of magnitude of the condition number of the preconditioned matrix as a function of the mesh steplength (Chan (1991)). In practise this property often holds even with , although an improvement in condition can result from increasing dscale slightly (Van der Vorst (1990)).
Some illustrations of the application of nag_sparse_sym_chol_fac (f11jac) to linear systems arising from the discretization of two-dimensional elliptic partial differential equations, and to random-valued randomly structured symmetric positive-definite linear systems, can be found in Salvini and Shaw (1995).
|
|
|
|
Examples
|
|
>
|
n := 7:
nnz := 16:
la := 32:
lfill := 0:
dtol := 0.0:
mic := "Nag_SparseSym_UnModFact":
dscale := 0:
pstrat := "Nag_SparseSym_MarkPiv":
comm := NAG:-Nag_Sparse_Comm():
a := Vector(la,[4,1,5,2,2,3,-1,1,4,1,-2,3,2,-1,-2,5],datatype=float[8]):
irow := Vector(la,[1,2,2,3,4,4,5,5,5,6,6,6,7,7,7,7],datatype=integer[kernelopts('wordsize')/8]):
icol := Vector(la,[1,1,2,3,2,4,1,4,5,2,5,6,1,2,3,7],datatype=integer[kernelopts('wordsize')/8]):
ipiv := Vector(n, datatype=integer[kernelopts('wordsize')/8]):
istr := Vector(n+1, datatype=integer[kernelopts('wordsize')/8]):
NAG:-f11jac(nnz, a, la, irow, icol, lfill, dtol, mic, dscale, ipiv, istr, nnzc, npivm, comm, 'n' = n, 'pstrat' = pstrat):
|
|
|
See Also
|
|
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 (ed O Axelsson and L Y Kolotilina) 1457 Springer–Verlag
f11 Chapter Introduction.
NAG Toolbox Overview.
NAG Web Site.
|
|