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.8.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 ILAENV2STAGE
340: DOUBLE PRECISION DLAMCH, DLANSY
341: EXTERNAL LSAME, DLAMCH, DLANSY, ILAENV2STAGE
342: * ..
343: * .. External Subroutines ..
344: EXTERNAL DCOPY, DLACPY, DORGTR, DORMTR, DSCAL, DSTEBZ,
345: $ DSTEIN, DSTEQR, DSTERF, DSWAP, XERBLA,
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 = ILAENV2STAGE( 1, 'DSYTRD_2STAGE', JOBZ,
397: $ N, -1, -1, -1 )
398: IB = ILAENV2STAGE( 2, 'DSYTRD_2STAGE', JOBZ,
399: $ N, KD, -1, -1 )
400: LHTRD = ILAENV2STAGE( 3, 'DSYTRD_2STAGE', JOBZ,
401: $ N, KD, IB, -1 )
402: LWTRD = ILAENV2STAGE( 4, 'DSYTRD_2STAGE', JOBZ,
403: $ N, KD, IB, -1 )
404: LWMIN = MAX( 8*N, 3*N + LHTRD + LWTRD )
405: WORK( 1 ) = LWMIN
406: END IF
407: *
408: IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY )
409: $ INFO = -17
410: END IF
411: *
412: IF( INFO.NE.0 ) THEN
413: CALL XERBLA( 'DSYEVX_2STAGE', -INFO )
414: RETURN
415: ELSE IF( LQUERY ) THEN
416: RETURN
417: END IF
418: *
419: * Quick return if possible
420: *
421: M = 0
422: IF( N.EQ.0 ) THEN
423: RETURN
424: END IF
425: *
426: IF( N.EQ.1 ) THEN
427: IF( ALLEIG .OR. INDEIG ) THEN
428: M = 1
429: W( 1 ) = A( 1, 1 )
430: ELSE
431: IF( VL.LT.A( 1, 1 ) .AND. VU.GE.A( 1, 1 ) ) THEN
432: M = 1
433: W( 1 ) = A( 1, 1 )
434: END IF
435: END IF
436: IF( WANTZ )
437: $ Z( 1, 1 ) = ONE
438: RETURN
439: END IF
440: *
441: * Get machine constants.
442: *
443: SAFMIN = DLAMCH( 'Safe minimum' )
444: EPS = DLAMCH( 'Precision' )
445: SMLNUM = SAFMIN / EPS
446: BIGNUM = ONE / SMLNUM
447: RMIN = SQRT( SMLNUM )
448: RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) )
449: *
450: * Scale matrix to allowable range, if necessary.
451: *
452: ISCALE = 0
453: ABSTLL = ABSTOL
454: IF( VALEIG ) THEN
455: VLL = VL
456: VUU = VU
457: END IF
458: ANRM = DLANSY( 'M', UPLO, N, A, LDA, WORK )
459: IF( ANRM.GT.ZERO .AND. ANRM.LT.RMIN ) THEN
460: ISCALE = 1
461: SIGMA = RMIN / ANRM
462: ELSE IF( ANRM.GT.RMAX ) THEN
463: ISCALE = 1
464: SIGMA = RMAX / ANRM
465: END IF
466: IF( ISCALE.EQ.1 ) THEN
467: IF( LOWER ) THEN
468: DO 10 J = 1, N
469: CALL DSCAL( N-J+1, SIGMA, A( J, J ), 1 )
470: 10 CONTINUE
471: ELSE
472: DO 20 J = 1, N
473: CALL DSCAL( J, SIGMA, A( 1, J ), 1 )
474: 20 CONTINUE
475: END IF
476: IF( ABSTOL.GT.0 )
477: $ ABSTLL = ABSTOL*SIGMA
478: IF( VALEIG ) THEN
479: VLL = VL*SIGMA
480: VUU = VU*SIGMA
481: END IF
482: END IF
483: *
484: * Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
485: *
486: INDTAU = 1
487: INDE = INDTAU + N
488: INDD = INDE + N
489: INDHOUS = INDD + N
490: INDWRK = INDHOUS + LHTRD
491: LLWORK = LWORK - INDWRK + 1
492: *
493: CALL DSYTRD_2STAGE( JOBZ, UPLO, N, A, LDA, WORK( INDD ),
494: $ WORK( INDE ), WORK( INDTAU ), WORK( INDHOUS ),
495: $ LHTRD, WORK( INDWRK ), LLWORK, IINFO )
496: *
497: * If all eigenvalues are desired and ABSTOL is less than or equal to
498: * zero, then call DSTERF or DORGTR and SSTEQR. If this fails for
499: * some eigenvalue, then try DSTEBZ.
500: *
501: TEST = .FALSE.
502: IF( INDEIG ) THEN
503: IF( IL.EQ.1 .AND. IU.EQ.N ) THEN
504: TEST = .TRUE.
505: END IF
506: END IF
507: IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN
508: CALL DCOPY( N, WORK( INDD ), 1, W, 1 )
509: INDEE = INDWRK + 2*N
510: IF( .NOT.WANTZ ) THEN
511: CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
512: CALL DSTERF( N, W, WORK( INDEE ), INFO )
513: ELSE
514: CALL DLACPY( 'A', N, N, A, LDA, Z, LDZ )
515: CALL DORGTR( UPLO, N, Z, LDZ, WORK( INDTAU ),
516: $ WORK( INDWRK ), LLWORK, IINFO )
517: CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
518: CALL DSTEQR( JOBZ, N, W, WORK( INDEE ), Z, LDZ,
519: $ WORK( INDWRK ), INFO )
520: IF( INFO.EQ.0 ) THEN
521: DO 30 I = 1, N
522: IFAIL( I ) = 0
523: 30 CONTINUE
524: END IF
525: END IF
526: IF( INFO.EQ.0 ) THEN
527: M = N
528: GO TO 40
529: END IF
530: INFO = 0
531: END IF
532: *
533: * Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
534: *
535: IF( WANTZ ) THEN
536: ORDER = 'B'
537: ELSE
538: ORDER = 'E'
539: END IF
540: INDIBL = 1
541: INDISP = INDIBL + N
542: INDIWO = INDISP + N
543: CALL DSTEBZ( RANGE, ORDER, N, VLL, VUU, IL, IU, ABSTLL,
544: $ WORK( INDD ), WORK( INDE ), M, NSPLIT, W,
545: $ IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWRK ),
546: $ IWORK( INDIWO ), INFO )
547: *
548: IF( WANTZ ) THEN
549: CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W,
550: $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
551: $ WORK( INDWRK ), IWORK( INDIWO ), IFAIL, INFO )
552: *
553: * Apply orthogonal matrix used in reduction to tridiagonal
554: * form to eigenvectors returned by DSTEIN.
555: *
556: INDWKN = INDE
557: LLWRKN = LWORK - INDWKN + 1
558: CALL DORMTR( 'L', UPLO, 'N', N, M, A, LDA, WORK( INDTAU ), Z,
559: $ LDZ, WORK( INDWKN ), LLWRKN, IINFO )
560: END IF
561: *
562: * If matrix was scaled, then rescale eigenvalues appropriately.
563: *
564: 40 CONTINUE
565: IF( ISCALE.EQ.1 ) THEN
566: IF( INFO.EQ.0 ) THEN
567: IMAX = M
568: ELSE
569: IMAX = INFO - 1
570: END IF
571: CALL DSCAL( IMAX, ONE / SIGMA, W, 1 )
572: END IF
573: *
574: * If eigenvalues are not in order, then sort them, along with
575: * eigenvectors.
576: *
577: IF( WANTZ ) THEN
578: DO 60 J = 1, M - 1
579: I = 0
580: TMP1 = W( J )
581: DO 50 JJ = J + 1, M
582: IF( W( JJ ).LT.TMP1 ) THEN
583: I = JJ
584: TMP1 = W( JJ )
585: END IF
586: 50 CONTINUE
587: *
588: IF( I.NE.0 ) THEN
589: ITMP1 = IWORK( INDIBL+I-1 )
590: W( I ) = W( J )
591: IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
592: W( J ) = TMP1
593: IWORK( INDIBL+J-1 ) = ITMP1
594: CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
595: IF( INFO.NE.0 ) THEN
596: ITMP1 = IFAIL( I )
597: IFAIL( I ) = IFAIL( J )
598: IFAIL( J ) = ITMP1
599: END IF
600: END IF
601: 60 CONTINUE
602: END IF
603: *
604: * Set WORK(1) to optimal workspace size.
605: *
606: WORK( 1 ) = LWMIN
607: *
608: RETURN
609: *
610: * End of DSYEVX_2STAGE
611: *
612: END
CVSweb interface <joel.bertrand@systella.fr>