00001 SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB, 00002 $ LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, 00003 $ LDZ, WORK, IWORK, IFAIL, INFO ) 00004 * 00005 * -- LAPACK driver routine (version 3.2) -- 00006 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00007 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00008 * November 2006 00009 * 00010 * .. Scalar Arguments .. 00011 CHARACTER JOBZ, RANGE, UPLO 00012 INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M, 00013 $ N 00014 DOUBLE PRECISION ABSTOL, VL, VU 00015 * .. 00016 * .. Array Arguments .. 00017 INTEGER IFAIL( * ), IWORK( * ) 00018 DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ), 00019 $ W( * ), WORK( * ), Z( LDZ, * ) 00020 * .. 00021 * 00022 * Purpose 00023 * ======= 00024 * 00025 * DSBGVX computes selected eigenvalues, and optionally, eigenvectors 00026 * of a real generalized symmetric-definite banded eigenproblem, of 00027 * the form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric 00028 * and banded, and B is also positive definite. Eigenvalues and 00029 * eigenvectors can be selected by specifying either all eigenvalues, 00030 * a range of values or a range of indices for the desired eigenvalues. 00031 * 00032 * Arguments 00033 * ========= 00034 * 00035 * JOBZ (input) CHARACTER*1 00036 * = 'N': Compute eigenvalues only; 00037 * = 'V': Compute eigenvalues and eigenvectors. 00038 * 00039 * RANGE (input) CHARACTER*1 00040 * = 'A': all eigenvalues will be found. 00041 * = 'V': all eigenvalues in the half-open interval (VL,VU] 00042 * will be found. 00043 * = 'I': the IL-th through IU-th eigenvalues will be found. 00044 * 00045 * UPLO (input) CHARACTER*1 00046 * = 'U': Upper triangles of A and B are stored; 00047 * = 'L': Lower triangles of A and B are stored. 00048 * 00049 * N (input) INTEGER 00050 * The order of the matrices A and B. N >= 0. 00051 * 00052 * KA (input) INTEGER 00053 * The number of superdiagonals of the matrix A if UPLO = 'U', 00054 * or the number of subdiagonals if UPLO = 'L'. KA >= 0. 00055 * 00056 * KB (input) INTEGER 00057 * The number of superdiagonals of the matrix B if UPLO = 'U', 00058 * or the number of subdiagonals if UPLO = 'L'. KB >= 0. 00059 * 00060 * AB (input/output) DOUBLE PRECISION array, dimension (LDAB, N) 00061 * On entry, the upper or lower triangle of the symmetric band 00062 * matrix A, stored in the first ka+1 rows of the array. The 00063 * j-th column of A is stored in the j-th column of the array AB 00064 * as follows: 00065 * if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j; 00066 * if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka). 00067 * 00068 * On exit, the contents of AB are destroyed. 00069 * 00070 * LDAB (input) INTEGER 00071 * The leading dimension of the array AB. LDAB >= KA+1. 00072 * 00073 * BB (input/output) DOUBLE PRECISION array, dimension (LDBB, N) 00074 * On entry, the upper or lower triangle of the symmetric band 00075 * matrix B, stored in the first kb+1 rows of the array. The 00076 * j-th column of B is stored in the j-th column of the array BB 00077 * as follows: 00078 * if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j; 00079 * if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb). 00080 * 00081 * On exit, the factor S from the split Cholesky factorization 00082 * B = S**T*S, as returned by DPBSTF. 00083 * 00084 * LDBB (input) INTEGER 00085 * The leading dimension of the array BB. LDBB >= KB+1. 00086 * 00087 * Q (output) DOUBLE PRECISION array, dimension (LDQ, N) 00088 * If JOBZ = 'V', the n-by-n matrix used in the reduction of 00089 * A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x, 00090 * and consequently C to tridiagonal form. 00091 * If JOBZ = 'N', the array Q is not referenced. 00092 * 00093 * LDQ (input) INTEGER 00094 * The leading dimension of the array Q. If JOBZ = 'N', 00095 * LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N). 00096 * 00097 * VL (input) DOUBLE PRECISION 00098 * VU (input) DOUBLE PRECISION 00099 * If RANGE='V', the lower and upper bounds of the interval to 00100 * be searched for eigenvalues. VL < VU. 00101 * Not referenced if RANGE = 'A' or 'I'. 00102 * 00103 * IL (input) INTEGER 00104 * IU (input) INTEGER 00105 * If RANGE='I', the indices (in ascending order) of the 00106 * smallest and largest eigenvalues to be returned. 00107 * 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. 00108 * Not referenced if RANGE = 'A' or 'V'. 00109 * 00110 * ABSTOL (input) DOUBLE PRECISION 00111 * The absolute error tolerance for the eigenvalues. 00112 * An approximate eigenvalue is accepted as converged 00113 * when it is determined to lie in an interval [a,b] 00114 * of width less than or equal to 00115 * 00116 * ABSTOL + EPS * max( |a|,|b| ) , 00117 * 00118 * where EPS is the machine precision. If ABSTOL is less than 00119 * or equal to zero, then EPS*|T| will be used in its place, 00120 * where |T| is the 1-norm of the tridiagonal matrix obtained 00121 * by reducing A to tridiagonal form. 00122 * 00123 * Eigenvalues will be computed most accurately when ABSTOL is 00124 * set to twice the underflow threshold 2*DLAMCH('S'), not zero. 00125 * If this routine returns with INFO>0, indicating that some 00126 * eigenvectors did not converge, try setting ABSTOL to 00127 * 2*DLAMCH('S'). 00128 * 00129 * M (output) INTEGER 00130 * The total number of eigenvalues found. 0 <= M <= N. 00131 * If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1. 00132 * 00133 * W (output) DOUBLE PRECISION array, dimension (N) 00134 * If INFO = 0, the eigenvalues in ascending order. 00135 * 00136 * Z (output) DOUBLE PRECISION array, dimension (LDZ, N) 00137 * If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of 00138 * eigenvectors, with the i-th column of Z holding the 00139 * eigenvector associated with W(i). The eigenvectors are 00140 * normalized so Z**T*B*Z = I. 00141 * If JOBZ = 'N', then Z is not referenced. 00142 * 00143 * LDZ (input) INTEGER 00144 * The leading dimension of the array Z. LDZ >= 1, and if 00145 * JOBZ = 'V', LDZ >= max(1,N). 00146 * 00147 * WORK (workspace/output) DOUBLE PRECISION array, dimension (7*N) 00148 * 00149 * IWORK (workspace/output) INTEGER array, dimension (5*N) 00150 * 00151 * IFAIL (output) INTEGER array, dimension (M) 00152 * If JOBZ = 'V', then if INFO = 0, the first M elements of 00153 * IFAIL are zero. If INFO > 0, then IFAIL contains the 00154 * indices of the eigenvalues that failed to converge. 00155 * If JOBZ = 'N', then IFAIL is not referenced. 00156 * 00157 * INFO (output) INTEGER 00158 * = 0 : successful exit 00159 * < 0 : if INFO = -i, the i-th argument had an illegal value 00160 * <= N: if INFO = i, then i eigenvectors failed to converge. 00161 * Their indices are stored in IFAIL. 00162 * > N : DPBSTF returned an error code; i.e., 00163 * if INFO = N + i, for 1 <= i <= N, then the leading 00164 * minor of order i of B is not positive definite. 00165 * The factorization of B could not be completed and 00166 * no eigenvalues or eigenvectors were computed. 00167 * 00168 * Further Details 00169 * =============== 00170 * 00171 * Based on contributions by 00172 * Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA 00173 * 00174 * ===================================================================== 00175 * 00176 * .. Parameters .. 00177 DOUBLE PRECISION ZERO, ONE 00178 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 00179 * .. 00180 * .. Local Scalars .. 00181 LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ 00182 CHARACTER ORDER, VECT 00183 INTEGER I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP, 00184 $ INDIWO, INDWRK, ITMP1, J, JJ, NSPLIT 00185 DOUBLE PRECISION TMP1 00186 * .. 00187 * .. External Functions .. 00188 LOGICAL LSAME 00189 EXTERNAL LSAME 00190 * .. 00191 * .. External Subroutines .. 00192 EXTERNAL DCOPY, DGEMV, DLACPY, DPBSTF, DSBGST, DSBTRD, 00193 $ DSTEBZ, DSTEIN, DSTEQR, DSTERF, DSWAP, XERBLA 00194 * .. 00195 * .. Intrinsic Functions .. 00196 INTRINSIC MIN 00197 * .. 00198 * .. Executable Statements .. 00199 * 00200 * Test the input parameters. 00201 * 00202 WANTZ = LSAME( JOBZ, 'V' ) 00203 UPPER = LSAME( UPLO, 'U' ) 00204 ALLEIG = LSAME( RANGE, 'A' ) 00205 VALEIG = LSAME( RANGE, 'V' ) 00206 INDEIG = LSAME( RANGE, 'I' ) 00207 * 00208 INFO = 0 00209 IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN 00210 INFO = -1 00211 ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN 00212 INFO = -2 00213 ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN 00214 INFO = -3 00215 ELSE IF( N.LT.0 ) THEN 00216 INFO = -4 00217 ELSE IF( KA.LT.0 ) THEN 00218 INFO = -5 00219 ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN 00220 INFO = -6 00221 ELSE IF( LDAB.LT.KA+1 ) THEN 00222 INFO = -8 00223 ELSE IF( LDBB.LT.KB+1 ) THEN 00224 INFO = -10 00225 ELSE IF( LDQ.LT.1 .OR. ( WANTZ .AND. LDQ.LT.N ) ) THEN 00226 INFO = -12 00227 ELSE 00228 IF( VALEIG ) THEN 00229 IF( N.GT.0 .AND. VU.LE.VL ) 00230 $ INFO = -14 00231 ELSE IF( INDEIG ) THEN 00232 IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN 00233 INFO = -15 00234 ELSE IF ( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN 00235 INFO = -16 00236 END IF 00237 END IF 00238 END IF 00239 IF( INFO.EQ.0) THEN 00240 IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN 00241 INFO = -21 00242 END IF 00243 END IF 00244 * 00245 IF( INFO.NE.0 ) THEN 00246 CALL XERBLA( 'DSBGVX', -INFO ) 00247 RETURN 00248 END IF 00249 * 00250 * Quick return if possible 00251 * 00252 M = 0 00253 IF( N.EQ.0 ) 00254 $ RETURN 00255 * 00256 * Form a split Cholesky factorization of B. 00257 * 00258 CALL DPBSTF( UPLO, N, KB, BB, LDBB, INFO ) 00259 IF( INFO.NE.0 ) THEN 00260 INFO = N + INFO 00261 RETURN 00262 END IF 00263 * 00264 * Transform problem to standard eigenvalue problem. 00265 * 00266 CALL DSBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ, 00267 $ WORK, IINFO ) 00268 * 00269 * Reduce symmetric band matrix to tridiagonal form. 00270 * 00271 INDD = 1 00272 INDE = INDD + N 00273 INDWRK = INDE + N 00274 IF( WANTZ ) THEN 00275 VECT = 'U' 00276 ELSE 00277 VECT = 'N' 00278 END IF 00279 CALL DSBTRD( VECT, UPLO, N, KA, AB, LDAB, WORK( INDD ), 00280 $ WORK( INDE ), Q, LDQ, WORK( INDWRK ), IINFO ) 00281 * 00282 * If all eigenvalues are desired and ABSTOL is less than or equal 00283 * to zero, then call DSTERF or SSTEQR. If this fails for some 00284 * eigenvalue, then try DSTEBZ. 00285 * 00286 TEST = .FALSE. 00287 IF( INDEIG ) THEN 00288 IF( IL.EQ.1 .AND. IU.EQ.N ) THEN 00289 TEST = .TRUE. 00290 END IF 00291 END IF 00292 IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN 00293 CALL DCOPY( N, WORK( INDD ), 1, W, 1 ) 00294 INDEE = INDWRK + 2*N 00295 CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 ) 00296 IF( .NOT.WANTZ ) THEN 00297 CALL DSTERF( N, W, WORK( INDEE ), INFO ) 00298 ELSE 00299 CALL DLACPY( 'A', N, N, Q, LDQ, Z, LDZ ) 00300 CALL DSTEQR( JOBZ, N, W, WORK( INDEE ), Z, LDZ, 00301 $ WORK( INDWRK ), INFO ) 00302 IF( INFO.EQ.0 ) THEN 00303 DO 10 I = 1, N 00304 IFAIL( I ) = 0 00305 10 CONTINUE 00306 END IF 00307 END IF 00308 IF( INFO.EQ.0 ) THEN 00309 M = N 00310 GO TO 30 00311 END IF 00312 INFO = 0 00313 END IF 00314 * 00315 * Otherwise, call DSTEBZ and, if eigenvectors are desired, 00316 * call DSTEIN. 00317 * 00318 IF( WANTZ ) THEN 00319 ORDER = 'B' 00320 ELSE 00321 ORDER = 'E' 00322 END IF 00323 INDIBL = 1 00324 INDISP = INDIBL + N 00325 INDIWO = INDISP + N 00326 CALL DSTEBZ( RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, 00327 $ WORK( INDD ), WORK( INDE ), M, NSPLIT, W, 00328 $ IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWRK ), 00329 $ IWORK( INDIWO ), INFO ) 00330 * 00331 IF( WANTZ ) THEN 00332 CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W, 00333 $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ, 00334 $ WORK( INDWRK ), IWORK( INDIWO ), IFAIL, INFO ) 00335 * 00336 * Apply transformation matrix used in reduction to tridiagonal 00337 * form to eigenvectors returned by DSTEIN. 00338 * 00339 DO 20 J = 1, M 00340 CALL DCOPY( N, Z( 1, J ), 1, WORK( 1 ), 1 ) 00341 CALL DGEMV( 'N', N, N, ONE, Q, LDQ, WORK, 1, ZERO, 00342 $ Z( 1, J ), 1 ) 00343 20 CONTINUE 00344 END IF 00345 * 00346 30 CONTINUE 00347 * 00348 * If eigenvalues are not in order, then sort them, along with 00349 * eigenvectors. 00350 * 00351 IF( WANTZ ) THEN 00352 DO 50 J = 1, M - 1 00353 I = 0 00354 TMP1 = W( J ) 00355 DO 40 JJ = J + 1, M 00356 IF( W( JJ ).LT.TMP1 ) THEN 00357 I = JJ 00358 TMP1 = W( JJ ) 00359 END IF 00360 40 CONTINUE 00361 * 00362 IF( I.NE.0 ) THEN 00363 ITMP1 = IWORK( INDIBL+I-1 ) 00364 W( I ) = W( J ) 00365 IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 ) 00366 W( J ) = TMP1 00367 IWORK( INDIBL+J-1 ) = ITMP1 00368 CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 ) 00369 IF( INFO.NE.0 ) THEN 00370 ITMP1 = IFAIL( I ) 00371 IFAIL( I ) = IFAIL( J ) 00372 IFAIL( J ) = ITMP1 00373 END IF 00374 END IF 00375 50 CONTINUE 00376 END IF 00377 * 00378 RETURN 00379 * 00380 * End of DSBGVX 00381 * 00382 END