LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsbevx_2stage.f
Go to the documentation of this file.
1*> \brief <b> DSBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* @precisions fortran d -> s
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> \htmlonly
11*> Download DSBEVX_2STAGE + dependencies
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbevx_2stage.f">
13*> [TGZ]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbevx_2stage.f">
15*> [ZIP]</a>
16*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbevx_2stage.f">
17*> [TXT]</a>
18*> \endhtmlonly
19*
20* Definition:
21* ===========
22*
23* SUBROUTINE DSBEVX_2STAGE( JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q,
24* LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
25* LDZ, WORK, LWORK, IWORK, IFAIL, INFO )
26*
27* IMPLICIT NONE
28*
29* .. Scalar Arguments ..
30* CHARACTER JOBZ, RANGE, UPLO
31* INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N, LWORK
32* DOUBLE PRECISION ABSTOL, VL, VU
33* ..
34* .. Array Arguments ..
35* INTEGER IFAIL( * ), IWORK( * )
36* DOUBLE PRECISION AB( LDAB, * ), Q( LDQ, * ), W( * ), WORK( * ),
37* $ Z( LDZ, * )
38* ..
39*
40*
41*> \par Purpose:
42* =============
43*>
44*> \verbatim
45*>
46*> DSBEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
47*> of a real symmetric band matrix A using the 2stage technique for
48*> the reduction to tridiagonal. Eigenvalues and eigenvectors can
49*> be selected by specifying either a range of values or a range of
50*> indices for the desired eigenvalues.
51*> \endverbatim
52*
53* Arguments:
54* ==========
55*
56*> \param[in] JOBZ
57*> \verbatim
58*> JOBZ is CHARACTER*1
59*> = 'N': Compute eigenvalues only;
60*> = 'V': Compute eigenvalues and eigenvectors.
61*> Not available in this release.
62*> \endverbatim
63*>
64*> \param[in] RANGE
65*> \verbatim
66*> RANGE is CHARACTER*1
67*> = 'A': all eigenvalues will be found;
68*> = 'V': all eigenvalues in the half-open interval (VL,VU]
69*> will be found;
70*> = 'I': the IL-th through IU-th eigenvalues will be found.
71*> \endverbatim
72*>
73*> \param[in] UPLO
74*> \verbatim
75*> UPLO is CHARACTER*1
76*> = 'U': Upper triangle of A is stored;
77*> = 'L': Lower triangle of A is stored.
78*> \endverbatim
79*>
80*> \param[in] N
81*> \verbatim
82*> N is INTEGER
83*> The order of the matrix A. N >= 0.
84*> \endverbatim
85*>
86*> \param[in] KD
87*> \verbatim
88*> KD is INTEGER
89*> The number of superdiagonals of the matrix A if UPLO = 'U',
90*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
91*> \endverbatim
92*>
93*> \param[in,out] AB
94*> \verbatim
95*> AB is DOUBLE PRECISION array, dimension (LDAB, N)
96*> On entry, the upper or lower triangle of the symmetric band
97*> matrix A, stored in the first KD+1 rows of the array. The
98*> j-th column of A is stored in the j-th column of the array AB
99*> as follows:
100*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
101*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
102*>
103*> On exit, AB is overwritten by values generated during the
104*> reduction to tridiagonal form. If UPLO = 'U', the first
105*> superdiagonal and the diagonal of the tridiagonal matrix T
106*> are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
107*> the diagonal and first subdiagonal of T are returned in the
108*> first two rows of AB.
109*> \endverbatim
110*>
111*> \param[in] LDAB
112*> \verbatim
113*> LDAB is INTEGER
114*> The leading dimension of the array AB. LDAB >= KD + 1.
115*> \endverbatim
116*>
117*> \param[out] Q
118*> \verbatim
119*> Q is DOUBLE PRECISION array, dimension (LDQ, N)
120*> If JOBZ = 'V', the N-by-N orthogonal matrix used in the
121*> reduction to tridiagonal form.
122*> If JOBZ = 'N', the array Q is not referenced.
123*> \endverbatim
124*>
125*> \param[in] LDQ
126*> \verbatim
127*> LDQ is INTEGER
128*> The leading dimension of the array Q. If JOBZ = 'V', then
129*> LDQ >= max(1,N).
130*> \endverbatim
131*>
132*> \param[in] VL
133*> \verbatim
134*> VL is DOUBLE PRECISION
135*> If RANGE='V', the lower 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] VU
141*> \verbatim
142*> VU is DOUBLE PRECISION
143*> If RANGE='V', the upper bound of the interval to
144*> be searched for eigenvalues. VL < VU.
145*> Not referenced if RANGE = 'A' or 'I'.
146*> \endverbatim
147*>
148*> \param[in] IL
149*> \verbatim
150*> IL is INTEGER
151*> If RANGE='I', the index of the
152*> smallest eigenvalue to be returned.
153*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
154*> Not referenced if RANGE = 'A' or 'V'.
155*> \endverbatim
156*>
157*> \param[in] IU
158*> \verbatim
159*> IU is INTEGER
160*> If RANGE='I', the index of the
161*> largest eigenvalue to be returned.
162*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
163*> Not referenced if RANGE = 'A' or 'V'.
164*> \endverbatim
165*>
166*> \param[in] ABSTOL
167*> \verbatim
168*> ABSTOL is DOUBLE PRECISION
169*> The absolute error tolerance for the eigenvalues.
170*> An approximate eigenvalue is accepted as converged
171*> when it is determined to lie in an interval [a,b]
172*> of width less than or equal to
173*>
174*> ABSTOL + EPS * max( |a|,|b| ) ,
175*>
176*> where EPS is the machine precision. If ABSTOL is less than
177*> or equal to zero, then EPS*|T| will be used in its place,
178*> where |T| is the 1-norm of the tridiagonal matrix obtained
179*> by reducing AB to tridiagonal form.
180*>
181*> Eigenvalues will be computed most accurately when ABSTOL is
182*> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
183*> If this routine returns with INFO>0, indicating that some
184*> eigenvectors did not converge, try setting ABSTOL to
185*> 2*DLAMCH('S').
186*>
187*> See "Computing Small Singular Values of Bidiagonal Matrices
188*> with Guaranteed High Relative Accuracy," by Demmel and
189*> Kahan, LAPACK Working Note #3.
190*> \endverbatim
191*>
192*> \param[out] M
193*> \verbatim
194*> M is INTEGER
195*> The total number of eigenvalues found. 0 <= M <= N.
196*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
197*> \endverbatim
198*>
199*> \param[out] W
200*> \verbatim
201*> W is DOUBLE PRECISION array, dimension (N)
202*> The first M elements contain the selected eigenvalues in
203*> ascending order.
204*> \endverbatim
205*>
206*> \param[out] Z
207*> \verbatim
208*> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
209*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
210*> contain the orthonormal eigenvectors of the matrix A
211*> corresponding to the selected eigenvalues, with the i-th
212*> column of Z holding the eigenvector associated with W(i).
213*> If an eigenvector fails to converge, then that column of Z
214*> contains the latest approximation to the eigenvector, and the
215*> index of the eigenvector is returned in IFAIL.
216*> If JOBZ = 'N', then Z is not referenced.
217*> Note: the user must ensure that at least max(1,M) columns are
218*> supplied in the array Z; if RANGE = 'V', the exact value of M
219*> is not known in advance and an upper bound must be used.
220*> \endverbatim
221*>
222*> \param[in] LDZ
223*> \verbatim
224*> LDZ is INTEGER
225*> The leading dimension of the array Z. LDZ >= 1, and if
226*> JOBZ = 'V', LDZ >= max(1,N).
227*> \endverbatim
228*>
229*> \param[out] WORK
230*> \verbatim
231*> WORK is DOUBLE PRECISION array, dimension (LWORK)
232*> \endverbatim
233*>
234*> \param[in] LWORK
235*> \verbatim
236*> LWORK is INTEGER
237*> The length of the array WORK. LWORK >= 1, when N <= 1;
238*> otherwise
239*> If JOBZ = 'N' and N > 1, LWORK must be queried.
240*> LWORK = MAX(1, 7*N, dimension) where
241*> dimension = (2KD+1)*N + KD*NTHREADS + 2*N
242*> where KD is the size of the band.
243*> NTHREADS is the number of threads used when
244*> openMP compilation is enabled, otherwise =1.
245*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
246*>
247*> If LWORK = -1, then a workspace query is assumed; the routine
248*> only calculates the optimal size of the WORK array, returns
249*> this value as the first entry of the WORK array, and no error
250*> message related to LWORK is issued by XERBLA.
251*> \endverbatim
252*>
253*> \param[out] IWORK
254*> \verbatim
255*> IWORK is INTEGER array, dimension (5*N)
256*> \endverbatim
257*>
258*> \param[out] IFAIL
259*> \verbatim
260*> IFAIL is INTEGER array, dimension (N)
261*> If JOBZ = 'V', then if INFO = 0, the first M elements of
262*> IFAIL are zero. If INFO > 0, then IFAIL contains the
263*> indices of the eigenvectors that failed to converge.
264*> If JOBZ = 'N', then IFAIL is not referenced.
265*> \endverbatim
266*>
267*> \param[out] INFO
268*> \verbatim
269*> INFO is INTEGER
270*> = 0: successful exit.
271*> < 0: if INFO = -i, the i-th argument had an illegal value.
272*> > 0: if INFO = i, then i eigenvectors failed to converge.
273*> Their indices are stored in array IFAIL.
274*> \endverbatim
275*
276* Authors:
277* ========
278*
279*> \author Univ. of Tennessee
280*> \author Univ. of California Berkeley
281*> \author Univ. of Colorado Denver
282*> \author NAG Ltd.
283*
284*> \ingroup hbevx_2stage
285*
286*> \par Further Details:
287* =====================
288*>
289*> \verbatim
290*>
291*> All details about the 2stage techniques are available in:
292*>
293*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
294*> Parallel reduction to condensed forms for symmetric eigenvalue problems
295*> using aggregated fine-grained and memory-aware kernels. In Proceedings
296*> of 2011 International Conference for High Performance Computing,
297*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
298*> Article 8 , 11 pages.
299*> http://doi.acm.org/10.1145/2063384.2063394
300*>
301*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
302*> An improved parallel singular value algorithm and its implementation
303*> for multicore hardware, In Proceedings of 2013 International Conference
304*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
305*> Denver, Colorado, USA, 2013.
306*> Article 90, 12 pages.
307*> http://doi.acm.org/10.1145/2503210.2503292
308*>
309*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
310*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
311*> calculations based on fine-grained memory aware tasks.
312*> International Journal of High Performance Computing Applications.
313*> Volume 28 Issue 2, Pages 196-209, May 2014.
314*> http://hpc.sagepub.com/content/28/2/196
315*>
316*> \endverbatim
317*
318* =====================================================================
319 SUBROUTINE dsbevx_2stage( JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q,
320 $ LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
321 $ LDZ, WORK, LWORK, IWORK, IFAIL, INFO )
322*
323 IMPLICIT NONE
324*
325* -- LAPACK driver routine --
326* -- LAPACK is a software package provided by Univ. of Tennessee, --
327* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
328*
329* .. Scalar Arguments ..
330 CHARACTER JOBZ, RANGE, UPLO
331 INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N, LWORK
332 DOUBLE PRECISION ABSTOL, VL, VU
333* ..
334* .. Array Arguments ..
335 INTEGER IFAIL( * ), IWORK( * )
336 DOUBLE PRECISION AB( LDAB, * ), Q( LDQ, * ), W( * ), WORK( * ),
337 $ z( ldz, * )
338* ..
339*
340* =====================================================================
341*
342* .. Parameters ..
343 DOUBLE PRECISION ZERO, ONE
344 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
345* ..
346* .. Local Scalars ..
347 LOGICAL ALLEIG, INDEIG, LOWER, TEST, VALEIG, WANTZ,
348 $ LQUERY
349 CHARACTER ORDER
350 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
351 $ indisp, indiwo, indwrk, iscale, itmp1, j, jj,
352 $ llwork, lwmin, lhtrd, lwtrd, ib, indhous,
353 $ nsplit
354 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
355 $ SIGMA, SMLNUM, TMP1, VLL, VUU
356* ..
357* .. External Functions ..
358 LOGICAL LSAME
359 INTEGER ILAENV2STAGE
360 DOUBLE PRECISION DLAMCH, DLANSB
361 EXTERNAL lsame, dlamch, dlansb, ilaenv2stage
362* ..
363* .. External Subroutines ..
364 EXTERNAL dcopy, dgemv, dlacpy, dlascl, dscal,
367* ..
368* .. Intrinsic Functions ..
369 INTRINSIC max, min, sqrt
370* ..
371* .. Executable Statements ..
372*
373* Test the input parameters.
374*
375 wantz = lsame( jobz, 'V' )
376 alleig = lsame( range, 'A' )
377 valeig = lsame( range, 'V' )
378 indeig = lsame( range, 'I' )
379 lower = lsame( uplo, 'L' )
380 lquery = ( lwork.EQ.-1 )
381*
382 info = 0
383 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
384 info = -1
385 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
386 info = -2
387 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
388 info = -3
389 ELSE IF( n.LT.0 ) THEN
390 info = -4
391 ELSE IF( kd.LT.0 ) THEN
392 info = -5
393 ELSE IF( ldab.LT.kd+1 ) THEN
394 info = -7
395 ELSE IF( wantz .AND. ldq.LT.max( 1, n ) ) THEN
396 info = -9
397 ELSE
398 IF( valeig ) THEN
399 IF( n.GT.0 .AND. vu.LE.vl )
400 $ info = -11
401 ELSE IF( indeig ) THEN
402 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
403 info = -12
404 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
405 info = -13
406 END IF
407 END IF
408 END IF
409 IF( info.EQ.0 ) THEN
410 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
411 $ info = -18
412 END IF
413*
414 IF( info.EQ.0 ) THEN
415 IF( n.LE.1 ) THEN
416 lwmin = 1
417 work( 1 ) = lwmin
418 ELSE
419 ib = ilaenv2stage( 2, 'DSYTRD_SB2ST', jobz,
420 $ n, kd, -1, -1 )
421 lhtrd = ilaenv2stage( 3, 'DSYTRD_SB2ST', jobz,
422 $ n, kd, ib, -1 )
423 lwtrd = ilaenv2stage( 4, 'DSYTRD_SB2ST', jobz,
424 $ n, kd, ib, -1 )
425 lwmin = 2*n + lhtrd + lwtrd
426 work( 1 ) = lwmin
427 ENDIF
428*
429 IF( lwork.LT.lwmin .AND. .NOT.lquery )
430 $ info = -20
431 END IF
432*
433 IF( info.NE.0 ) THEN
434 CALL xerbla( 'DSBEVX_2STAGE ', -info )
435 RETURN
436 ELSE IF( lquery ) THEN
437 RETURN
438 END IF
439*
440* Quick return if possible
441*
442 m = 0
443 IF( n.EQ.0 )
444 $ RETURN
445*
446 IF( n.EQ.1 ) THEN
447 m = 1
448 IF( lower ) THEN
449 tmp1 = ab( 1, 1 )
450 ELSE
451 tmp1 = ab( kd+1, 1 )
452 END IF
453 IF( valeig ) THEN
454 IF( .NOT.( vl.LT.tmp1 .AND. vu.GE.tmp1 ) )
455 $ m = 0
456 END IF
457 IF( m.EQ.1 ) THEN
458 w( 1 ) = tmp1
459 IF( wantz )
460 $ z( 1, 1 ) = one
461 END IF
462 RETURN
463 END IF
464*
465* Get machine constants.
466*
467 safmin = dlamch( 'Safe minimum' )
468 eps = dlamch( 'Precision' )
469 smlnum = safmin / eps
470 bignum = one / smlnum
471 rmin = sqrt( smlnum )
472 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
473*
474* Scale matrix to allowable range, if necessary.
475*
476 iscale = 0
477 abstll = abstol
478 IF( valeig ) THEN
479 vll = vl
480 vuu = vu
481 ELSE
482 vll = zero
483 vuu = zero
484 END IF
485 anrm = dlansb( 'M', uplo, n, kd, ab, ldab, work )
486 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
487 iscale = 1
488 sigma = rmin / anrm
489 ELSE IF( anrm.GT.rmax ) THEN
490 iscale = 1
491 sigma = rmax / anrm
492 END IF
493 IF( iscale.EQ.1 ) THEN
494 IF( lower ) THEN
495 CALL dlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
496 ELSE
497 CALL dlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
498 END IF
499 IF( abstol.GT.0 )
500 $ abstll = abstol*sigma
501 IF( valeig ) THEN
502 vll = vl*sigma
503 vuu = vu*sigma
504 END IF
505 END IF
506*
507* Call DSYTRD_SB2ST to reduce symmetric band matrix to tridiagonal form.
508*
509 indd = 1
510 inde = indd + n
511 indhous = inde + n
512 indwrk = indhous + lhtrd
513 llwork = lwork - indwrk + 1
514*
515 CALL dsytrd_sb2st( "N", jobz, uplo, n, kd, ab, ldab, work( indd ),
516 $ work( inde ), work( indhous ), lhtrd,
517 $ work( indwrk ), llwork, iinfo )
518*
519* If all eigenvalues are desired and ABSTOL is less than or equal
520* to zero, then call DSTERF or SSTEQR. If this fails for some
521* eigenvalue, then try DSTEBZ.
522*
523 test = .false.
524 IF (indeig) THEN
525 IF (il.EQ.1 .AND. iu.EQ.n) THEN
526 test = .true.
527 END IF
528 END IF
529 IF ((alleig .OR. test) .AND. (abstol.LE.zero)) THEN
530 CALL dcopy( n, work( indd ), 1, w, 1 )
531 indee = indwrk + 2*n
532 IF( .NOT.wantz ) THEN
533 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
534 CALL dsterf( n, w, work( indee ), info )
535 ELSE
536 CALL dlacpy( 'A', n, n, q, ldq, z, ldz )
537 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
538 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
539 $ work( indwrk ), info )
540 IF( info.EQ.0 ) THEN
541 DO 10 i = 1, n
542 ifail( i ) = 0
543 10 CONTINUE
544 END IF
545 END IF
546 IF( info.EQ.0 ) THEN
547 m = n
548 GO TO 30
549 END IF
550 info = 0
551 END IF
552*
553* Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
554*
555 IF( wantz ) THEN
556 order = 'B'
557 ELSE
558 order = 'E'
559 END IF
560 indibl = 1
561 indisp = indibl + n
562 indiwo = indisp + n
563 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
564 $ work( indd ), work( inde ), m, nsplit, w,
565 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
566 $ iwork( indiwo ), info )
567*
568 IF( wantz ) THEN
569 CALL dstein( n, work( indd ), work( inde ), m, w,
570 $ iwork( indibl ), iwork( indisp ), z, ldz,
571 $ work( indwrk ), iwork( indiwo ), ifail, info )
572*
573* Apply orthogonal matrix used in reduction to tridiagonal
574* form to eigenvectors returned by DSTEIN.
575*
576 DO 20 j = 1, m
577 CALL dcopy( n, z( 1, j ), 1, work( 1 ), 1 )
578 CALL dgemv( 'N', n, n, one, q, ldq, work, 1, zero,
579 $ z( 1, j ), 1 )
580 20 CONTINUE
581 END IF
582*
583* If matrix was scaled, then rescale eigenvalues appropriately.
584*
585 30 CONTINUE
586 IF( iscale.EQ.1 ) THEN
587 IF( info.EQ.0 ) THEN
588 imax = m
589 ELSE
590 imax = info - 1
591 END IF
592 CALL dscal( imax, one / sigma, w, 1 )
593 END IF
594*
595* If eigenvalues are not in order, then sort them, along with
596* eigenvectors.
597*
598 IF( wantz ) THEN
599 DO 50 j = 1, m - 1
600 i = 0
601 tmp1 = w( j )
602 DO 40 jj = j + 1, m
603 IF( w( jj ).LT.tmp1 ) THEN
604 i = jj
605 tmp1 = w( jj )
606 END IF
607 40 CONTINUE
608*
609 IF( i.NE.0 ) THEN
610 itmp1 = iwork( indibl+i-1 )
611 w( i ) = w( j )
612 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
613 w( j ) = tmp1
614 iwork( indibl+j-1 ) = itmp1
615 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
616 IF( info.NE.0 ) THEN
617 itmp1 = ifail( i )
618 ifail( i ) = ifail( j )
619 ifail( j ) = itmp1
620 END IF
621 END IF
622 50 CONTINUE
623 END IF
624*
625* Set WORK(1) to optimal workspace size.
626*
627 work( 1 ) = lwmin
628*
629 RETURN
630*
631* End of DSBEVX_2STAGE
632*
633 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_2stage(jobz, range, uplo, n, kd, ab, ldab, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, iwork, ifail, info)
DSBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
subroutine dsytrd_sb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
DSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric tridiagonal form T
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
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:143
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:273
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:174
subroutine dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:131
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82