NAG FL Interface
f11gdf (real_​symm_​basic_​setup)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f11gdf is a setup routine, the first in a suite of three routines for the iterative solution of a symmetric system of simultaneous linear equations. f11gdf must be called before the iterative solver, f11gef. The third routine in the suite, f11gff, can be used to return additional information about the computation.
These three routines are suitable for the solution of large sparse symmetric systems of equations.

2 Specification

Fortran Interface
Integer, Intent (In) :: iterm, n, maxitn, maxits, monit, lwork
Integer, Intent (Inout) :: ifail
Integer, Intent (Out) :: lwreq
Real (Kind=nag_wp), Intent (In) :: tol, anorm, sigmax, sigtol
Real (Kind=nag_wp), Intent (Out) :: work(lwork)
Character (*), Intent (In) :: method
Character (1), Intent (In) :: precon, sigcmp, norm, weight
C Header Interface
#include <nag.h>
void  f11gdf_ (const char *method, const char *precon, const char *sigcmp, const char *norm, const char *weight, const Integer *iterm, const Integer *n, const double *tol, const Integer *maxitn, const double *anorm, const double *sigmax, const double *sigtol, const Integer *maxits, const Integer *monit, Integer *lwreq, double work[], const Integer *lwork, Integer *ifail, const Charlen length_method, const Charlen length_precon, const Charlen length_sigcmp, const Charlen length_norm, const Charlen length_weight)
The routine may be called by the names f11gdf or nagf_sparse_real_symm_basic_setup.

3 Description

The suite consisting of the routines f11gdf, f11gef and f11gff is designed to solve the symmetric system of simultaneous linear equations Ax=b of order n, where n is large and the matrix of the coefficients A is sparse.
f11gdf is a setup routine which must be called before f11gef, the iterative solver. The third routine in the suite, f11gff can be used to return additional information about the computation. One of the following methods can be used:
  1. 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 A 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. 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 A is positive definite.
  3. 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 r0=b-Ax0, where x0 is an initial estimate for the solution (often x0=0), and generate an orthogonal basis for the Krylov subspace span{Akr0}, for k=0,1,, by means of three-term recurrence relations (see Golub and Van Loan (1996)). A sequence of symmetric tridiagonal matrices {Tk} is also generated. Here and in the following, the index k 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 {xk} is thus generated such that the sequence of the norms of the residuals {rk} converges to a required tolerance. Note that, in general, the convergence is not monotonic.
In exact arithmetic, after n iterations, this process is equivalent to an orthogonal reduction of A to symmetric tridiagonal form, Tn=QTAQ; the solution xn 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 LDLT (Cholesky) factorization whereas the Lanczos method (SYMMLQ) uses an LQ 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
A¯x¯=b¯, (1)
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 span{A¯kr¯0}, for k=0,1,, 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 M=EET, where M is nonsingular, and such that A¯=E-1AE-TIn in (1), where In is the identity matrix of order n. Also, we can define r¯=E-1r and x¯=ETx. These are formal definitions, used only in the design of the algorithms; in practice, only the means to compute the matrix-vector products v=Au and to solve the preconditioning equations Mv=u are required, that is, explicit information about M, E or their inverses is not required at any stage.
The first termination criterion
rkp τ ( bp + Ap ×xkp) (2)
is available for both conjugate gradient and Lanczos (SYMMLQ) methods. In (2), p=1, or 2 and τ denotes a user-specified tolerance subject to max(10,n)ετ<1, where ε is the machine precision. Facilities are provided for the estimation of the norm of the matrix of the coefficients A1=A, when this is not known in advance, used in (2), by applying Higham's method (see Higham (1988)). Note that A2 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 κ(A) which is not easily obtainable.
The second termination criterion
r¯k2 τ max(1.0, b2 / r02 ) ( r¯02 + σ1 (A¯) ×Δx¯k2) (3)
is available only for the Lanczos method (SYMMLQ). In (3), σ1(A¯)=A¯2 is the largest singular value of the (preconditioned) iteration matrix A¯. 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 σ1(A¯) is not supplied, facilities are provided for its estimation by σ1(A¯)maxkσ1(Tk). The interlacing property σ1(Tk-1)σ1(Tk) and Gerschgorin's theorem provide lower and upper bounds from which σ1(Tk) can be easily computed by bisection. Alternatively, the less expensive estimate σ1(A¯)maxkTk1 can be used, where σ1(A¯)Tk1 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 x0x, 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 w of user-specified weights can be used in the computation of the vector norms in termination criterion (2), i.e., vp(w)=v (w) p, where (v (w) )i=wi vi, for i=1,2,,n. Note that the use of weights increases the computational costs.
The MINRES algorithm terminates when the norm of the residual of the preconditioned system F, F2τ× A¯2× xk2 , where A¯ is the preconditioned matrix.
The termination criteria discussed are not robust in the presence of a non-trivial nullspace of A, i.e., when A is singular. It is then possible for xkp to grow without limit, spuriously satisfying the termination criterion. If singularity is suspected, more robust routines can be found in Chapter E04.
The sequence of calls to the routines comprising the suite is enforced: first, the setup routine f11gdf must be called, followed by the solver f11gef. The diagnostic routine f11gff can be called either when f11gef is carrying out a monitoring step or after f11gef has completed its tasks. Incorrect sequencing will raise an error condition.

4 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

5 Arguments

1: method Character(*) Input
On entry: the iterative method to be used.
method='CG'
Conjugate gradient method (CG).
method='SYMMLQ'
Lanczos method (SYMMLQ).
method='MINRES'
Minimum residual method (MINRES).
Constraint: method='CG', 'SYMMLQ' or 'MINRES'.
2: precon Character(1) Input
On entry: determines whether preconditioning is used.
precon='N'
No preconditioning.
precon='P'
Preconditioning.
Constraint: precon='N' or 'P'.
3: sigcmp Character(1) Input
On entry: determines whether an estimate of σ1(A¯)=E-1AE-T2, 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 {Tk} generated during the iteration. Note that A¯=A when a preconditioner is not used.
If sigmax>0.0 (see sigmax), i.e., when σ1(A¯) is supplied, the value of sigcmp is ignored.
sigcmp='S'
σ1(A¯) is to be computed using the bisection method.
sigcmp='N'
The bisection method is not used.
If the termination criterion (3) is used, requiring σ1(A¯), an inexpensive estimate is computed and used (see Section 3).
It is not used if method='MINRES'.
Suggested value: sigcmp='N'.
Constraint: sigcmp='S' or 'N'.
4: norm Character(1) Input
On entry: if method='CG' or 'SYMMLQ', norm defines the matrix and vector norm to be used in the termination criteria.
norm='1'
Use the l1 norm.
norm='I'
Use the l norm.
norm='2'
Use the l2 norm.
It has no effect if method='MINRES'.
Suggested values:
  • if iterm=1, norm='I';
  • if iterm=2, norm='2'.
Constraints:
  • if iterm=1, norm='1', 'I' or '2';
  • if iterm=2, norm='2'.
5: weight Character(1) Input
On entry: specifies whether a vector w of user-supplied weights is to be used in the vector norms used in the computation of termination criterion (2) (iterm=1): vp(w)=v (w) p, where vi (w) =wi vi, for i=1,2,,n. The suffix p=1,2, denotes the vector norm used, as specified by the argument norm. Note that weights cannot be used when iterm=2, i.e., when criterion (3) is used.
weight='W'
User-supplied weights are to be used and must be supplied on initial entry to f11gef.
weight='N'
All weights are implicitly set equal to one. Weights do not need to be supplied on initial entry to f11gef.
It has no effect if method='MINRES'.
Suggested value: weight='N'.
Constraints:
  • if iterm=1, weight='W' or 'N';
  • if iterm=2, weight='N'.
6: iterm Integer Input
On entry: defines the termination criterion to be used.
iterm=1
Use the termination criterion defined in (2) (both conjugate gradient and Lanczos (SYMMLQ) methods).
iterm=2
Use the termination criterion defined in (3) (Lanczos method (SYMMLQ) only).
It has no effect if method='MINRES'.
Suggested value: iterm=1.
Constraints:
  • if method='CG', iterm=1;
  • if method='SYMMLQ', iterm=1 or 2.
7: n Integer Input
On entry: n, the order of the matrix A.
Constraint: n>0.
8: tol Real (Kind=nag_wp) Input
On entry: the tolerance τ for the termination criterion.
If tol0.0, τ=max(ε,nε) is used, where ε is the machine precision.
Otherwise τ=max(tol,10ε,nε) is used.
Constraint: tol<1.0.
9: maxitn Integer Input
On entry: the maximum number of iterations.
Constraint: maxitn>0.
10: anorm Real (Kind=nag_wp) Input
On entry: if anorm>0.0, the value of Ap to be used in the termination criterion (2) (iterm=1).
If anorm0.0, iterm=1 and norm='1' or 'I', A1=A is estimated internally by f11gef.
If iterm=2, anorm is not referenced.
It has no effect if method='MINRES'.
Constraint: if iterm=1 and norm=2, anorm>0.0.
11: sigmax Real (Kind=nag_wp) Input
On entry: if sigmax>0.0, the value of σ1(A¯)=E-1AE-T2.
If sigmax0.0, σ1(A¯) is estimated by f11gef when either sigcmp='S' or termination criterion (3) (iterm=2) is employed, though it will be used only in the latter case.
Otherwise, or if method='MINRES', sigmax is not referenced.
12: sigtol Real (Kind=nag_wp) Input
On entry: the tolerance used in assessing the convergence of the estimate of σ1(A¯)=A¯2 when the bisection method is used.
If sigtol0.0, the default value sigtol=0.01 is used. The actual value used is max(sigtol,ε).
If sigcmp='N' or sigmax>0.0, sigtol is not referenced.
It has no effect if method='MINRES'.
Suggested value: sigtol=0.01 should be sufficient in most cases.
Constraint: if sigcmp='S' and sigmax0.0, sigtol<1.0.
13: maxits Integer Input
On entry: the maximum iteration number k=maxits for which σ1(Tk) is computed by bisection (see also Section 3). If sigcmp='N' or sigmax>0.0, or if method='MINRES', maxits is not referenced.
Suggested value: maxits=min(10,n) when sigtol is of the order of its default value (0.01).
Constraint: if sigcmp='S' and sigmax0.0, 1maxitsmaxitn.
14: monit Integer Input
On entry: if monit>0, the frequency at which a monitoring step is executed by f11gef: the current solution and residual iterates will be returned by f11gef and a call to f11gff made possible every monit iterations, starting from the (monit)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.
Constraint: monitmaxitn.
15: lwreq Integer Output
On exit: the minimum amount of workspace required by f11gef. (See also Section 5 in f11gef.)
16: work(lwork) Real (Kind=nag_wp) array Communication Array
On exit: the array work is initialized by f11gdf. It must not be modified before calling the next routine in the suite, namely f11gef.
17: lwork Integer Input
On entry: the dimension of the array work as declared in the (sub)program from which f11gdf is called.
Constraint: lwork120.
Note:  although the minimum value of lwork ensures the correct functioning of f11gdf, a larger value is required by the other routines in the suite, namely f11gef and f11gff. The required value is as follows:
Method Requirements
CG lwork=120+5n+p.
SYMMLQ lwork=120+6n+p,
MINRES lwork=120+9n,
where
  • p=2*(maxits+1), when an estimate of σ1(A) (sigmax) is computed;
  • p=0, otherwise.
18: ifail Integer Input/Output
On entry: ifail must be set to 0, −1 or 1 to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of 0 causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of −1 means that an error message is printed while a value of 1 means that it is not.
If halting is not appropriate, the value −1 or 1 is recommended. If message printing is undesirable, then the value 1 is recommended. Otherwise, the value 0 is recommended. 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
f11gdf has been called out of sequence: either f11gdf has been called twice or f11gef has not terminated its current task.
ifail=-1
On entry, method=value.
Constraint: method='CG', 'SYMMLQ' or 'MINRES'.
ifail=-2
On entry, precon=value.
Constraint: precon='N' or 'P'.
ifail=-3
On entry, sigcmp=value.
Constraint: sigcmp='S' or 'N'.
ifail=-4
On entry, norm=value.
Constraint: norm='1', 'I' or '2'.
ifail=-5
On entry, weight=value.
Constraint: weight='N' or 'W'.
ifail=-6
On entry, iterm=value.
Constraint: iterm=1 or 2.
On entry, iterm=2 and method='CG'.
Constraint: if iterm=2, method'CG'.
On entry, iterm=2 and norm=value.
Constraint: if iterm=2, norm='2'.
On entry, iterm=2 and weight=value.
Constraint: if iterm=2, weight='N'.
ifail=-7
On entry, n=value.
Constraint: n>0.
ifail=-8
On entry, tol=value.
Constraint: tol<1.0.
ifail=-9
On entry, maxitn=value.
Constraint: maxitn>0.
ifail=-10
On entry, iterm=1, norm='2' and anorm=value.
Constraint: if iterm=1 and norm='2', anorm>0.0.
ifail=-12
On entry, sigcmp='S', sigmax0.0 and sigtol=value.
Constraint: if sigcmp='S' and sigmax0.0, sigtol<1.0.
ifail=-13
On entry, sigcmp='S', sigmax0.0, maxits=value and maxitn=value.
Constraint: if sigcmp='S' and sigmax0.0, maxitsmaxitn.
On entry, sigcmp='S', sigmax0.0 and maxits=value.
Constraint: if sigcmp='S' and sigmax0.0, maxits1.
ifail=-14
On entry, monit=value and maxitn=value.
Constraint: monitmaxitn.
ifail=-17
On entry, lwork=value.
Constraint: lwork120.
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

Not applicable.

8 Parallelism and Performance

Background information to multithreading can be found in the Multithreading documentation.
f11gdf is not threaded in any implementation.

9 Further Comments

When σ1(A¯) is not supplied (sigmax0.0) but it is required, it is estimated by f11gef using either of the two methods described in Section 3, as specified by the argument sigcmp. In particular, if sigcmp='S', then the computation of σ1(A¯) is deemed to have converged when the differences between three successive values of σ1(Tk) differ, in a relative sense, by less than the tolerance sigtol, i.e., when
max( |σ1(k)-σ1(k-1)| σ1(k) , ​ |σ1(k)-σ1(k-2)| σ1(k) ) sigtol.  
The computation of σ1(A¯) is also terminated when the iteration count exceeds the maximum value allowed, i.e., kmaxits.
Bisection is increasingly expensive with increasing iteration count. A reasonably large value of sigtol, of the order of the suggested value, is recommended and an excessive value of maxits should be avoided. Under these conditions, σ1(A¯) usually converges within very few iterations.

10 Example

This example solves a symmetric system of simultaneous linear equations using the conjugate gradient method, where the matrix of the coefficients A, has a random sparsity pattern. An incomplete Cholesky preconditioner is used (f11jaf and f11jbf).

10.1 Program Text

Program Text (f11gdfe.f90)

10.2 Program Data

Program Data (f11gdfe.d)

10.3 Program Results

Program Results (f11gdfe.r)