1: *> \brief <b> DSYEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices</b>
2: *
3: * =========== DOCUMENTATION ===========
4: *
5: * Online html documentation available at
6: * http://www.netlib.org/lapack/explore-html/
7: *
8: *> \htmlonly
9: *> Download DSYEVR + dependencies
10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyevr.f">
11: *> [TGZ]</a>
12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyevr.f">
13: *> [ZIP]</a>
14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyevr.f">
15: *> [TXT]</a>
16: *> \endhtmlonly
17: *
18: * Definition:
19: * ===========
20: *
21: * SUBROUTINE DSYEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
22: * ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK,
23: * IWORK, LIWORK, INFO )
24: *
25: * .. Scalar Arguments ..
26: * CHARACTER JOBZ, RANGE, UPLO
27: * INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
28: * DOUBLE PRECISION ABSTOL, VL, VU
29: * ..
30: * .. Array Arguments ..
31: * INTEGER ISUPPZ( * ), IWORK( * )
32: * DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
33: * ..
34: *
35: *
36: *> \par Purpose:
37: * =============
38: *>
39: *> \verbatim
40: *>
41: *> DSYEVR computes selected eigenvalues and, optionally, eigenvectors
42: *> of a real symmetric matrix A. Eigenvalues and eigenvectors can be
43: *> selected by specifying either a range of values or a range of
44: *> indices for the desired eigenvalues.
45: *>
46: *> DSYEVR first reduces the matrix A to tridiagonal form T with a call
47: *> to DSYTRD. Then, whenever possible, DSYEVR calls DSTEMR to compute
48: *> the eigenspectrum using Relatively Robust Representations. DSTEMR
49: *> computes eigenvalues by the dqds algorithm, while orthogonal
50: *> eigenvectors are computed from various "good" L D L^T representations
51: *> (also known as Relatively Robust Representations). Gram-Schmidt
52: *> orthogonalization is avoided as far as possible. More specifically,
53: *> the various steps of the algorithm are as follows.
54: *>
55: *> For each unreduced block (submatrix) of T,
56: *> (a) Compute T - sigma I = L D L^T, so that L and D
57: *> define all the wanted eigenvalues to high relative accuracy.
58: *> This means that small relative changes in the entries of D and L
59: *> cause only small relative changes in the eigenvalues and
60: *> eigenvectors. The standard (unfactored) representation of the
61: *> tridiagonal matrix T does not have this property in general.
62: *> (b) Compute the eigenvalues to suitable accuracy.
63: *> If the eigenvectors are desired, the algorithm attains full
64: *> accuracy of the computed eigenvalues only right before
65: *> the corresponding vectors have to be computed, see steps c) and d).
66: *> (c) For each cluster of close eigenvalues, select a new
67: *> shift close to the cluster, find a new factorization, and refine
68: *> the shifted eigenvalues to suitable accuracy.
69: *> (d) For each eigenvalue with a large enough relative separation compute
70: *> the corresponding eigenvector by forming a rank revealing twisted
71: *> factorization. Go back to (c) for any clusters that remain.
72: *>
73: *> The desired accuracy of the output can be specified by the input
74: *> parameter ABSTOL.
75: *>
76: *> For more details, see DSTEMR's documentation and:
77: *> - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
78: *> to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
79: *> Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
80: *> - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
81: *> Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
82: *> 2004. Also LAPACK Working Note 154.
83: *> - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
84: *> tridiagonal eigenvalue/eigenvector problem",
85: *> Computer Science Division Technical Report No. UCB/CSD-97-971,
86: *> UC Berkeley, May 1997.
87: *>
88: *>
89: *> Note 1 : DSYEVR calls DSTEMR when the full spectrum is requested
90: *> on machines which conform to the ieee-754 floating point standard.
91: *> DSYEVR calls DSTEBZ and SSTEIN on non-ieee machines and
92: *> when partial spectrum requests are made.
93: *>
94: *> Normal execution of DSTEMR may create NaNs and infinities and
95: *> hence may abort due to a floating point exception in environments
96: *> which do not handle NaNs and infinities in the ieee standard default
97: *> manner.
98: *> \endverbatim
99: *
100: * Arguments:
101: * ==========
102: *
103: *> \param[in] JOBZ
104: *> \verbatim
105: *> JOBZ is CHARACTER*1
106: *> = 'N': Compute eigenvalues only;
107: *> = 'V': Compute eigenvalues and eigenvectors.
108: *> \endverbatim
109: *>
110: *> \param[in] RANGE
111: *> \verbatim
112: *> RANGE is CHARACTER*1
113: *> = 'A': all eigenvalues will be found.
114: *> = 'V': all eigenvalues in the half-open interval (VL,VU]
115: *> will be found.
116: *> = 'I': the IL-th through IU-th eigenvalues will be found.
117: *> For RANGE = 'V' or 'I' and IU - IL < N - 1, DSTEBZ and
118: *> DSTEIN are called
119: *> \endverbatim
120: *>
121: *> \param[in] UPLO
122: *> \verbatim
123: *> UPLO is CHARACTER*1
124: *> = 'U': Upper triangle of A is stored;
125: *> = 'L': Lower triangle of A is stored.
126: *> \endverbatim
127: *>
128: *> \param[in] N
129: *> \verbatim
130: *> N is INTEGER
131: *> The order of the matrix A. N >= 0.
132: *> \endverbatim
133: *>
134: *> \param[in,out] A
135: *> \verbatim
136: *> A is DOUBLE PRECISION array, dimension (LDA, N)
137: *> On entry, the symmetric matrix A. If UPLO = 'U', the
138: *> leading N-by-N upper triangular part of A contains the
139: *> upper triangular part of the matrix A. If UPLO = 'L',
140: *> the leading N-by-N lower triangular part of A contains
141: *> the lower triangular part of the matrix A.
142: *> On exit, the lower triangle (if UPLO='L') or the upper
143: *> triangle (if UPLO='U') of A, including the diagonal, is
144: *> destroyed.
145: *> \endverbatim
146: *>
147: *> \param[in] LDA
148: *> \verbatim
149: *> LDA is INTEGER
150: *> The leading dimension of the array A. LDA >= max(1,N).
151: *> \endverbatim
152: *>
153: *> \param[in] VL
154: *> \verbatim
155: *> VL is DOUBLE PRECISION
156: *> \endverbatim
157: *>
158: *> \param[in] VU
159: *> \verbatim
160: *> VU is DOUBLE PRECISION
161: *> If RANGE='V', the lower and upper bounds of the interval to
162: *> be searched for eigenvalues. VL < VU.
163: *> Not referenced if RANGE = 'A' or 'I'.
164: *> \endverbatim
165: *>
166: *> \param[in] IL
167: *> \verbatim
168: *> IL is INTEGER
169: *> \endverbatim
170: *>
171: *> \param[in] IU
172: *> \verbatim
173: *> IU is INTEGER
174: *> If RANGE='I', the indices (in ascending order) of the
175: *> smallest and largest eigenvalues to be returned.
176: *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
177: *> Not referenced if RANGE = 'A' or 'V'.
178: *> \endverbatim
179: *>
180: *> \param[in] ABSTOL
181: *> \verbatim
182: *> ABSTOL is DOUBLE PRECISION
183: *> The absolute error tolerance for the eigenvalues.
184: *> An approximate eigenvalue is accepted as converged
185: *> when it is determined to lie in an interval [a,b]
186: *> of width less than or equal to
187: *>
188: *> ABSTOL + EPS * max( |a|,|b| ) ,
189: *>
190: *> where EPS is the machine precision. If ABSTOL is less than
191: *> or equal to zero, then EPS*|T| will be used in its place,
192: *> where |T| is the 1-norm of the tridiagonal matrix obtained
193: *> by reducing A to tridiagonal form.
194: *>
195: *> See "Computing Small Singular Values of Bidiagonal Matrices
196: *> with Guaranteed High Relative Accuracy," by Demmel and
197: *> Kahan, LAPACK Working Note #3.
198: *>
199: *> If high relative accuracy is important, set ABSTOL to
200: *> DLAMCH( 'Safe minimum' ). Doing so will guarantee that
201: *> eigenvalues are computed to high relative accuracy when
202: *> possible in future releases. The current code does not
203: *> make any guarantees about high relative accuracy, but
204: *> future releases will. See J. Barlow and J. Demmel,
205: *> "Computing Accurate Eigensystems of Scaled Diagonally
206: *> Dominant Matrices", LAPACK Working Note #7, for a discussion
207: *> of which matrices define their eigenvalues to high relative
208: *> accuracy.
209: *> \endverbatim
210: *>
211: *> \param[out] M
212: *> \verbatim
213: *> M is INTEGER
214: *> The total number of eigenvalues found. 0 <= M <= N.
215: *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
216: *> \endverbatim
217: *>
218: *> \param[out] W
219: *> \verbatim
220: *> W is DOUBLE PRECISION array, dimension (N)
221: *> The first M elements contain the selected eigenvalues in
222: *> ascending order.
223: *> \endverbatim
224: *>
225: *> \param[out] Z
226: *> \verbatim
227: *> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
228: *> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
229: *> contain the orthonormal eigenvectors of the matrix A
230: *> corresponding to the selected eigenvalues, with the i-th
231: *> column of Z holding the eigenvector associated with W(i).
232: *> If JOBZ = 'N', then Z is not referenced.
233: *> Note: the user must ensure that at least max(1,M) columns are
234: *> supplied in the array Z; if RANGE = 'V', the exact value of M
235: *> is not known in advance and an upper bound must be used.
236: *> Supplying N columns is always safe.
237: *> \endverbatim
238: *>
239: *> \param[in] LDZ
240: *> \verbatim
241: *> LDZ is INTEGER
242: *> The leading dimension of the array Z. LDZ >= 1, and if
243: *> JOBZ = 'V', LDZ >= max(1,N).
244: *> \endverbatim
245: *>
246: *> \param[out] ISUPPZ
247: *> \verbatim
248: *> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
249: *> The support of the eigenvectors in Z, i.e., the indices
250: *> indicating the nonzero elements in Z. The i-th eigenvector
251: *> is nonzero only in elements ISUPPZ( 2*i-1 ) through
252: *> ISUPPZ( 2*i ).
253: *> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
254: *> \endverbatim
255: *>
256: *> \param[out] WORK
257: *> \verbatim
258: *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
259: *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
260: *> \endverbatim
261: *>
262: *> \param[in] LWORK
263: *> \verbatim
264: *> LWORK is INTEGER
265: *> The dimension of the array WORK. LWORK >= max(1,26*N).
266: *> For optimal efficiency, LWORK >= (NB+6)*N,
267: *> where NB is the max of the blocksize for DSYTRD and DORMTR
268: *> returned by ILAENV.
269: *>
270: *> If LWORK = -1, then a workspace query is assumed; the routine
271: *> only calculates the optimal size of the WORK array, returns
272: *> this value as the first entry of the WORK array, and no error
273: *> message related to LWORK is issued by XERBLA.
274: *> \endverbatim
275: *>
276: *> \param[out] IWORK
277: *> \verbatim
278: *> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
279: *> On exit, if INFO = 0, IWORK(1) returns the optimal LWORK.
280: *> \endverbatim
281: *>
282: *> \param[in] LIWORK
283: *> \verbatim
284: *> LIWORK is INTEGER
285: *> The dimension of the array IWORK. LIWORK >= max(1,10*N).
286: *>
287: *> If LIWORK = -1, then a workspace query is assumed; the
288: *> routine only calculates the optimal size of the IWORK array,
289: *> returns this value as the first entry of the IWORK array, and
290: *> no error message related to LIWORK is issued by XERBLA.
291: *> \endverbatim
292: *>
293: *> \param[out] INFO
294: *> \verbatim
295: *> INFO is INTEGER
296: *> = 0: successful exit
297: *> < 0: if INFO = -i, the i-th argument had an illegal value
298: *> > 0: Internal error
299: *> \endverbatim
300: *
301: * Authors:
302: * ========
303: *
304: *> \author Univ. of Tennessee
305: *> \author Univ. of California Berkeley
306: *> \author Univ. of Colorado Denver
307: *> \author NAG Ltd.
308: *
309: *> \date November 2011
310: *
311: *> \ingroup doubleSYeigen
312: *
313: *> \par Contributors:
314: * ==================
315: *>
316: *> Inderjit Dhillon, IBM Almaden, USA \n
317: *> Osni Marques, LBNL/NERSC, USA \n
318: *> Ken Stanley, Computer Science Division, University of
319: *> California at Berkeley, USA \n
320: *> Jason Riedy, Computer Science Division, University of
321: *> California at Berkeley, USA \n
322: *>
323: * =====================================================================
324: SUBROUTINE DSYEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
325: $ ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK,
326: $ IWORK, LIWORK, INFO )
327: *
328: * -- LAPACK driver routine (version 3.4.0) --
329: * -- LAPACK is a software package provided by Univ. of Tennessee, --
330: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
331: * November 2011
332: *
333: * .. Scalar Arguments ..
334: CHARACTER JOBZ, RANGE, UPLO
335: INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
336: DOUBLE PRECISION ABSTOL, VL, VU
337: * ..
338: * .. Array Arguments ..
339: INTEGER ISUPPZ( * ), IWORK( * )
340: DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
341: * ..
342: *
343: * =====================================================================
344: *
345: * .. Parameters ..
346: DOUBLE PRECISION ZERO, ONE, TWO
347: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0, TWO = 2.0D+0 )
348: * ..
349: * .. Local Scalars ..
350: LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, VALEIG, WANTZ,
351: $ TRYRAC
352: CHARACTER ORDER
353: INTEGER I, IEEEOK, IINFO, IMAX, INDD, INDDD, INDE,
354: $ INDEE, INDIBL, INDIFL, INDISP, INDIWO, INDTAU,
355: $ INDWK, INDWKN, ISCALE, J, JJ, LIWMIN,
356: $ LLWORK, LLWRKN, LWKOPT, LWMIN, NB, 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, DLANSY
364: EXTERNAL LSAME, ILAENV, DLAMCH, DLANSY
365: * ..
366: * .. External Subroutines ..
367: EXTERNAL DCOPY, DORMTR, DSCAL, DSTEBZ, DSTEMR, DSTEIN,
368: $ DSTERF, DSWAP, DSYTRD, XERBLA
369: * ..
370: * .. Intrinsic Functions ..
371: INTRINSIC MAX, MIN, SQRT
372: * ..
373: * .. Executable Statements ..
374: *
375: * Test the input parameters.
376: *
377: IEEEOK = ILAENV( 10, 'DSYEVR', 'N', 1, 2, 3, 4 )
378: *
379: LOWER = LSAME( UPLO, 'L' )
380: WANTZ = LSAME( JOBZ, 'V' )
381: ALLEIG = LSAME( RANGE, 'A' )
382: VALEIG = LSAME( RANGE, 'V' )
383: INDEIG = LSAME( RANGE, 'I' )
384: *
385: LQUERY = ( ( LWORK.EQ.-1 ) .OR. ( LIWORK.EQ.-1 ) )
386: *
387: LWMIN = MAX( 1, 26*N )
388: LIWMIN = MAX( 1, 10*N )
389: *
390: INFO = 0
391: IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
392: INFO = -1
393: ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
394: INFO = -2
395: ELSE IF( .NOT.( LOWER .OR. LSAME( UPLO, 'U' ) ) ) THEN
396: INFO = -3
397: ELSE IF( N.LT.0 ) THEN
398: INFO = -4
399: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
400: INFO = -6
401: ELSE
402: IF( VALEIG ) THEN
403: IF( N.GT.0 .AND. VU.LE.VL )
404: $ INFO = -8
405: ELSE IF( INDEIG ) THEN
406: IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
407: INFO = -9
408: ELSE IF( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
409: INFO = -10
410: END IF
411: END IF
412: END IF
413: IF( INFO.EQ.0 ) THEN
414: IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
415: INFO = -15
416: ELSE IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
417: INFO = -18
418: ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN
419: INFO = -20
420: END IF
421: END IF
422: *
423: IF( INFO.EQ.0 ) THEN
424: NB = ILAENV( 1, 'DSYTRD', UPLO, N, -1, -1, -1 )
425: NB = MAX( NB, ILAENV( 1, 'DORMTR', UPLO, N, -1, -1, -1 ) )
426: LWKOPT = MAX( ( NB+1 )*N, LWMIN )
427: WORK( 1 ) = LWKOPT
428: IWORK( 1 ) = LIWMIN
429: END IF
430: *
431: IF( INFO.NE.0 ) THEN
432: CALL XERBLA( 'DSYEVR', -INFO )
433: RETURN
434: ELSE IF( LQUERY ) THEN
435: RETURN
436: END IF
437: *
438: * Quick return if possible
439: *
440: M = 0
441: IF( N.EQ.0 ) THEN
442: WORK( 1 ) = 1
443: RETURN
444: END IF
445: *
446: IF( N.EQ.1 ) THEN
447: WORK( 1 ) = 7
448: IF( ALLEIG .OR. INDEIG ) THEN
449: M = 1
450: W( 1 ) = A( 1, 1 )
451: ELSE
452: IF( VL.LT.A( 1, 1 ) .AND. VU.GE.A( 1, 1 ) ) THEN
453: M = 1
454: W( 1 ) = A( 1, 1 )
455: END IF
456: END IF
457: IF( WANTZ ) THEN
458: Z( 1, 1 ) = ONE
459: ISUPPZ( 1 ) = 1
460: ISUPPZ( 2 ) = 1
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: END IF
482: ANRM = DLANSY( 'M', UPLO, N, A, LDA, WORK )
483: IF( ANRM.GT.ZERO .AND. ANRM.LT.RMIN ) THEN
484: ISCALE = 1
485: SIGMA = RMIN / ANRM
486: ELSE IF( ANRM.GT.RMAX ) THEN
487: ISCALE = 1
488: SIGMA = RMAX / ANRM
489: END IF
490: IF( ISCALE.EQ.1 ) THEN
491: IF( LOWER ) THEN
492: DO 10 J = 1, N
493: CALL DSCAL( N-J+1, SIGMA, A( J, J ), 1 )
494: 10 CONTINUE
495: ELSE
496: DO 20 J = 1, N
497: CALL DSCAL( J, SIGMA, A( 1, J ), 1 )
498: 20 CONTINUE
499: END IF
500: IF( ABSTOL.GT.0 )
501: $ ABSTLL = ABSTOL*SIGMA
502: IF( VALEIG ) THEN
503: VLL = VL*SIGMA
504: VUU = VU*SIGMA
505: END IF
506: END IF
507:
508: * Initialize indices into workspaces. Note: The IWORK indices are
509: * used only if DSTERF or DSTEMR fail.
510:
511: * WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
512: * elementary reflectors used in DSYTRD.
513: INDTAU = 1
514: * WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
515: INDD = INDTAU + N
516: * WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
517: * tridiagonal matrix from DSYTRD.
518: INDE = INDD + N
519: * WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
520: * -written by DSTEMR (the DSTERF path copies the diagonal to W).
521: INDDD = INDE + N
522: * WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
523: * -written while computing the eigenvalues in DSTERF and DSTEMR.
524: INDEE = INDDD + N
525: * INDWK is the starting offset of the left-over workspace, and
526: * LLWORK is the remaining workspace size.
527: INDWK = INDEE + N
528: LLWORK = LWORK - INDWK + 1
529:
530: * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
531: * stores the block indices of each of the M<=N eigenvalues.
532: INDIBL = 1
533: * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
534: * stores the starting and finishing indices of each block.
535: INDISP = INDIBL + N
536: * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
537: * that corresponding to eigenvectors that fail to converge in
538: * DSTEIN. This information is discarded; if any fail, the driver
539: * returns INFO > 0.
540: INDIFL = INDISP + N
541: * INDIWO is the offset of the remaining integer workspace.
542: INDIWO = INDISP + N
543:
544: *
545: * Call DSYTRD to reduce symmetric matrix to tridiagonal form.
546: *
547: CALL DSYTRD( UPLO, N, A, LDA, WORK( INDD ), WORK( INDE ),
548: $ WORK( INDTAU ), WORK( INDWK ), LLWORK, IINFO )
549: *
550: * If all eigenvalues are desired
551: * then call DSTERF or DSTEMR and DORMTR.
552: *
553: IF( ( ALLEIG .OR. ( INDEIG .AND. IL.EQ.1 .AND. IU.EQ.N ) ) .AND.
554: $ IEEEOK.EQ.1 ) THEN
555: IF( .NOT.WANTZ ) THEN
556: CALL DCOPY( N, WORK( INDD ), 1, W, 1 )
557: CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
558: CALL DSTERF( N, W, WORK( INDEE ), INFO )
559: ELSE
560: CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
561: CALL DCOPY( N, WORK( INDD ), 1, WORK( INDDD ), 1 )
562: *
563: IF (ABSTOL .LE. TWO*N*EPS) THEN
564: TRYRAC = .TRUE.
565: ELSE
566: TRYRAC = .FALSE.
567: END IF
568: CALL DSTEMR( JOBZ, 'A', N, WORK( INDDD ), WORK( INDEE ),
569: $ VL, VU, IL, IU, M, W, Z, LDZ, N, ISUPPZ,
570: $ TRYRAC, WORK( INDWK ), LWORK, IWORK, LIWORK,
571: $ INFO )
572: *
573: *
574: *
575: * Apply orthogonal matrix used in reduction to tridiagonal
576: * form to eigenvectors returned by DSTEIN.
577: *
578: IF( WANTZ .AND. INFO.EQ.0 ) THEN
579: INDWKN = INDE
580: LLWRKN = LWORK - INDWKN + 1
581: CALL DORMTR( 'L', UPLO, 'N', N, M, A, LDA,
582: $ WORK( INDTAU ), Z, LDZ, WORK( INDWKN ),
583: $ LLWRKN, IINFO )
584: END IF
585: END IF
586: *
587: *
588: IF( INFO.EQ.0 ) THEN
589: * Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
590: * undefined.
591: M = N
592: GO TO 30
593: END IF
594: INFO = 0
595: END IF
596: *
597: * Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
598: * Also call DSTEBZ and DSTEIN if DSTEMR fails.
599: *
600: IF( WANTZ ) THEN
601: ORDER = 'B'
602: ELSE
603: ORDER = 'E'
604: END IF
605:
606: CALL DSTEBZ( RANGE, ORDER, N, VLL, VUU, IL, IU, ABSTLL,
607: $ WORK( INDD ), WORK( INDE ), M, NSPLIT, W,
608: $ IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWK ),
609: $ IWORK( INDIWO ), INFO )
610: *
611: IF( WANTZ ) THEN
612: CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W,
613: $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
614: $ WORK( INDWK ), IWORK( INDIWO ), IWORK( INDIFL ),
615: $ INFO )
616: *
617: * Apply orthogonal matrix used in reduction to tridiagonal
618: * form to eigenvectors returned by DSTEIN.
619: *
620: INDWKN = INDE
621: LLWRKN = LWORK - INDWKN + 1
622: CALL DORMTR( 'L', UPLO, 'N', N, M, A, LDA, WORK( INDTAU ), Z,
623: $ LDZ, WORK( INDWKN ), LLWRKN, IINFO )
624: END IF
625: *
626: * If matrix was scaled, then rescale eigenvalues appropriately.
627: *
628: * Jump here if DSTEMR/DSTEIN succeeded.
629: 30 CONTINUE
630: IF( ISCALE.EQ.1 ) THEN
631: IF( INFO.EQ.0 ) THEN
632: IMAX = M
633: ELSE
634: IMAX = INFO - 1
635: END IF
636: CALL DSCAL( IMAX, ONE / SIGMA, W, 1 )
637: END IF
638: *
639: * If eigenvalues are not in order, then sort them, along with
640: * eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
641: * It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
642: * not return this detailed information to the user.
643: *
644: IF( WANTZ ) THEN
645: DO 50 J = 1, M - 1
646: I = 0
647: TMP1 = W( J )
648: DO 40 JJ = J + 1, M
649: IF( W( JJ ).LT.TMP1 ) THEN
650: I = JJ
651: TMP1 = W( JJ )
652: END IF
653: 40 CONTINUE
654: *
655: IF( I.NE.0 ) THEN
656: W( I ) = W( J )
657: W( J ) = TMP1
658: CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
659: END IF
660: 50 CONTINUE
661: END IF
662: *
663: * Set WORK(1) to optimal workspace size.
664: *
665: WORK( 1 ) = LWKOPT
666: IWORK( 1 ) = LIWMIN
667: *
668: RETURN
669: *
670: * End of DSYEVR
671: *
672: END
CVSweb interface <joel.bertrand@systella.fr>