001:       SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
002:      $                   LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
003:      $                   LDZ, WORK, IWORK, IFAIL, INFO )
004: *
005: *  -- LAPACK driver routine (version 3.2) --
006: *     Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
007: *     November 2006
008: *
009: *     .. Scalar Arguments ..
010:       CHARACTER          JOBZ, RANGE, UPLO
011:       INTEGER            IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
012:      $                   N
013:       DOUBLE PRECISION   ABSTOL, VL, VU
014: *     ..
015: *     .. Array Arguments ..
016:       INTEGER            IFAIL( * ), IWORK( * )
017:       DOUBLE PRECISION   AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
018:      $                   W( * ), WORK( * ), Z( LDZ, * )
019: *     ..
020: *
021: *  Purpose
022: *  =======
023: *
024: *  DSBGVX computes selected eigenvalues, and optionally, eigenvectors
025: *  of a real generalized symmetric-definite banded eigenproblem, of
026: *  the form A*x=(lambda)*B*x.  Here A and B are assumed to be symmetric
027: *  and banded, and B is also positive definite.  Eigenvalues and
028: *  eigenvectors can be selected by specifying either all eigenvalues,
029: *  a range of values or a range of indices for the desired eigenvalues.
030: *
031: *  Arguments
032: *  =========
033: *
034: *  JOBZ    (input) CHARACTER*1
035: *          = 'N':  Compute eigenvalues only;
036: *          = 'V':  Compute eigenvalues and eigenvectors.
037: *
038: *  RANGE   (input) CHARACTER*1
039: *          = 'A': all eigenvalues will be found.
040: *          = 'V': all eigenvalues in the half-open interval (VL,VU]
041: *                 will be found.
042: *          = 'I': the IL-th through IU-th eigenvalues will be found.
043: *
044: *  UPLO    (input) CHARACTER*1
045: *          = 'U':  Upper triangles of A and B are stored;
046: *          = 'L':  Lower triangles of A and B are stored.
047: *
048: *  N       (input) INTEGER
049: *          The order of the matrices A and B.  N >= 0.
050: *
051: *  KA      (input) INTEGER
052: *          The number of superdiagonals of the matrix A if UPLO = 'U',
053: *          or the number of subdiagonals if UPLO = 'L'.  KA >= 0.
054: *
055: *  KB      (input) INTEGER
056: *          The number of superdiagonals of the matrix B if UPLO = 'U',
057: *          or the number of subdiagonals if UPLO = 'L'.  KB >= 0.
058: *
059: *  AB      (input/output) DOUBLE PRECISION array, dimension (LDAB, N)
060: *          On entry, the upper or lower triangle of the symmetric band
061: *          matrix A, stored in the first ka+1 rows of the array.  The
062: *          j-th column of A is stored in the j-th column of the array AB
063: *          as follows:
064: *          if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
065: *          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+ka).
066: *
067: *          On exit, the contents of AB are destroyed.
068: *
069: *  LDAB    (input) INTEGER
070: *          The leading dimension of the array AB.  LDAB >= KA+1.
071: *
072: *  BB      (input/output) DOUBLE PRECISION array, dimension (LDBB, N)
073: *          On entry, the upper or lower triangle of the symmetric band
074: *          matrix B, stored in the first kb+1 rows of the array.  The
075: *          j-th column of B is stored in the j-th column of the array BB
076: *          as follows:
077: *          if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
078: *          if UPLO = 'L', BB(1+i-j,j)    = B(i,j) for j<=i<=min(n,j+kb).
079: *
080: *          On exit, the factor S from the split Cholesky factorization
081: *          B = S**T*S, as returned by DPBSTF.
082: *
083: *  LDBB    (input) INTEGER
084: *          The leading dimension of the array BB.  LDBB >= KB+1.
085: *
086: *  Q       (output) DOUBLE PRECISION array, dimension (LDQ, N)
087: *          If JOBZ = 'V', the n-by-n matrix used in the reduction of
088: *          A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
089: *          and consequently C to tridiagonal form.
090: *          If JOBZ = 'N', the array Q is not referenced.
091: *
092: *  LDQ     (input) INTEGER
093: *          The leading dimension of the array Q.  If JOBZ = 'N',
094: *          LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).
095: *
096: *  VL      (input) DOUBLE PRECISION
097: *  VU      (input) DOUBLE PRECISION
098: *          If RANGE='V', the lower and upper bounds of the interval to
099: *          be searched for eigenvalues. VL < VU.
100: *          Not referenced if RANGE = 'A' or 'I'.
101: *
102: *  IL      (input) INTEGER
103: *  IU      (input) INTEGER
104: *          If RANGE='I', the indices (in ascending order) of the
105: *          smallest and largest eigenvalues to be returned.
106: *          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
107: *          Not referenced if RANGE = 'A' or 'V'.
108: *
109: *  ABSTOL  (input) DOUBLE PRECISION
110: *          The absolute error tolerance for the eigenvalues.
111: *          An approximate eigenvalue is accepted as converged
112: *          when it is determined to lie in an interval [a,b]
113: *          of width less than or equal to
114: *
115: *                  ABSTOL + EPS *   max( |a|,|b| ) ,
116: *
117: *          where EPS is the machine precision.  If ABSTOL is less than
118: *          or equal to zero, then  EPS*|T|  will be used in its place,
119: *          where |T| is the 1-norm of the tridiagonal matrix obtained
120: *          by reducing A to tridiagonal form.
121: *
122: *          Eigenvalues will be computed most accurately when ABSTOL is
123: *          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
124: *          If this routine returns with INFO>0, indicating that some
125: *          eigenvectors did not converge, try setting ABSTOL to
126: *          2*DLAMCH('S').
127: *
128: *  M       (output) INTEGER
129: *          The total number of eigenvalues found.  0 <= M <= N.
130: *          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
131: *
132: *  W       (output) DOUBLE PRECISION array, dimension (N)
133: *          If INFO = 0, the eigenvalues in ascending order.
134: *
135: *  Z       (output) DOUBLE PRECISION array, dimension (LDZ, N)
136: *          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
137: *          eigenvectors, with the i-th column of Z holding the
138: *          eigenvector associated with W(i).  The eigenvectors are
139: *          normalized so Z**T*B*Z = I.
140: *          If JOBZ = 'N', then Z is not referenced.
141: *
142: *  LDZ     (input) INTEGER
143: *          The leading dimension of the array Z.  LDZ >= 1, and if
144: *          JOBZ = 'V', LDZ >= max(1,N).
145: *
146: *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (7*N)
147: *
148: *  IWORK   (workspace/output) INTEGER array, dimension (5*N)
149: *
150: *  IFAIL   (output) INTEGER array, dimension (M)
151: *          If JOBZ = 'V', then if INFO = 0, the first M elements of
152: *          IFAIL are zero.  If INFO > 0, then IFAIL contains the
153: *          indices of the eigenvalues that failed to converge.
154: *          If JOBZ = 'N', then IFAIL is not referenced.
155: *
156: *  INFO    (output) INTEGER
157: *          = 0 : successful exit
158: *          < 0 : if INFO = -i, the i-th argument had an illegal value
159: *          <= N: if INFO = i, then i eigenvectors failed to converge.
160: *                  Their indices are stored in IFAIL.
161: *          > N : DPBSTF returned an error code; i.e.,
162: *                if INFO = N + i, for 1 <= i <= N, then the leading
163: *                minor of order i of B is not positive definite.
164: *                The factorization of B could not be completed and
165: *                no eigenvalues or eigenvectors were computed.
166: *
167: *  Further Details
168: *  ===============
169: *
170: *  Based on contributions by
171: *     Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
172: *
173: *  =====================================================================
174: *
175: *     .. Parameters ..
176:       DOUBLE PRECISION   ZERO, ONE
177:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
178: *     ..
179: *     .. Local Scalars ..
180:       LOGICAL            ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
181:       CHARACTER          ORDER, VECT
182:       INTEGER            I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP,
183:      $                   INDIWO, INDWRK, ITMP1, J, JJ, NSPLIT
184:       DOUBLE PRECISION   TMP1
185: *     ..
186: *     .. External Functions ..
187:       LOGICAL            LSAME
188:       EXTERNAL           LSAME
189: *     ..
190: *     .. External Subroutines ..
191:       EXTERNAL           DCOPY, DGEMV, DLACPY, DPBSTF, DSBGST, DSBTRD,
192:      $                   DSTEBZ, DSTEIN, DSTEQR, DSTERF, DSWAP, XERBLA
193: *     ..
194: *     .. Intrinsic Functions ..
195:       INTRINSIC          MIN
196: *     ..
197: *     .. Executable Statements ..
198: *
199: *     Test the input parameters.
200: *
201:       WANTZ = LSAME( JOBZ, 'V' )
202:       UPPER = LSAME( UPLO, 'U' )
203:       ALLEIG = LSAME( RANGE, 'A' )
204:       VALEIG = LSAME( RANGE, 'V' )
205:       INDEIG = LSAME( RANGE, 'I' )
206: *
207:       INFO = 0
208:       IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
209:          INFO = -1
210:       ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
211:          INFO = -2
212:       ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
213:          INFO = -3
214:       ELSE IF( N.LT.0 ) THEN
215:          INFO = -4
216:       ELSE IF( KA.LT.0 ) THEN
217:          INFO = -5
218:       ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN
219:          INFO = -6
220:       ELSE IF( LDAB.LT.KA+1 ) THEN
221:          INFO = -8
222:       ELSE IF( LDBB.LT.KB+1 ) THEN
223:          INFO = -10
224:       ELSE IF( LDQ.LT.1 .OR. ( WANTZ .AND. LDQ.LT.N ) ) THEN
225:          INFO = -12
226:       ELSE
227:          IF( VALEIG ) THEN
228:             IF( N.GT.0 .AND. VU.LE.VL )
229:      $         INFO = -14
230:          ELSE IF( INDEIG ) THEN
231:             IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
232:                INFO = -15
233:             ELSE IF ( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
234:                INFO = -16
235:             END IF
236:          END IF
237:       END IF
238:       IF( INFO.EQ.0) THEN
239:          IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
240:             INFO = -21
241:          END IF
242:       END IF
243: *
244:       IF( INFO.NE.0 ) THEN
245:          CALL XERBLA( 'DSBGVX', -INFO )
246:          RETURN
247:       END IF
248: *
249: *     Quick return if possible
250: *
251:       M = 0
252:       IF( N.EQ.0 )
253:      $   RETURN
254: *
255: *     Form a split Cholesky factorization of B.
256: *
257:       CALL DPBSTF( UPLO, N, KB, BB, LDBB, INFO )
258:       IF( INFO.NE.0 ) THEN
259:          INFO = N + INFO
260:          RETURN
261:       END IF
262: *
263: *     Transform problem to standard eigenvalue problem.
264: *
265:       CALL DSBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ,
266:      $             WORK, IINFO )
267: *
268: *     Reduce symmetric band matrix to tridiagonal form.
269: *
270:       INDD = 1
271:       INDE = INDD + N
272:       INDWRK = INDE + N
273:       IF( WANTZ ) THEN
274:          VECT = 'U'
275:       ELSE
276:          VECT = 'N'
277:       END IF
278:       CALL DSBTRD( VECT, UPLO, N, KA, AB, LDAB, WORK( INDD ),
279:      $             WORK( INDE ), Q, LDQ, WORK( INDWRK ), IINFO )
280: *
281: *     If all eigenvalues are desired and ABSTOL is less than or equal
282: *     to zero, then call DSTERF or SSTEQR.  If this fails for some
283: *     eigenvalue, then try DSTEBZ.
284: *
285:       TEST = .FALSE.
286:       IF( INDEIG ) THEN
287:          IF( IL.EQ.1 .AND. IU.EQ.N ) THEN
288:             TEST = .TRUE.
289:          END IF
290:       END IF
291:       IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN
292:          CALL DCOPY( N, WORK( INDD ), 1, W, 1 )
293:          INDEE = INDWRK + 2*N
294:          CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
295:          IF( .NOT.WANTZ ) THEN
296:             CALL DSTERF( N, W, WORK( INDEE ), INFO )
297:          ELSE
298:             CALL DLACPY( 'A', N, N, Q, LDQ, Z, LDZ )
299:             CALL DSTEQR( JOBZ, N, W, WORK( INDEE ), Z, LDZ,
300:      $                   WORK( INDWRK ), INFO )
301:             IF( INFO.EQ.0 ) THEN
302:                DO 10 I = 1, N
303:                   IFAIL( I ) = 0
304:    10          CONTINUE
305:             END IF
306:          END IF
307:          IF( INFO.EQ.0 ) THEN
308:             M = N
309:             GO TO 30
310:          END IF
311:          INFO = 0
312:       END IF
313: *
314: *     Otherwise, call DSTEBZ and, if eigenvectors are desired,
315: *     call DSTEIN.
316: *
317:       IF( WANTZ ) THEN
318:          ORDER = 'B'
319:       ELSE
320:          ORDER = 'E'
321:       END IF
322:       INDIBL = 1
323:       INDISP = INDIBL + N
324:       INDIWO = INDISP + N
325:       CALL DSTEBZ( RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL,
326:      $             WORK( INDD ), WORK( INDE ), M, NSPLIT, W,
327:      $             IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWRK ),
328:      $             IWORK( INDIWO ), INFO )
329: *
330:       IF( WANTZ ) THEN
331:          CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W,
332:      $                IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
333:      $                WORK( INDWRK ), IWORK( INDIWO ), IFAIL, INFO )
334: *
335: *        Apply transformation matrix used in reduction to tridiagonal
336: *        form to eigenvectors returned by DSTEIN.
337: *
338:          DO 20 J = 1, M
339:             CALL DCOPY( N, Z( 1, J ), 1, WORK( 1 ), 1 )
340:             CALL DGEMV( 'N', N, N, ONE, Q, LDQ, WORK, 1, ZERO,
341:      $                  Z( 1, J ), 1 )
342:    20    CONTINUE
343:       END IF
344: *
345:    30 CONTINUE
346: *
347: *     If eigenvalues are not in order, then sort them, along with
348: *     eigenvectors.
349: *
350:       IF( WANTZ ) THEN
351:          DO 50 J = 1, M - 1
352:             I = 0
353:             TMP1 = W( J )
354:             DO 40 JJ = J + 1, M
355:                IF( W( JJ ).LT.TMP1 ) THEN
356:                   I = JJ
357:                   TMP1 = W( JJ )
358:                END IF
359:    40       CONTINUE
360: *
361:             IF( I.NE.0 ) THEN
362:                ITMP1 = IWORK( INDIBL+I-1 )
363:                W( I ) = W( J )
364:                IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
365:                W( J ) = TMP1
366:                IWORK( INDIBL+J-1 ) = ITMP1
367:                CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
368:                IF( INFO.NE.0 ) THEN
369:                   ITMP1 = IFAIL( I )
370:                   IFAIL( I ) = IFAIL( J )
371:                   IFAIL( J ) = ITMP1
372:                END IF
373:             END IF
374:    50    CONTINUE
375:       END IF
376: *
377:       RETURN
378: *
379: *     End of DSBGVX
380: *
381:       END
382: