LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
dsbevx.f
Go to the documentation of this file.
1*> \brief <b> DSBEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DSBEVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbevx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbevx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbevx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DSBEVX( JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL,
20* VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK,
21* IFAIL, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER JOBZ, RANGE, UPLO
25* INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N
26* DOUBLE PRECISION ABSTOL, VL, VU
27* ..
28* .. Array Arguments ..
29* INTEGER IFAIL( * ), IWORK( * )
30* DOUBLE PRECISION AB( LDAB, * ), Q( LDQ, * ), W( * ), WORK( * ),
31* $ Z( LDZ, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> DSBEVX computes selected eigenvalues and, optionally, eigenvectors
41*> of a real symmetric band matrix A. Eigenvalues and eigenvectors can
42*> be selected by specifying either a range of values or a range of
43*> indices for the desired eigenvalues.
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] JOBZ
50*> \verbatim
51*> JOBZ is CHARACTER*1
52*> = 'N': Compute eigenvalues only;
53*> = 'V': Compute eigenvalues and eigenvectors.
54*> \endverbatim
55*>
56*> \param[in] RANGE
57*> \verbatim
58*> RANGE is CHARACTER*1
59*> = 'A': all eigenvalues will be found;
60*> = 'V': all eigenvalues in the half-open interval (VL,VU]
61*> will be found;
62*> = 'I': the IL-th through IU-th eigenvalues will be found.
63*> \endverbatim
64*>
65*> \param[in] UPLO
66*> \verbatim
67*> UPLO is CHARACTER*1
68*> = 'U': Upper triangle of A is stored;
69*> = 'L': Lower triangle of A is stored.
70*> \endverbatim
71*>
72*> \param[in] N
73*> \verbatim
74*> N is INTEGER
75*> The order of the matrix A. N >= 0.
76*> \endverbatim
77*>
78*> \param[in] KD
79*> \verbatim
80*> KD is INTEGER
81*> The number of superdiagonals of the matrix A if UPLO = 'U',
82*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
83*> \endverbatim
84*>
85*> \param[in,out] AB
86*> \verbatim
87*> AB is DOUBLE PRECISION array, dimension (LDAB, N)
88*> On entry, the upper or lower triangle of the symmetric band
89*> matrix A, stored in the first KD+1 rows of the array. The
90*> j-th column of A is stored in the j-th column of the array AB
91*> as follows:
92*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
93*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
94*>
95*> On exit, AB is overwritten by values generated during the
96*> reduction to tridiagonal form. If UPLO = 'U', the first
97*> superdiagonal and the diagonal of the tridiagonal matrix T
98*> are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
99*> the diagonal and first subdiagonal of T are returned in the
100*> first two rows of AB.
101*> \endverbatim
102*>
103*> \param[in] LDAB
104*> \verbatim
105*> LDAB is INTEGER
106*> The leading dimension of the array AB. LDAB >= KD + 1.
107*> \endverbatim
108*>
109*> \param[out] Q
110*> \verbatim
111*> Q is DOUBLE PRECISION array, dimension (LDQ, N)
112*> If JOBZ = 'V', the N-by-N orthogonal matrix used in the
113*> reduction to tridiagonal form.
114*> If JOBZ = 'N', the array Q is not referenced.
115*> \endverbatim
116*>
117*> \param[in] LDQ
118*> \verbatim
119*> LDQ is INTEGER
120*> The leading dimension of the array Q. If JOBZ = 'V', then
121*> LDQ >= max(1,N).
122*> \endverbatim
123*>
124*> \param[in] VL
125*> \verbatim
126*> VL is DOUBLE PRECISION
127*> If RANGE='V', the lower bound of the interval to
128*> be searched for eigenvalues. VL < VU.
129*> Not referenced if RANGE = 'A' or 'I'.
130*> \endverbatim
131*>
132*> \param[in] VU
133*> \verbatim
134*> VU is DOUBLE PRECISION
135*> If RANGE='V', the upper bound of the interval to
136*> be searched for eigenvalues. VL < VU.
137*> Not referenced if RANGE = 'A' or 'I'.
138*> \endverbatim
139*>
140*> \param[in] IL
141*> \verbatim
142*> IL is INTEGER
143*> If RANGE='I', the index of the
144*> smallest eigenvalue to be returned.
145*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
146*> Not referenced if RANGE = 'A' or 'V'.
147*> \endverbatim
148*>
149*> \param[in] IU
150*> \verbatim
151*> IU is INTEGER
152*> If RANGE='I', the index of the
153*> largest eigenvalue to be returned.
154*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
155*> Not referenced if RANGE = 'A' or 'V'.
156*> \endverbatim
157*>
158*> \param[in] ABSTOL
159*> \verbatim
160*> ABSTOL is DOUBLE PRECISION
161*> The absolute error tolerance for the eigenvalues.
162*> An approximate eigenvalue is accepted as converged
163*> when it is determined to lie in an interval [a,b]
164*> of width less than or equal to
165*>
166*> ABSTOL + EPS * max( |a|,|b| ) ,
167*>
168*> where EPS is the machine precision. If ABSTOL is less than
169*> or equal to zero, then EPS*|T| will be used in its place,
170*> where |T| is the 1-norm of the tridiagonal matrix obtained
171*> by reducing AB to tridiagonal form.
172*>
173*> Eigenvalues will be computed most accurately when ABSTOL is
174*> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
175*> If this routine returns with INFO>0, indicating that some
176*> eigenvectors did not converge, try setting ABSTOL to
177*> 2*DLAMCH('S').
178*>
179*> See "Computing Small Singular Values of Bidiagonal Matrices
180*> with Guaranteed High Relative Accuracy," by Demmel and
181*> Kahan, LAPACK Working Note #3.
182*> \endverbatim
183*>
184*> \param[out] M
185*> \verbatim
186*> M is INTEGER
187*> The total number of eigenvalues found. 0 <= M <= N.
188*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
189*> \endverbatim
190*>
191*> \param[out] W
192*> \verbatim
193*> W is DOUBLE PRECISION array, dimension (N)
194*> The first M elements contain the selected eigenvalues in
195*> ascending order.
196*> \endverbatim
197*>
198*> \param[out] Z
199*> \verbatim
200*> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
201*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
202*> contain the orthonormal eigenvectors of the matrix A
203*> corresponding to the selected eigenvalues, with the i-th
204*> column of Z holding the eigenvector associated with W(i).
205*> If an eigenvector fails to converge, then that column of Z
206*> contains the latest approximation to the eigenvector, and the
207*> index of the eigenvector is returned in IFAIL.
208*> If JOBZ = 'N', then Z is not referenced.
209*> Note: the user must ensure that at least max(1,M) columns are
210*> supplied in the array Z; if RANGE = 'V', the exact value of M
211*> is not known in advance and an upper bound must be used.
212*> \endverbatim
213*>
214*> \param[in] LDZ
215*> \verbatim
216*> LDZ is INTEGER
217*> The leading dimension of the array Z. LDZ >= 1, and if
218*> JOBZ = 'V', LDZ >= max(1,N).
219*> \endverbatim
220*>
221*> \param[out] WORK
222*> \verbatim
223*> WORK is DOUBLE PRECISION array, dimension (7*N)
224*> \endverbatim
225*>
226*> \param[out] IWORK
227*> \verbatim
228*> IWORK is INTEGER array, dimension (5*N)
229*> \endverbatim
230*>
231*> \param[out] IFAIL
232*> \verbatim
233*> IFAIL is INTEGER array, dimension (N)
234*> If JOBZ = 'V', then if INFO = 0, the first M elements of
235*> IFAIL are zero. If INFO > 0, then IFAIL contains the
236*> indices of the eigenvectors that failed to converge.
237*> If JOBZ = 'N', then IFAIL is not referenced.
238*> \endverbatim
239*>
240*> \param[out] INFO
241*> \verbatim
242*> INFO is INTEGER
243*> = 0: successful exit.
244*> < 0: if INFO = -i, the i-th argument had an illegal value.
245*> > 0: if INFO = i, then i eigenvectors failed to converge.
246*> Their indices are stored in array IFAIL.
247*> \endverbatim
248*
249* Authors:
250* ========
251*
252*> \author Univ. of Tennessee
253*> \author Univ. of California Berkeley
254*> \author Univ. of Colorado Denver
255*> \author NAG Ltd.
256*
257*> \ingroup hbevx
258*
259* =====================================================================
260 SUBROUTINE dsbevx( JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ,
261 $ VL,
262 $ VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK,
263 $ IFAIL, INFO )
264*
265* -- LAPACK driver routine --
266* -- LAPACK is a software package provided by Univ. of Tennessee, --
267* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
268*
269* .. Scalar Arguments ..
270 CHARACTER JOBZ, RANGE, UPLO
271 INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N
272 DOUBLE PRECISION ABSTOL, VL, VU
273* ..
274* .. Array Arguments ..
275 INTEGER IFAIL( * ), IWORK( * )
276 DOUBLE PRECISION AB( LDAB, * ), Q( LDQ, * ), W( * ), WORK( * ),
277 $ Z( LDZ, * )
278* ..
279*
280* =====================================================================
281*
282* .. Parameters ..
283 DOUBLE PRECISION ZERO, ONE
284 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
285* ..
286* .. Local Scalars ..
287 LOGICAL ALLEIG, INDEIG, LOWER, TEST, VALEIG, WANTZ
288 CHARACTER ORDER
289 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
290 $ indisp, indiwo, indwrk, iscale, itmp1, j, jj,
291 $ nsplit
292 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
293 $ SIGMA, SMLNUM, TMP1, VLL, VUU
294* ..
295* .. External Functions ..
296 LOGICAL LSAME
297 DOUBLE PRECISION DLAMCH, DLANSB
298 EXTERNAL LSAME, DLAMCH, DLANSB
299* ..
300* .. External Subroutines ..
301 EXTERNAL dcopy, dgemv, dlacpy, dlascl, dsbtrd,
302 $ dscal,
304* ..
305* .. Intrinsic Functions ..
306 INTRINSIC max, min, sqrt
307* ..
308* .. Executable Statements ..
309*
310* Test the input parameters.
311*
312 wantz = lsame( jobz, 'V' )
313 alleig = lsame( range, 'A' )
314 valeig = lsame( range, 'V' )
315 indeig = lsame( range, 'I' )
316 lower = lsame( uplo, 'L' )
317*
318 info = 0
319 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
320 info = -1
321 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
322 info = -2
323 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
324 info = -3
325 ELSE IF( n.LT.0 ) THEN
326 info = -4
327 ELSE IF( kd.LT.0 ) THEN
328 info = -5
329 ELSE IF( ldab.LT.kd+1 ) THEN
330 info = -7
331 ELSE IF( wantz .AND. ldq.LT.max( 1, n ) ) THEN
332 info = -9
333 ELSE
334 IF( valeig ) THEN
335 IF( n.GT.0 .AND. vu.LE.vl )
336 $ info = -11
337 ELSE IF( indeig ) THEN
338 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
339 info = -12
340 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
341 info = -13
342 END IF
343 END IF
344 END IF
345 IF( info.EQ.0 ) THEN
346 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
347 $ info = -18
348 END IF
349*
350 IF( info.NE.0 ) THEN
351 CALL xerbla( 'DSBEVX', -info )
352 RETURN
353 END IF
354*
355* Quick return if possible
356*
357 m = 0
358 IF( n.EQ.0 )
359 $ RETURN
360*
361 IF( n.EQ.1 ) THEN
362 m = 1
363 IF( lower ) THEN
364 tmp1 = ab( 1, 1 )
365 ELSE
366 tmp1 = ab( kd+1, 1 )
367 END IF
368 IF( valeig ) THEN
369 IF( .NOT.( vl.LT.tmp1 .AND. vu.GE.tmp1 ) )
370 $ m = 0
371 END IF
372 IF( m.EQ.1 ) THEN
373 w( 1 ) = tmp1
374 IF( wantz )
375 $ z( 1, 1 ) = one
376 END IF
377 RETURN
378 END IF
379*
380* Get machine constants.
381*
382 safmin = dlamch( 'Safe minimum' )
383 eps = dlamch( 'Precision' )
384 smlnum = safmin / eps
385 bignum = one / smlnum
386 rmin = sqrt( smlnum )
387 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
388*
389* Scale matrix to allowable range, if necessary.
390*
391 iscale = 0
392 abstll = abstol
393 IF( valeig ) THEN
394 vll = vl
395 vuu = vu
396 ELSE
397 vll = zero
398 vuu = zero
399 END IF
400 anrm = dlansb( 'M', uplo, n, kd, ab, ldab, work )
401 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
402 iscale = 1
403 sigma = rmin / anrm
404 ELSE IF( anrm.GT.rmax ) THEN
405 iscale = 1
406 sigma = rmax / anrm
407 END IF
408 IF( iscale.EQ.1 ) THEN
409 IF( lower ) THEN
410 CALL dlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab,
411 $ info )
412 ELSE
413 CALL dlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab,
414 $ info )
415 END IF
416 IF( abstol.GT.0 )
417 $ abstll = abstol*sigma
418 IF( valeig ) THEN
419 vll = vl*sigma
420 vuu = vu*sigma
421 END IF
422 END IF
423*
424* Call DSBTRD to reduce symmetric band matrix to tridiagonal form.
425*
426 indd = 1
427 inde = indd + n
428 indwrk = inde + n
429 CALL dsbtrd( jobz, uplo, n, kd, ab, ldab, work( indd ),
430 $ work( inde ), q, ldq, work( indwrk ), iinfo )
431*
432* If all eigenvalues are desired and ABSTOL is less than or equal
433* to zero, then call DSTERF or SSTEQR. If this fails for some
434* eigenvalue, then try DSTEBZ.
435*
436 test = .false.
437 IF (indeig) THEN
438 IF (il.EQ.1 .AND. iu.EQ.n) THEN
439 test = .true.
440 END IF
441 END IF
442 IF ((alleig .OR. test) .AND. (abstol.LE.zero)) THEN
443 CALL dcopy( n, work( indd ), 1, w, 1 )
444 indee = indwrk + 2*n
445 IF( .NOT.wantz ) THEN
446 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
447 CALL dsterf( n, w, work( indee ), info )
448 ELSE
449 CALL dlacpy( 'A', n, n, q, ldq, z, ldz )
450 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
451 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
452 $ work( indwrk ), info )
453 IF( info.EQ.0 ) THEN
454 DO 10 i = 1, n
455 ifail( i ) = 0
456 10 CONTINUE
457 END IF
458 END IF
459 IF( info.EQ.0 ) THEN
460 m = n
461 GO TO 30
462 END IF
463 info = 0
464 END IF
465*
466* Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
467*
468 IF( wantz ) THEN
469 order = 'B'
470 ELSE
471 order = 'E'
472 END IF
473 indibl = 1
474 indisp = indibl + n
475 indiwo = indisp + n
476 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
477 $ work( indd ), work( inde ), m, nsplit, w,
478 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
479 $ iwork( indiwo ), info )
480*
481 IF( wantz ) THEN
482 CALL dstein( n, work( indd ), work( inde ), m, w,
483 $ iwork( indibl ), iwork( indisp ), z, ldz,
484 $ work( indwrk ), iwork( indiwo ), ifail, info )
485*
486* Apply orthogonal matrix used in reduction to tridiagonal
487* form to eigenvectors returned by DSTEIN.
488*
489 DO 20 j = 1, m
490 CALL dcopy( n, z( 1, j ), 1, work( 1 ), 1 )
491 CALL dgemv( 'N', n, n, one, q, ldq, work, 1, zero,
492 $ z( 1, j ), 1 )
493 20 CONTINUE
494 END IF
495*
496* If matrix was scaled, then rescale eigenvalues appropriately.
497*
498 30 CONTINUE
499 IF( iscale.EQ.1 ) THEN
500 IF( info.EQ.0 ) THEN
501 imax = m
502 ELSE
503 imax = info - 1
504 END IF
505 CALL dscal( imax, one / sigma, w, 1 )
506 END IF
507*
508* If eigenvalues are not in order, then sort them, along with
509* eigenvectors.
510*
511 IF( wantz ) THEN
512 DO 50 j = 1, m - 1
513 i = 0
514 tmp1 = w( j )
515 DO 40 jj = j + 1, m
516 IF( w( jj ).LT.tmp1 ) THEN
517 i = jj
518 tmp1 = w( jj )
519 END IF
520 40 CONTINUE
521*
522 IF( i.NE.0 ) THEN
523 itmp1 = iwork( indibl+i-1 )
524 w( i ) = w( j )
525 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
526 w( j ) = tmp1
527 iwork( indibl+j-1 ) = itmp1
528 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
529 IF( info.NE.0 ) THEN
530 itmp1 = ifail( i )
531 ifail( i ) = ifail( j )
532 ifail( j ) = itmp1
533 END IF
534 END IF
535 50 CONTINUE
536 END IF
537*
538 RETURN
539*
540* End of DSBEVX
541*
542 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
subroutine dsbevx(jobz, range, uplo, n, kd, ab, ldab, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, iwork, ifail, info)
DSBEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition dsbevx.f:264
subroutine dsbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
DSBTRD
Definition dsbtrd.f:161
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:142
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
DSTEBZ
Definition dstebz.f:272
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:172
subroutine dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:129
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82