NAG FL Interface
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.
Settings help

FL Name Style:


FL Specification Language:


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 <nag.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)
The routine may be called by the names f12fgf or nagf_sparseig_real_symm_band_solve.

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):

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, Philadelphia

5 Arguments

1: kl Integer Input
On entry: the number of subdiagonals of the matrices A and B.
Constraint: kl0.
2: ku Integer Input
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: ab(ldab,*) Real (Kind=nag_wp) array Input
Note: the second dimension of the array ab must be at least max(1,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 Integer Input
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: mb(ldmb,*) Real (Kind=nag_wp) array Input
Note: the second dimension of the array mb must be at least max(1,n) if a Generalized problem has be selected (see f12fdf); otherwise mb is not referenced. (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 Integer Input
On entry: the first dimension of the array mb as declared in the (sub)program from which f12fgf is called.
Constraint: if the optional parameter Generalized has been selected, 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 Integer Output
On exit: the number of converged eigenvalues.
9: d(*) Real (Kind=nag_wp) array Output
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: z(ldz,*) Real (Kind=nag_wp) array Output
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 Integer Input
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) array Input/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: v(ldv,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array v must be at least max(1,ncv) (see f12fff).
On exit: if the option Vectors=SCHUR or RITZ (see f12fdf) 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 Integer Input
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) array Communication Array
Note: the actual argument supplied must be the array comm supplied to the initialization routine f12fdf.
On initial entry: must remain unchanged from the prior call to f12fdf and f12fff.
On exit: contains no useful information.
16: icomm(*) Integer array Communication Array
Note: the actual argument supplied must be the array icomm supplied to the initialization routine f12fff.
On initial entry: must remain unchanged from the prior call to f12fbf and f12fdf.
On exit: contains no useful information.
17: 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
On entry, kl=value.
Constraint: kl 0.
ifail=2
On entry, ku=value.
Constraint: ku 0.
ifail=3
On entry, ldab=value, 2×kl+ku+1=value.
Constraint: ldab2×kl+ku+1.
ifail=4
The maximum number of iterations 0, the option Iteration Limit has been set to value.
ifail=5
The options Generalized and Regular are incompatible.
ifail=6
Eigenvalues from both ends of the spectrum were requested, but the number of eigenvalues (nev in f12fff) requested is one.
ifail=7
The option Initial Residual was selected but the starting vector held in resid is zero.
ifail=8
On entry, ldz=value, n=value (see n in f12fff).
Constraint: ldzmax(1,n).
ifail=9
On entry, Vectors = Select , 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. The size of the current Lanczos factorization = value.
ifail=12
Error in internal call to compute eigenvalues and corresponding error bounds of the current upper Hessenberg matrix. Please contact NAG.
ifail=13
During calculation of a tridiagonal form, there was a failure to compute value eigenvalues in a total of value iterations.
ifail=14
Failure during internal factorization of banded matrix. Please contact NAG.
ifail=15
Failure during internal solution of banded system. Please contact NAG.
ifail=16
The maximum number of iterations has been reached: there have been value iterations. There are value converged eigenvalues.
ifail=17
No shifts could be applied during a cycle of the implicitly restarted Lanczos iteration.
ifail=18
Either an initial call to the setup routine has not been made or the communication arrays have become corrupted.
ifail=19
On entry, ldmb=value, 2×kl+ku+1=value.
Constraint: ldmb2×kl+ku+1.
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

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)