Annotation of rpl/lapack/lapack/dsbgvx.f, revision 1.8
1.8 ! bertrand 1: *> \brief \b DSBGST
! 2: *
! 3: * =========== DOCUMENTATION ===========
! 4: *
! 5: * Online html documentation available at
! 6: * http://www.netlib.org/lapack/explore-html/
! 7: *
! 8: *> \htmlonly
! 9: *> Download DSBGVX + dependencies
! 10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbgvx.f">
! 11: *> [TGZ]</a>
! 12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbgvx.f">
! 13: *> [ZIP]</a>
! 14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbgvx.f">
! 15: *> [TXT]</a>
! 16: *> \endhtmlonly
! 17: *
! 18: * Definition:
! 19: * ===========
! 20: *
! 21: * SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
! 22: * LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
! 23: * LDZ, WORK, IWORK, IFAIL, INFO )
! 24: *
! 25: * .. Scalar Arguments ..
! 26: * CHARACTER JOBZ, RANGE, UPLO
! 27: * INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
! 28: * $ N
! 29: * DOUBLE PRECISION ABSTOL, VL, VU
! 30: * ..
! 31: * .. Array Arguments ..
! 32: * INTEGER IFAIL( * ), IWORK( * )
! 33: * DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
! 34: * $ W( * ), WORK( * ), Z( LDZ, * )
! 35: * ..
! 36: *
! 37: *
! 38: *> \par Purpose:
! 39: * =============
! 40: *>
! 41: *> \verbatim
! 42: *>
! 43: *> DSBGVX computes selected eigenvalues, and optionally, eigenvectors
! 44: *> of a real generalized symmetric-definite banded eigenproblem, of
! 45: *> the form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric
! 46: *> and banded, and B is also positive definite. Eigenvalues and
! 47: *> eigenvectors can be selected by specifying either all eigenvalues,
! 48: *> a range of values or a range of indices for the desired eigenvalues.
! 49: *> \endverbatim
! 50: *
! 51: * Arguments:
! 52: * ==========
! 53: *
! 54: *> \param[in] JOBZ
! 55: *> \verbatim
! 56: *> JOBZ is CHARACTER*1
! 57: *> = 'N': Compute eigenvalues only;
! 58: *> = 'V': Compute eigenvalues and eigenvectors.
! 59: *> \endverbatim
! 60: *>
! 61: *> \param[in] RANGE
! 62: *> \verbatim
! 63: *> RANGE is CHARACTER*1
! 64: *> = 'A': all eigenvalues will be found.
! 65: *> = 'V': all eigenvalues in the half-open interval (VL,VU]
! 66: *> will be found.
! 67: *> = 'I': the IL-th through IU-th eigenvalues will be found.
! 68: *> \endverbatim
! 69: *>
! 70: *> \param[in] UPLO
! 71: *> \verbatim
! 72: *> UPLO is CHARACTER*1
! 73: *> = 'U': Upper triangles of A and B are stored;
! 74: *> = 'L': Lower triangles of A and B are stored.
! 75: *> \endverbatim
! 76: *>
! 77: *> \param[in] N
! 78: *> \verbatim
! 79: *> N is INTEGER
! 80: *> The order of the matrices A and B. N >= 0.
! 81: *> \endverbatim
! 82: *>
! 83: *> \param[in] KA
! 84: *> \verbatim
! 85: *> KA is INTEGER
! 86: *> The number of superdiagonals of the matrix A if UPLO = 'U',
! 87: *> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
! 88: *> \endverbatim
! 89: *>
! 90: *> \param[in] KB
! 91: *> \verbatim
! 92: *> KB is INTEGER
! 93: *> The number of superdiagonals of the matrix B if UPLO = 'U',
! 94: *> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
! 95: *> \endverbatim
! 96: *>
! 97: *> \param[in,out] AB
! 98: *> \verbatim
! 99: *> AB is DOUBLE PRECISION array, dimension (LDAB, N)
! 100: *> On entry, the upper or lower triangle of the symmetric band
! 101: *> matrix A, stored in the first ka+1 rows of the array. The
! 102: *> j-th column of A is stored in the j-th column of the array AB
! 103: *> as follows:
! 104: *> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
! 105: *> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
! 106: *>
! 107: *> On exit, the contents of AB are destroyed.
! 108: *> \endverbatim
! 109: *>
! 110: *> \param[in] LDAB
! 111: *> \verbatim
! 112: *> LDAB is INTEGER
! 113: *> The leading dimension of the array AB. LDAB >= KA+1.
! 114: *> \endverbatim
! 115: *>
! 116: *> \param[in,out] BB
! 117: *> \verbatim
! 118: *> BB is DOUBLE PRECISION array, dimension (LDBB, N)
! 119: *> On entry, the upper or lower triangle of the symmetric band
! 120: *> matrix B, stored in the first kb+1 rows of the array. The
! 121: *> j-th column of B is stored in the j-th column of the array BB
! 122: *> as follows:
! 123: *> if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
! 124: *> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
! 125: *>
! 126: *> On exit, the factor S from the split Cholesky factorization
! 127: *> B = S**T*S, as returned by DPBSTF.
! 128: *> \endverbatim
! 129: *>
! 130: *> \param[in] LDBB
! 131: *> \verbatim
! 132: *> LDBB is INTEGER
! 133: *> The leading dimension of the array BB. LDBB >= KB+1.
! 134: *> \endverbatim
! 135: *>
! 136: *> \param[out] Q
! 137: *> \verbatim
! 138: *> Q is DOUBLE PRECISION array, dimension (LDQ, N)
! 139: *> If JOBZ = 'V', the n-by-n matrix used in the reduction of
! 140: *> A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
! 141: *> and consequently C to tridiagonal form.
! 142: *> If JOBZ = 'N', the array Q is not referenced.
! 143: *> \endverbatim
! 144: *>
! 145: *> \param[in] LDQ
! 146: *> \verbatim
! 147: *> LDQ is INTEGER
! 148: *> The leading dimension of the array Q. If JOBZ = 'N',
! 149: *> LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).
! 150: *> \endverbatim
! 151: *>
! 152: *> \param[in] VL
! 153: *> \verbatim
! 154: *> VL is DOUBLE PRECISION
! 155: *> \endverbatim
! 156: *>
! 157: *> \param[in] VU
! 158: *> \verbatim
! 159: *> VU is DOUBLE PRECISION
! 160: *>
! 161: *> If RANGE='V', the lower and upper bounds of the interval to
! 162: *> be searched for eigenvalues. VL < VU.
! 163: *> Not referenced if RANGE = 'A' or 'I'.
! 164: *> \endverbatim
! 165: *>
! 166: *> \param[in] IL
! 167: *> \verbatim
! 168: *> IL is INTEGER
! 169: *> \endverbatim
! 170: *>
! 171: *> \param[in] IU
! 172: *> \verbatim
! 173: *> IU is INTEGER
! 174: *>
! 175: *> If RANGE='I', the indices (in ascending order) of the
! 176: *> smallest and largest eigenvalues to be returned.
! 177: *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
! 178: *> Not referenced if RANGE = 'A' or 'V'.
! 179: *> \endverbatim
! 180: *>
! 181: *> \param[in] ABSTOL
! 182: *> \verbatim
! 183: *> ABSTOL is DOUBLE PRECISION
! 184: *> The absolute error tolerance for the eigenvalues.
! 185: *> An approximate eigenvalue is accepted as converged
! 186: *> when it is determined to lie in an interval [a,b]
! 187: *> of width less than or equal to
! 188: *>
! 189: *> ABSTOL + EPS * max( |a|,|b| ) ,
! 190: *>
! 191: *> where EPS is the machine precision. If ABSTOL is less than
! 192: *> or equal to zero, then EPS*|T| will be used in its place,
! 193: *> where |T| is the 1-norm of the tridiagonal matrix obtained
! 194: *> by reducing A to tridiagonal form.
! 195: *>
! 196: *> Eigenvalues will be computed most accurately when ABSTOL is
! 197: *> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
! 198: *> If this routine returns with INFO>0, indicating that some
! 199: *> eigenvectors did not converge, try setting ABSTOL to
! 200: *> 2*DLAMCH('S').
! 201: *> \endverbatim
! 202: *>
! 203: *> \param[out] M
! 204: *> \verbatim
! 205: *> M is INTEGER
! 206: *> The total number of eigenvalues found. 0 <= M <= N.
! 207: *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
! 208: *> \endverbatim
! 209: *>
! 210: *> \param[out] W
! 211: *> \verbatim
! 212: *> W is DOUBLE PRECISION array, dimension (N)
! 213: *> If INFO = 0, the eigenvalues in ascending order.
! 214: *> \endverbatim
! 215: *>
! 216: *> \param[out] Z
! 217: *> \verbatim
! 218: *> Z is DOUBLE PRECISION array, dimension (LDZ, N)
! 219: *> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
! 220: *> eigenvectors, with the i-th column of Z holding the
! 221: *> eigenvector associated with W(i). The eigenvectors are
! 222: *> normalized so Z**T*B*Z = I.
! 223: *> If JOBZ = 'N', then Z is not referenced.
! 224: *> \endverbatim
! 225: *>
! 226: *> \param[in] LDZ
! 227: *> \verbatim
! 228: *> LDZ is INTEGER
! 229: *> The leading dimension of the array Z. LDZ >= 1, and if
! 230: *> JOBZ = 'V', LDZ >= max(1,N).
! 231: *> \endverbatim
! 232: *>
! 233: *> \param[out] WORK
! 234: *> \verbatim
! 235: *> WORK is DOUBLE PRECISION array, dimension (7*N)
! 236: *> \endverbatim
! 237: *>
! 238: *> \param[out] IWORK
! 239: *> \verbatim
! 240: *> IWORK is INTEGER array, dimension (5*N)
! 241: *> \endverbatim
! 242: *>
! 243: *> \param[out] IFAIL
! 244: *> \verbatim
! 245: *> IFAIL is INTEGER array, dimension (M)
! 246: *> If JOBZ = 'V', then if INFO = 0, the first M elements of
! 247: *> IFAIL are zero. If INFO > 0, then IFAIL contains the
! 248: *> indices of the eigenvalues that failed to converge.
! 249: *> If JOBZ = 'N', then IFAIL is not referenced.
! 250: *> \endverbatim
! 251: *>
! 252: *> \param[out] INFO
! 253: *> \verbatim
! 254: *> INFO is INTEGER
! 255: *> = 0 : successful exit
! 256: *> < 0 : if INFO = -i, the i-th argument had an illegal value
! 257: *> <= N: if INFO = i, then i eigenvectors failed to converge.
! 258: *> Their indices are stored in IFAIL.
! 259: *> > N : DPBSTF returned an error code; i.e.,
! 260: *> if INFO = N + i, for 1 <= i <= N, then the leading
! 261: *> minor of order i of B is not positive definite.
! 262: *> The factorization of B could not be completed and
! 263: *> no eigenvalues or eigenvectors were computed.
! 264: *> \endverbatim
! 265: *
! 266: * Authors:
! 267: * ========
! 268: *
! 269: *> \author Univ. of Tennessee
! 270: *> \author Univ. of California Berkeley
! 271: *> \author Univ. of Colorado Denver
! 272: *> \author NAG Ltd.
! 273: *
! 274: *> \date November 2011
! 275: *
! 276: *> \ingroup doubleOTHEReigen
! 277: *
! 278: *> \par Contributors:
! 279: * ==================
! 280: *>
! 281: *> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
! 282: *
! 283: * =====================================================================
1.1 bertrand 284: SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
285: $ LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
286: $ LDZ, WORK, IWORK, IFAIL, INFO )
287: *
1.8 ! bertrand 288: * -- LAPACK driver routine (version 3.4.0) --
1.1 bertrand 289: * -- LAPACK is a software package provided by Univ. of Tennessee, --
290: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.8 ! bertrand 291: * November 2011
1.1 bertrand 292: *
293: * .. Scalar Arguments ..
294: CHARACTER JOBZ, RANGE, UPLO
295: INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
296: $ N
297: DOUBLE PRECISION ABSTOL, VL, VU
298: * ..
299: * .. Array Arguments ..
300: INTEGER IFAIL( * ), IWORK( * )
301: DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
302: $ W( * ), WORK( * ), Z( LDZ, * )
303: * ..
304: *
305: * =====================================================================
306: *
307: * .. Parameters ..
308: DOUBLE PRECISION ZERO, ONE
309: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
310: * ..
311: * .. Local Scalars ..
312: LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
313: CHARACTER ORDER, VECT
314: INTEGER I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP,
315: $ INDIWO, INDWRK, ITMP1, J, JJ, NSPLIT
316: DOUBLE PRECISION TMP1
317: * ..
318: * .. External Functions ..
319: LOGICAL LSAME
320: EXTERNAL LSAME
321: * ..
322: * .. External Subroutines ..
323: EXTERNAL DCOPY, DGEMV, DLACPY, DPBSTF, DSBGST, DSBTRD,
324: $ DSTEBZ, DSTEIN, DSTEQR, DSTERF, DSWAP, XERBLA
325: * ..
326: * .. Intrinsic Functions ..
327: INTRINSIC MIN
328: * ..
329: * .. Executable Statements ..
330: *
331: * Test the input parameters.
332: *
333: WANTZ = LSAME( JOBZ, 'V' )
334: UPPER = LSAME( UPLO, 'U' )
335: ALLEIG = LSAME( RANGE, 'A' )
336: VALEIG = LSAME( RANGE, 'V' )
337: INDEIG = LSAME( RANGE, 'I' )
338: *
339: INFO = 0
340: IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
341: INFO = -1
342: ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
343: INFO = -2
344: ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
345: INFO = -3
346: ELSE IF( N.LT.0 ) THEN
347: INFO = -4
348: ELSE IF( KA.LT.0 ) THEN
349: INFO = -5
350: ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN
351: INFO = -6
352: ELSE IF( LDAB.LT.KA+1 ) THEN
353: INFO = -8
354: ELSE IF( LDBB.LT.KB+1 ) THEN
355: INFO = -10
356: ELSE IF( LDQ.LT.1 .OR. ( WANTZ .AND. LDQ.LT.N ) ) THEN
357: INFO = -12
358: ELSE
359: IF( VALEIG ) THEN
360: IF( N.GT.0 .AND. VU.LE.VL )
361: $ INFO = -14
362: ELSE IF( INDEIG ) THEN
363: IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
364: INFO = -15
365: ELSE IF ( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
366: INFO = -16
367: END IF
368: END IF
369: END IF
370: IF( INFO.EQ.0) THEN
371: IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
372: INFO = -21
373: END IF
374: END IF
375: *
376: IF( INFO.NE.0 ) THEN
377: CALL XERBLA( 'DSBGVX', -INFO )
378: RETURN
379: END IF
380: *
381: * Quick return if possible
382: *
383: M = 0
384: IF( N.EQ.0 )
385: $ RETURN
386: *
387: * Form a split Cholesky factorization of B.
388: *
389: CALL DPBSTF( UPLO, N, KB, BB, LDBB, INFO )
390: IF( INFO.NE.0 ) THEN
391: INFO = N + INFO
392: RETURN
393: END IF
394: *
395: * Transform problem to standard eigenvalue problem.
396: *
397: CALL DSBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ,
398: $ WORK, IINFO )
399: *
400: * Reduce symmetric band matrix to tridiagonal form.
401: *
402: INDD = 1
403: INDE = INDD + N
404: INDWRK = INDE + N
405: IF( WANTZ ) THEN
406: VECT = 'U'
407: ELSE
408: VECT = 'N'
409: END IF
410: CALL DSBTRD( VECT, UPLO, N, KA, AB, LDAB, WORK( INDD ),
411: $ WORK( INDE ), Q, LDQ, WORK( INDWRK ), IINFO )
412: *
413: * If all eigenvalues are desired and ABSTOL is less than or equal
414: * to zero, then call DSTERF or SSTEQR. If this fails for some
415: * eigenvalue, then try DSTEBZ.
416: *
417: TEST = .FALSE.
418: IF( INDEIG ) THEN
419: IF( IL.EQ.1 .AND. IU.EQ.N ) THEN
420: TEST = .TRUE.
421: END IF
422: END IF
423: IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN
424: CALL DCOPY( N, WORK( INDD ), 1, W, 1 )
425: INDEE = INDWRK + 2*N
426: CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
427: IF( .NOT.WANTZ ) THEN
428: CALL DSTERF( N, W, WORK( INDEE ), INFO )
429: ELSE
430: CALL DLACPY( 'A', N, N, Q, LDQ, Z, LDZ )
431: CALL DSTEQR( JOBZ, N, W, WORK( INDEE ), Z, LDZ,
432: $ WORK( INDWRK ), INFO )
433: IF( INFO.EQ.0 ) THEN
434: DO 10 I = 1, N
435: IFAIL( I ) = 0
436: 10 CONTINUE
437: END IF
438: END IF
439: IF( INFO.EQ.0 ) THEN
440: M = N
441: GO TO 30
442: END IF
443: INFO = 0
444: END IF
445: *
446: * Otherwise, call DSTEBZ and, if eigenvectors are desired,
447: * call DSTEIN.
448: *
449: IF( WANTZ ) THEN
450: ORDER = 'B'
451: ELSE
452: ORDER = 'E'
453: END IF
454: INDIBL = 1
455: INDISP = INDIBL + N
456: INDIWO = INDISP + N
457: CALL DSTEBZ( RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL,
458: $ WORK( INDD ), WORK( INDE ), M, NSPLIT, W,
459: $ IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWRK ),
460: $ IWORK( INDIWO ), INFO )
461: *
462: IF( WANTZ ) THEN
463: CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W,
464: $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
465: $ WORK( INDWRK ), IWORK( INDIWO ), IFAIL, INFO )
466: *
467: * Apply transformation matrix used in reduction to tridiagonal
468: * form to eigenvectors returned by DSTEIN.
469: *
470: DO 20 J = 1, M
471: CALL DCOPY( N, Z( 1, J ), 1, WORK( 1 ), 1 )
472: CALL DGEMV( 'N', N, N, ONE, Q, LDQ, WORK, 1, ZERO,
473: $ Z( 1, J ), 1 )
474: 20 CONTINUE
475: END IF
476: *
477: 30 CONTINUE
478: *
479: * If eigenvalues are not in order, then sort them, along with
480: * eigenvectors.
481: *
482: IF( WANTZ ) THEN
483: DO 50 J = 1, M - 1
484: I = 0
485: TMP1 = W( J )
486: DO 40 JJ = J + 1, M
487: IF( W( JJ ).LT.TMP1 ) THEN
488: I = JJ
489: TMP1 = W( JJ )
490: END IF
491: 40 CONTINUE
492: *
493: IF( I.NE.0 ) THEN
494: ITMP1 = IWORK( INDIBL+I-1 )
495: W( I ) = W( J )
496: IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
497: W( J ) = TMP1
498: IWORK( INDIBL+J-1 ) = ITMP1
499: CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
500: IF( INFO.NE.0 ) THEN
501: ITMP1 = IFAIL( I )
502: IFAIL( I ) = IFAIL( J )
503: IFAIL( J ) = ITMP1
504: END IF
505: END IF
506: 50 CONTINUE
507: END IF
508: *
509: RETURN
510: *
511: * End of DSBGVX
512: *
513: END
CVSweb interface <joel.bertrand@systella.fr>