naginterfaces.library.sparse.real_​symm_​solve_​ichol

naginterfaces.library.sparse.real_symm_solve_ichol(method, nnz, a, irow, icol, ipiv, istr, b, tol, maxitn, x)[source]

real_symm_solve_ichol solves a real sparse symmetric system of linear equations, represented in symmetric coordinate storage format, using a conjugate gradient or Lanczos method, with incomplete Cholesky preconditioning.

For full information please refer to the NAG Library document for f11jc

https://support.nag.com/numeric/nl/nagdoc_30.2/flhtml/f11/f11jcf.html

Parameters
methodstr

Specifies the iterative method to be used.

Conjugate gradient method.

Lanczos method (SYMMLQ).

nnzint

The number of nonzero elements in the lower triangular part of the matrix . This must be the same value as was supplied in the preceding call to real_symm_precon_ichol().

afloat, array-like, shape

The values returned in the array by a previous call to real_symm_precon_ichol().

irowint, array-like, shape

The values returned in arrays , , and by a previous call to real_symm_precon_ichol().

icolint, array-like, shape

The values returned in arrays , , and by a previous call to real_symm_precon_ichol().

ipivint, array-like, shape

The values returned in arrays , , and by a previous call to real_symm_precon_ichol().

istrint, array-like, shape

The values returned in arrays , , and by a previous call to real_symm_precon_ichol().

bfloat, array-like, shape

The right-hand side vector .

tolfloat

The required tolerance. Let denote the approximate solution at iteration , and the corresponding residual. The algorithm is considered to have converged at iteration if

If , is used, where is the machine precision. Otherwise is used.

maxitnint

The maximum number of iterations allowed.

xfloat, array-like, shape

An initial approximation to the solution vector .

Returns
xfloat, ndarray, shape

An improved approximation to the solution vector .

rnormfloat

The final value of the residual norm , where is the output value of .

itnint

The number of iterations carried out.

Raises
NagValueError
(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: or .

(errno )

On entry, and .

Constraint: .

(errno )

On entry, and .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, the location () is a duplicate: .

(errno )

On entry, is out of order: .

(errno )

On entry, , , .

Constraint: and .

(errno )

On entry, , and .

Constraint: and .

(errno )

The SCS representation of the preconditioner is invalid.

(errno )

The solution has not converged after iterations.

(errno )

The preconditioner appears not to be positive definite. The computation cannot continue.

(errno )

The matrix of the coefficients appears not to be positive definite. The computation cannot continue.

(errno )

A serious error has occurred in an internal call: . Check all function calls and array sizes. Seek expert help.

(errno )

A serious error has occurred in an internal call: . Check all function calls and array sizes. Seek expert help.

Warns
NagAlgorithmicWarning
(errno )

The required accuracy could not be obtained. However a reasonable accuracy has been achieved.

Notes

In the NAG Library the traditional C interface for this routine uses a different algorithmic base. Please contact NAG if you have any questions about compatibility.

real_symm_solve_ichol solves a real sparse symmetric linear system of equations

using a preconditioned conjugate gradient method (see Meijerink and Van der Vorst (1977)), or a preconditioned Lanczos method based on the algorithm SYMMLQ (see Paige and Saunders (1975)). The conjugate gradient method is more efficient if is positive definite, but may fail to converge for indefinite matrices. In this case the Lanczos method should be used instead. For further details see Barrett et al. (1994).

real_symm_solve_ichol uses the incomplete Cholesky factorization determined by real_symm_precon_ichol() as the preconditioning matrix. A call to real_symm_solve_ichol must always be preceded by a call to real_symm_precon_ichol(). Alternative preconditioners for the same storage scheme are available by calling real_symm_solve_jacssor().

The matrix , and the preconditioning matrix , are represented in symmetric coordinate storage (SCS) format (see the F11 Introduction) in the arrays , and , as returned from real_symm_precon_ichol(). The array holds the nonzero entries in the lower triangular parts of these matrices, while and hold the corresponding row and column indices.

References

Barrett, R, Berry, M, Chan, T F, Demmel, J, Donato, J, Dongarra, J, Eijkhout, V, Pozo, R, Romine, C and Van der Vorst, H, 1994, Templates for the Solution of Linear Systems: Building Blocks for Iterative Methods, SIAM, Philadelphia

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

Paige, C C and Saunders, M A, 1975, Solution of sparse indefinite systems of linear equations, SIAM J. Numer. Anal. (12), 617–629

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