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