naginterfaces.library.sparse.real_​symm_​basic_​setup

naginterfaces.library.sparse.real_symm_basic_setup(method, precon, n, tol, maxitn, anorm, sigmax, maxits, monit, sigcmp='N', norm=None, weight='N', iterm=1, sigtol=0.01)[source]

real_symm_basic_setup is a setup function, the first in a suite of three functions for the iterative solution of a symmetric system of simultaneous linear equations. real_symm_basic_setup must be called before the iterative solver, real_symm_basic_solver(). The third function in the suite, real_symm_basic_diag(), can be used to return additional information about the computation.

These three functions are suitable for the solution of large sparse symmetric systems of equations.

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

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

Parameters
methodstr

The iterative method to be used.

Conjugate gradient method (CG).

Lanczos method (SYMMLQ).

Minimum residual method (MINRES).

preconstr, length 1

Determines whether preconditioning is used.

No preconditioning.

Preconditioning.

nint

, the order of the matrix .

tolfloat

The tolerance for the termination criterion.

If , is used, where is the machine precision.

Otherwise is used.

maxitnint

The maximum number of iterations.

anormfloat

If , the value of to be used in the termination criterion (2) ().

If , and or , is estimated internally by real_symm_basic_solver().

If , is not referenced.

It has no effect if .

sigmaxfloat

If , the value of .

If , is estimated by real_symm_basic_solver() when either or termination criterion (3) () is employed, though it will be used only in the latter case.

Otherwise, or if , is not referenced.

maxitsint

The maximum iteration number for which is computed by bisection (see also Notes). If or , or if , is not referenced.

Suggested value: when is of the order of its default value .

monitint

If , the frequency at which a monitoring step is executed by real_symm_basic_solver(): the current solution and residual iterates will be returned by real_symm_basic_solver() and a call to real_symm_basic_diag() made possible every iterations, starting from the ()th. Otherwise, no monitoring takes place.

There are some additional computational costs involved in monitoring the solution and residual vectors when the Lanczos method (SYMMLQ) is used.

sigcmpstr, length 1, optional

Determines whether an estimate of , the largest singular value of the preconditioned matrix of the coefficients, is to be computed using the bisection method on the sequence of tridiagonal matrices generated during the iteration. Note that when a preconditioner is not used.

If (see ), i.e., when is supplied, the value of is ignored.

is to be computed using the bisection method.

The bisection method is not used.

If the termination criterion (3) is used, requiring , an inexpensive estimate is computed and used (see Notes).

It is not used if .

normNone or str, length 1, optional

Note: if this argument is None then a default value will be used, determined as follows: if : ; otherwise: .

If or , defines the matrix and vector norm to be used in the termination criteria.

Use the norm.

Use the norm.

Use the norm.

It has no effect if .

weightstr, length 1, optional

Specifies whether a vector of user-supplied weights is to be used in the vector norms used in the computation of termination criterion (2) (): , where , for . The suffix denotes the vector norm used, as specified by the argument . Note that weights cannot be used when , i.e., when criterion (3) is used.

User-supplied weights are to be used and must be supplied on initial entry to real_symm_basic_solver().

All weights are implicitly set equal to one. Weights do not need to be supplied on initial entry to real_symm_basic_solver().

It has no effect if .

itermint, optional

Defines the termination criterion to be used.

Use the termination criterion defined in (2) (both conjugate gradient and Lanczos (SYMMLQ) methods).

Use the termination criterion defined in (3) (Lanczos method (SYMMLQ) only).

It has no effect if .

sigtolfloat, optional

The tolerance used in assessing the convergence of the estimate of when the bisection method is used.

If , the default value is used.

The actual value used is .

If or , is not referenced.

It has no effect if .

Returns
commdict, communication object

Communication structure.

Raises
NagValueError
(errno )

On entry, and .

Constraint: .

(errno )

On entry, , , and .

Constraint: if and , .

(errno )

On entry, , and .

Constraint: if and , .

(errno )

On entry, , and .

Constraint: if and , .

(errno )

On entry, , and .

Constraint: if and , .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, and .

Constraint: if , .

(errno )

On entry, and .

Constraint: if , .

(errno )

On entry, and .

Constraint: if , .

(errno )

On entry, .

Constraint: or .

(errno )

On entry, .

Constraint: or .

(errno )

On entry, .

Constraint: , or .

(errno )

On entry, .

Constraint: or .

(errno )

On entry, .

Constraint: or .

(errno )

On entry, .

Constraint: , or .

(errno )

real_symm_basic_setup has been called out of sequence: either real_symm_basic_setup has been called twice or real_symm_basic_solver() has not terminated its current task.

Notes

The suite consisting of the functions real_symm_basic_setup, real_symm_basic_solver() and real_symm_basic_diag() is designed to solve the symmetric system of simultaneous linear equations of order , where is large and the matrix of the coefficients is sparse.

real_symm_basic_setup is a setup function which must be called before real_symm_basic_solver(), the iterative solver. The third function in the suite, real_symm_basic_diag() can be used to return additional information about the computation. One of the following methods can be used:

  1. Conjugate Gradient Method (CG)

    For this method (see Hestenes and Stiefel (1952), Golub and Van Loan (1996), Barrett et al. (1994) and Dias da Cunha and Hopkins (1994)), the matrix should ideally be positive definite. The application of the Conjugate Gradient method to indefinite matrices may lead to failure or to lack of convergence.

  2. Lanczos Method (SYMMLQ)

    This method, based upon the algorithm SYMMLQ (see Paige and Saunders (1975) and Barrett et al. (1994)), is suitable for both positive definite and indefinite matrices. It is more robust than the Conjugate Gradient method but less efficient when is positive definite.

  3. Minimum Residual Method (MINRES)

    This method may be used when the matrix is indefinite. It seeks to reduce the norm of the residual at each iteration and often takes fewer iterations than the other methods. It does however require slightly more memory.

The CG and SYMMLQ methods start from the residual , where is an initial estimate for the solution (often ), and generate an orthogonal basis for the Krylov subspace , for , by means of three-term recurrence relations (see Golub and Van Loan (1996)). A sequence of symmetric tridiagonal matrices is also generated. Here and in the following, the index denotes the iteration count. The resulting symmetric tridiagonal systems of equations are usually more easily solved than the original problem. A sequence of solution iterates is thus generated such that the sequence of the norms of the residuals converges to a required tolerance. Note that, in general, the convergence is not monotonic.

In exact arithmetic, after iterations, this process is equivalent to an orthogonal reduction of to symmetric tridiagonal form, ; the solution would thus achieve exact convergence. In finite-precision arithmetic, cancellation and round-off errors accumulate causing loss of orthogonality. These methods must, therefore, be viewed as genuinely iterative methods, able to converge to a solution within a prescribed tolerance.

The orthogonal basis is not formed explicitly in either method. The basic difference between the Conjugate Gradient and Lanczos methods lies in the method of solution of the resulting symmetric tridiagonal systems of equations: the conjugate gradient method is equivalent to carrying out an (Cholesky) factorization whereas the Lanczos method (SYMMLQ) uses an factorization.

Faster convergence for all the methods can be achieved using a preconditioner (see Golub and Van Loan (1996) and Barrett et al. (1994)). A preconditioner maps the original system of equations onto a different system, say

with, hopefully, better characteristics with respect to its speed of convergence: for example, the condition number of the matrix of the coefficients can be improved or eigenvalues in its spectrum can be made to coalesce. An orthogonal basis for the Krylov subspace , for , is generated and the solution proceeds as outlined above. The algorithms used are such that the solution and residual iterates of the original system are produced, not their preconditioned counterparts. Note that an unsuitable preconditioner or no preconditioning at all may result in a very slow rate, or lack, of convergence. However, preconditioning involves a trade-off between the reduction in the number of iterations required for convergence and the additional computational costs per iteration. Also, setting up a preconditioner may involve non-negligible overheads.

A preconditioner must be symmetric and positive definite, i.e., representable by , where is nonsingular, and such that in (1), where is the identity matrix of order . Also, we can define and . These are formal definitions, used only in the design of the algorithms; in practice, only the means to compute the matrix-vector products and to solve the preconditioning equations are required, that is, explicit information about , or their inverses is not required at any stage.

The first termination criterion

is available for both conjugate gradient and Lanczos (SYMMLQ) methods. In (2), or and denotes a user-specified tolerance subject to , where is the machine precision. Facilities are provided for the estimation of the norm of the matrix of the coefficients , when this is not known in advance, used in (2), by applying Higham’s method (see Higham (1988)). Note that cannot be estimated internally. This criterion uses an error bound derived from backward error analysis to ensure that the computed solution is the exact solution of a problem as close to the original as the termination tolerance requires. Termination criteria employing bounds derived from forward error analysis could be used, but any such criteria would require information about the condition number which is not easily obtainable.

The second termination criterion

is available only for the Lanczos method (SYMMLQ). In (3), is the largest singular value of the (preconditioned) iteration matrix . This termination criterion monitors the progress of the solution of the preconditioned system of equations and is less expensive to apply than criterion (2). When is not supplied, facilities are provided for its estimation by . The interlacing property and Gerschgorin’s theorem provide lower and upper bounds from which can be easily computed by bisection. Alternatively, the less expensive estimate can be used, where by Gerschgorin’s theorem. Note that only order of magnitude estimates are required by the termination criterion.

Termination criterion (2) is the recommended choice, despite its (small) additional costs per iteration when using the Lanczos method (SYMMLQ). Also, if the norm of the initial estimate is much larger than the norm of the solution, that is, if , a dramatic loss of significant digits could result in complete lack of convergence. The use of criterion (2) will enable the detection of such a situation, and the iteration will be restarted at a suitable point. No such restart facilities are provided for criterion (3).

Optionally, a vector of user-specified weights can be used in the computation of the vector norms in termination criterion (2), i.e., , where , for . Note that the use of weights increases the computational costs.

The MINRES algorithm terminates when the norm of the residual of the preconditioned system , , where is the preconditioned matrix.

The termination criteria discussed are not robust in the presence of a non-trivial nullspace of , i.e., when is singular. It is then possible for to grow without limit, spuriously satisfying the termination criterion. If singularity is suspected, more robust functions can be found in submodule opt.

The sequence of calls to the functions comprising the suite is enforced: first, the setup function real_symm_basic_setup must be called, followed by the solver real_symm_basic_solver(). The diagnostic function real_symm_basic_diag() can be called either when real_symm_basic_solver() is carrying out a monitoring step or after real_symm_basic_solver() has completed its tasks. Incorrect sequencing will raise an error condition.

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

Dias da Cunha, R and Hopkins, T, 1994, PIM 1.1 — the parallel iterative method package for systems of linear equations user’s guide — Fortran 77 version, Technical Report, Computing Laboratory, University of Kent at Canterbury, Kent, UK

Golub, G H and Van Loan, C F, 1996, Matrix Computations, (3rd Edition), Johns Hopkins University Press, Baltimore

Hestenes, M and Stiefel, E, 1952, Methods of conjugate gradients for solving linear systems, J. Res. Nat. Bur. Stand. (49), 409–436

Higham, N J, 1988, FORTRAN codes for estimating the one-norm of a real or complex matrix, with applications to condition estimation, ACM Trans. Math. Software (14), 381–396

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