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