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