NAG CL Interface
f12fgc (real_​symm_​band_​solve)

Note: this function 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 function f12fdc need not be called. If, however, you wish to reset some or all of the settings please refer to Section 11 in f12fdc for a detailed description of the specification of the optional parameters.
Settings help

CL Name Style:


1 Purpose

f12fgc is the main solver function in a suite of functions which includes f12fdc and f12ffc. f12fgc must be called following an initial call to f12ffc and following any calls to f12fdc.
f12fgc 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

#include <nag.h>
void  f12fgc (Integer kl, Integer ku, const double ab[], const double mb[], double sigma, Integer *nconv, double d[], double z[], double resid[], double v[], double comm[], Integer icomm[], NagError *fail)
The function may be called by the names: f12fgc, nag_sparseig_real_symm_band_solve or nag_real_symm_banded_sparse_eigensystem_sol.

3 Description

The suite of functions 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 function f12ffc, f12fgc 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.4.2 in the F07 Chapter Introduction for details on this storage format.
f12fgc is based on the banded driver functions 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 functions 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 function.
f12fgc, is a general purpose direct communication function that must be called following initialization by f12ffc. f12fgc uses options, set either by default or explicitly by calling f12fdc, 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[dim] const double Input
Note: the dimension, dim, of the array ab must be at least max(1, n×(2×kl+ku+1) )  (see f12ffc).
On entry: must contain the matrix A in LAPACK column-ordered banded storage format for nonsymmetric matrices (see Section 3.4.4 in the F07 Chapter Introduction).
4: mb[dim] const double Input
Note: the dimension, dim, of the array mb must be at least max(1,n×(2×kl+ku+1)) (see f12ffc).
On entry: must contain the matrix B in LAPACK column-ordered banded storage format for nonsymmetric matrices (see Section 3.4.4 in the F07 Chapter Introduction).
5: sigma double Input
On entry: if one of the Shifted Inverse (see f12fdc) modes has been selected then sigma contains the real shift used; otherwise sigma is not referenced.
6: nconv Integer * Output
On exit: the number of converged eigenvalues.
7: d[dim] double Output
Note: the dimension, dim, of the array d must be at least ncv (see f12ffc).
On exit: the first nconv locations of the array d contain the converged approximate eigenvalues.
8: z[n×(nev+1)] double Output
On exit: if the default option Vectors=RITZ (see f12fdc) 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 at locations z[i-1×n+j-1], for j=1,2,,n.
9: resid[dim] double Input/Output
Note: the dimension, dim, of the array resid must be at least n (see f12ffc).
On entry: need not be set unless the option Initial Residual has been set in a prior call to f12fdc in which case resid must contain an initial residual vector.
On exit: contains the final residual vector.
10: v[n×ncv] double Output
On exit: if the option Vectors=SCHUR or RITZ (see f12fdc) and z does not equal v then the first nconv sections of v, of length n, will contain approximate Schur vectors that span the desired invariant subspace.
The jth Schur vector is stored in locations v[n×(j-1)+i-1], for j=1,2,,nconv and i=1,2,,n.
11: comm[dim] double Communication Array
Note: the actual argument supplied must be the array comm supplied to the initialization routine f12fdc.
On initial entry: must remain unchanged from the prior call to f12fdc and f12ffc.
On exit: contains no useful information.
12: icomm[dim] Integer Communication Array
Note: the actual argument supplied must be the array icomm supplied to the initialization routine f12ffc.
On initial entry: must remain unchanged from the prior call to f12fbc and f12fdc.
On exit: contains no useful information.
13: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.1.2 in the Introduction to the NAG Library CL Interface for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_BOTH_ENDS_1
Eigenvalues from both ends of the spectrum were requested, but the number of eigenvalues (nev in f12ffc) requested is one.
NE_INITIALIZATION
Either an initial call to the setup function has not been made or the communication arrays have become corrupted.
NE_INT
On entry, kl=value.
Constraint: kl0.
On entry, ku=value.
Constraint: ku0.
The maximum number of iterations 0, the option Iteration Limit has been set to value.
NE_INT_2
The maximum number of iterations has been reached: there have been value iterations. There are value converged eigenvalues.
NE_INTERNAL_EIGVAL_FAIL
Error in internal call to compute eigenvalues and corresponding error bounds of the current upper Hessenberg matrix. Please contact NAG.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
See Section 7.5 in the Introduction to the NAG Library CL Interface for further information.
NE_INVALID_OPTION
On entry, Vectors = Select , but this is not yet implemented.
NE_MAX_ITER
During calculation of a tridiagonal form, there was a failure to compute value eigenvalues in a total of value iterations.
NE_NO_LANCZOS_FAC
Could not build a Lanczos factorization. The size of the current Lanczos factorization = value.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library CL Interface for further information.
NE_NO_SHIFTS_APPLIED
No shifts could be applied during a cycle of the implicitly restarted Lanczos iteration.
NE_OPT_INCOMPAT
The options Generalized and Regular are incompatible.
NE_REAL_BAND_FAC
Failure during internal factorization of banded matrix. Please contact NAG.
NE_REAL_BAND_SOL
Failure during internal solution of banded system. Please contact NAG.
NE_ZERO_EIGS_FOUND
The number of eigenvalues found to sufficient accuracy is zero.
NE_ZERO_INIT_RESID
The option Initial Residual was selected but the starting vector held in resid is zero.

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 X02AJC.

8 Parallelism and Performance

Background information to multithreading can be found in the Multithreading documentation.
f12fgc is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f12fgc 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 function. 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 (f12fgce.c)

10.2 Program Data

Program Data (f12fgce.d)

10.3 Program Results

Program Results (f12fgce.r)