LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
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 *> \date June 2016
285 *
286 *> \ingroup doubleOTHEReigen
287 *
288 *> \par Further Details:
289 * =====================
290 *>
291 *> \verbatim
292 *>
293 *> All details about the 2stage techniques are available in:
294 *>
295 *> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
296 *> Parallel reduction to condensed forms for symmetric eigenvalue problems
297 *> using aggregated fine-grained and memory-aware kernels. In Proceedings
298 *> of 2011 International Conference for High Performance Computing,
299 *> Networking, Storage and Analysis (SC '11), New York, NY, USA,
300 *> Article 8 , 11 pages.
301 *> http://doi.acm.org/10.1145/2063384.2063394
302 *>
303 *> A. Haidar, J. Kurzak, P. Luszczek, 2013.
304 *> An improved parallel singular value algorithm and its implementation
305 *> for multicore hardware, In Proceedings of 2013 International Conference
306 *> for High Performance Computing, Networking, Storage and Analysis (SC '13).
307 *> Denver, Colorado, USA, 2013.
308 *> Article 90, 12 pages.
309 *> http://doi.acm.org/10.1145/2503210.2503292
310 *>
311 *> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
312 *> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
313 *> calculations based on fine-grained memory aware tasks.
314 *> International Journal of High Performance Computing Applications.
315 *> Volume 28 Issue 2, Pages 196-209, May 2014.
316 *> http://hpc.sagepub.com/content/28/2/196
317 *>
318 *> \endverbatim
319 *
320 * =====================================================================
321  SUBROUTINE dsbevx_2stage( JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q,
322  $ ldq, vl, vu, il, iu, abstol, m, w, z,
323  $ ldz, work, lwork, iwork, ifail, info )
324 *
325  IMPLICIT NONE
326 *
327 * -- LAPACK driver routine (version 3.7.0) --
328 * -- LAPACK is a software package provided by Univ. of Tennessee, --
329 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
330 * June 2016
331 *
332 * .. Scalar Arguments ..
333  CHARACTER JOBZ, RANGE, UPLO
334  INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N, LWORK
335  DOUBLE PRECISION ABSTOL, VL, VU
336 * ..
337 * .. Array Arguments ..
338  INTEGER IFAIL( * ), IWORK( * )
339  DOUBLE PRECISION AB( ldab, * ), Q( ldq, * ), W( * ), WORK( * ),
340  $ z( ldz, * )
341 * ..
342 *
343 * =====================================================================
344 *
345 * .. Parameters ..
346  DOUBLE PRECISION ZERO, ONE
347  parameter ( zero = 0.0d0, one = 1.0d0 )
348 * ..
349 * .. Local Scalars ..
350  LOGICAL ALLEIG, INDEIG, LOWER, TEST, VALEIG, WANTZ,
351  $ lquery
352  CHARACTER ORDER
353  INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
354  $ indisp, indiwo, indwrk, iscale, itmp1, j, jj,
355  $ llwork, lwmin, lhtrd, lwtrd, ib, indhous,
356  $ nsplit
357  DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
358  $ sigma, smlnum, tmp1, vll, vuu
359 * ..
360 * .. External Functions ..
361  LOGICAL LSAME
362  INTEGER ILAENV
363  DOUBLE PRECISION DLAMCH, DLANSB
364  EXTERNAL lsame, dlamch, dlansb, ilaenv
365 * ..
366 * .. External Subroutines ..
367  EXTERNAL dcopy, dgemv, dlacpy, dlascl, dscal,
369  $ dsytrd_sb2st
370 * ..
371 * .. Intrinsic Functions ..
372  INTRINSIC max, min, sqrt
373 * ..
374 * .. Executable Statements ..
375 *
376 * Test the input parameters.
377 *
378  wantz = lsame( jobz, 'V' )
379  alleig = lsame( range, 'A' )
380  valeig = lsame( range, 'V' )
381  indeig = lsame( range, 'I' )
382  lower = lsame( uplo, 'L' )
383  lquery = ( lwork.EQ.-1 )
384 *
385  info = 0
386  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
387  info = -1
388  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
389  info = -2
390  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
391  info = -3
392  ELSE IF( n.LT.0 ) THEN
393  info = -4
394  ELSE IF( kd.LT.0 ) THEN
395  info = -5
396  ELSE IF( ldab.LT.kd+1 ) THEN
397  info = -7
398  ELSE IF( wantz .AND. ldq.LT.max( 1, n ) ) THEN
399  info = -9
400  ELSE
401  IF( valeig ) THEN
402  IF( n.GT.0 .AND. vu.LE.vl )
403  $ info = -11
404  ELSE IF( indeig ) THEN
405  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
406  info = -12
407  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
408  info = -13
409  END IF
410  END IF
411  END IF
412  IF( info.EQ.0 ) THEN
413  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
414  $ info = -18
415  END IF
416 *
417  IF( info.EQ.0 ) THEN
418  IF( n.LE.1 ) THEN
419  lwmin = 1
420  work( 1 ) = lwmin
421  ELSE
422  ib = ilaenv( 18, 'DSYTRD_SB2ST', jobz, n, kd, -1, -1 )
423  lhtrd = ilaenv( 19, 'DSYTRD_SB2ST', jobz, n, kd, ib, -1 )
424  lwtrd = ilaenv( 20, 'DSYTRD_SB2ST', jobz, 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 dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:133
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
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:145
subroutine dstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEIN
Definition: dstein.f:176
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
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...