NAG Library Routine Document

f12fgf  (real_symm_band_solve)

Note: this routine uses optional parameters to define choices in the problem specification. If you wish to use default settings for all of the optional parameters, then the option setting routine f12fdf need not be called. If, however, you wish to reset some or all of the settings please refer to Section 11 in f12fdf for a detailed description of the specification of the optional parameters.

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

f12fgf is the main solver routine in a suite of routines which includes f12fdf and f12fff. f12fgf must be called following an initial call to f12fff and following any calls to f12fdf.
f12fgf returns approximations to selected eigenvalues, and (optionally) the corresponding eigenvectors, of a standard or generalized eigenvalue problem defined by real banded symmetric matrices. The banded matrix must be stored using the LAPACK storage format for real banded nonsymmetric matrices.

2
Specification

Fortran Interface
Subroutine f12fgf ( kl, ku, ab, ldab, mb, ldmb, sigma, nconv, d, z, ldz, resid, v, ldv, comm, icomm, ifail)
Integer, Intent (In):: kl, ku, ldab, ldmb, ldz, ldv
Integer, Intent (Inout):: icomm(*), ifail
Integer, Intent (Out):: nconv
Real (Kind=nag_wp), Intent (In):: ab(ldab,*), mb(ldmb,*), sigma
Real (Kind=nag_wp), Intent (Inout):: d(*), z(ldz,*), resid(*), v(ldv,*), comm(*)
C Header Interface
#include nagmk26.h
void  f12fgf_ ( const Integer *kl, const Integer *ku, const double ab[], const Integer *ldab, const double mb[], const Integer *ldmb, const double *sigma, Integer *nconv, double d[], double z[], const Integer *ldz, double resid[], double v[], const Integer *ldv, double comm[], Integer icomm[], Integer *ifail)

3
Description

The suite of routines is designed to calculate some of the eigenvalues, λ , (and optionally the corresponding eigenvectors, x ) of a standard eigenvalue problem Ax = λx , or of a generalized eigenvalue problem Ax = λBx  of order n , where n  is large and the coefficient matrices A  and B  are banded, real and symmetric.
Following a call to the initialization routine f12fff, f12fgf returns the converged approximations to eigenvalues and (optionally) the corresponding approximate eigenvectors and/or an orthonormal basis for the associated approximate invariant subspace. The eigenvalues (and eigenvectors) are selected from those of a standard or generalized eigenvalue problem defined by real banded symmetric matrices. There is negligible additional computational cost to obtain eigenvectors; an orthonormal basis is always computed, but there is an additional storage cost if both are requested.
The banded matrices A  and B must be stored using the LAPACK storage format for banded nonsymmetric matrices; please refer to Section 3.3.2 in the F07 Chapter Introduction for details on this storage format.
f12fgf is based on the banded driver routines dsbdr1 to dsbdr6 from the ARPACK package, which uses the Implicitly Restarted Lanczos iteration method. The method is described in Lehoucq and Sorensen (1996) and Lehoucq (2001) while its use within the ARPACK software is described in great detail in Lehoucq et al. (1998). This suite of routines offers the same functionality as the ARPACK banded driver software for real symmetric problems, but the interface design is quite different in order to make the option setting clearer and to combine the different drivers into a general purpose routine.
f12fgf, is a general purpose direct communication routine that must be called following initialization by f12fff. f12fgf uses options, set either by default or explicitly by calling f12fdf, to return the converged approximations to selected eigenvalues and (optionally):
the corresponding approximate eigenvectors;
an orthonormal basis for the associated approximate invariant subspace;
both.

4
References

Lehoucq R B (2001) Implicitly restarted Arnoldi methods and subspace iteration SIAM Journal on Matrix Analysis and Applications 23 551–562
Lehoucq R B and Scott J A (1996) An evaluation of software for computing eigenvalues of sparse nonsymmetric matrices Preprint MCS-P547-1195 Argonne National Laboratory
Lehoucq R B and Sorensen D C (1996) Deflation techniques for an implicitly restarted Arnoldi iteration SIAM Journal on Matrix Analysis and Applications 17 789–821
Lehoucq R B, Sorensen D C and Yang C (1998) ARPACK Users' Guide: Solution of Large-scale Eigenvalue Problems with Implicitly Restarted Arnoldi Methods SIAM, Philidelphia

5
Arguments

1:     kl – IntegerInput
On entry: the number of subdiagonals of the matrices A  and B.
Constraint: kl0.
2:     ku – IntegerInput
On entry: the number of superdiagonals of the matrices A  and B. Since A and B are symmetric, the normal case is ku=kl.
Constraint: ku0.
3:     abldab* – Real (Kind=nag_wp) arrayInput
Note: the second dimension of the array ab must be at least max1,n  (see f12fff).
On entry: must contain the matrix A in LAPACK banded storage format for nonsymmetric matrices (see Section 3.3.4 in the F07 Chapter Introduction).
4:     ldab – IntegerInput
On entry: the first dimension of the array ab as declared in the (sub)program from which f12fgf is called.
Constraint: ldab 2×kl+ku+1.
5:     mbldmb* – Real (Kind=nag_wp) arrayInput
Note: the second dimension of the array mb must be at least max1,n (see f12fff).
On entry: must contain the matrix B in LAPACK banded storage format for nonsymmetric matrices (see Section 3.3.4 in the F07 Chapter Introduction).
6:     ldmb – IntegerInput
On entry: the first dimension of the array mb as declared in the (sub)program from which f12fgf is called.
Constraint: ldmb2×kl+ku+1.
7:     sigma – Real (Kind=nag_wp)Input
On entry: if one of the Shifted Inverse (see f12fdf) modes has been selected then sigma contains the real shift used; otherwise sigma is not referenced.
8:     nconv – IntegerOutput
On exit: the number of converged eigenvalues.
9:     d* – Real (Kind=nag_wp) arrayOutput
Note: the dimension of the array d must be at least ncv (see f12fff).
On exit: the first nconv locations of the array d contain the converged approximate eigenvalues.
10:   zldz* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array z must be at least nev+1 if the default option Vectors=RITZ has been selected and at least 1 if the option Vectors=NONE or SCHUR has been selected (see f12fff).
On exit: if the default option Vectors=RITZ (see f12fdf) has been selected then z contains the final set of eigenvectors corresponding to the eigenvalues held in d. The real eigenvector associated with eigenvalue i-1, for i=1,2,,nconv, is stored in the ith column of z.
11:   ldz – IntegerInput
On entry: the first dimension of the array z as declared in the (sub)program from which f12fgf is called.
Constraints:
  • if the default option Vectors=Ritz has been selected, ldzn;
  • if the option Vectors=None or Schur has been selected, ldz1.
12:   resid* – Real (Kind=nag_wp) arrayInput/Output
Note: the dimension of the array resid must be at least n (see f12fff).
On entry: need not be set unless the option Initial Residual has been set in a prior call to f12fdf in which case resid must contain an initial residual vector.
On exit: contains the final residual vector.
13:   vldv* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array v must be at least max1,ncv (see f12fff).
On exit: if the option Vectors (see f12fdf) has been set to Schur or Ritz and a separate array z has been passed then the first nconv×n elements of v will contain approximate Schur vectors that span the desired invariant subspace.
The jth Schur vector is stored in the ith column of v.
14:   ldv – IntegerInput
On entry: the first dimension of the array v as declared in the (sub)program from which f12fgf is called.
Constraint: ldvn.
15:   comm* – Real (Kind=nag_wp) arrayCommunication Array
Note: the dimension of the array comm must be at least max1,lcomm (see f12fff).
On initial entry: must remain unchanged from the prior call to f12fdf and f12fff.
On exit: contains no useful information.
16:   icomm* – Integer arrayCommunication Array
Note: the dimension of the array icomm must be at least max1,licomm (see f12fff).
On initial entry: must remain unchanged from the prior call to f12fbf and f12fdf.
On exit: contains no useful information.
17:   ifail – IntegerInput/Output
On entry: ifail must be set to 0, -1​ or ​1. If you are unfamiliar with this argument you should refer to Section 3.4 in How to Use the NAG Library and its Documentation 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, kl < 0 .
ifail=2
On entry, ku < 0 .
ifail=3
On entry, ldab < 2×kl+ku+1 .
ifail=4
Iteration Limit < 0 .
ifail=5
The options Generalized and Regular are incompatible.
ifail=6
Eigenvalues from Both Ends of the spectrum were requested, but only one eigenvalue (nev) is requested.
ifail=7
The Initial Residual was selected but the starting vector held in resid is zero.
ifail=8
On entry, ldz < max1,n  or ldz < 1  when no vectors are required.
ifail=9
On entry, the option Vectors = Select  was selected, but this is not yet implemented.
ifail=10
The number of eigenvalues found to sufficient accuracy is zero.
ifail=11
Could not build a Lanczos factorization. Consider changing ncv or nev in the initialization routine (see Section 5 in f12faf for details of these arguments).
ifail=12
Unexpected error in internal call to compute eigenvalues and corresponding error bounds of the current symmetric tridiagonal matrix. Please contact NAG.
ifail=13
Unexpected error during calculation of a real Schur form: there was a failure to compute all the converged eigenvalues. Please contact NAG.
ifail=14
Failure during internal factorization of real banded matrix. Please contact NAG.
ifail=15
Failure during internal solution of real banded system. Please contact NAG.
ifail=16
The maximum number of iterations has been reached. Some Ritz values may have converged; nconv returns the number of converged values.
ifail=17
No shifts could be applied during a cycle of the implicitly restarted Lanczos iteration. One possibility is to increase the size of ncv relative to nev (see Section 5 in f12fff for details of these arguments).
ifail=18
Either an initial call to the setup routine has not been made or the communication arrays have become corrupted.
ifail=19
The routine was unable to dynamically allocate sufficient internal workspace. Please contact NAG.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 3.9 in How to Use the NAG Library and its Documentation for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 3.8 in How to Use the NAG Library and its Documentation for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 3.7 in How to Use the NAG Library and its Documentation for further information.

7
Accuracy

The relative accuracy of a Ritz value, λ , is considered acceptable if its Ritz estimate Tolerance × λ . The default Tolerance used is the machine precision given by x02ajf.

8
Parallelism and Performance

f12fgf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f12fgf 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

None.

10
Example

This example solves Ax = λx  in regular mode, where A  is obtained from the standard central difference discretization of the two-dimensional convection-diffusion operator d2u dx2 + d2u dy2 = ρ du dx  on the unit square with zero Dirichlet boundary conditions. A  is stored in LAPACK banded storage format.

10.1
Program Text

Program Text (f12fgfe.f90)

10.2
Program Data

Program Data (f12fgfe.d)

10.3
Program Results

Program Results (f12fgfe.r)

© The Numerical Algorithms Group Ltd, Oxford, UK. 2017