1: *> \brief \b DSTEMR
2: *
3: * =========== DOCUMENTATION ===========
4: *
5: * Online html documentation available at
6: * http://www.netlib.org/lapack/explore-html/
7: *
8: *> \htmlonly
9: *> Download DSTEMR + dependencies
10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dstemr.f">
11: *> [TGZ]</a>
12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dstemr.f">
13: *> [ZIP]</a>
14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dstemr.f">
15: *> [TXT]</a>
16: *> \endhtmlonly
17: *
18: * Definition:
19: * ===========
20: *
21: * SUBROUTINE DSTEMR( JOBZ, RANGE, N, D, E, VL, VU, IL, IU,
22: * M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK,
23: * IWORK, LIWORK, INFO )
24: *
25: * .. Scalar Arguments ..
26: * CHARACTER JOBZ, RANGE
27: * LOGICAL TRYRAC
28: * INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
29: * DOUBLE PRECISION VL, VU
30: * ..
31: * .. Array Arguments ..
32: * INTEGER ISUPPZ( * ), IWORK( * )
33: * DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * )
34: * DOUBLE PRECISION Z( LDZ, * )
35: * ..
36: *
37: *
38: *> \par Purpose:
39: * =============
40: *>
41: *> \verbatim
42: *>
43: *> DSTEMR computes selected eigenvalues and, optionally, eigenvectors
44: *> of a real symmetric tridiagonal matrix T. Any such unreduced matrix has
45: *> a well defined set of pairwise different real eigenvalues, the corresponding
46: *> real eigenvectors are pairwise orthogonal.
47: *>
48: *> The spectrum may be computed either completely or partially by specifying
49: *> either an interval (VL,VU] or a range of indices IL:IU for the desired
50: *> eigenvalues.
51: *>
52: *> Depending on the number of desired eigenvalues, these are computed either
53: *> by bisection or the dqds algorithm. Numerically orthogonal eigenvectors are
54: *> computed by the use of various suitable L D L^T factorizations near clusters
55: *> of close eigenvalues (referred to as RRRs, Relatively Robust
56: *> Representations). An informal sketch of the algorithm follows.
57: *>
58: *> For each unreduced block (submatrix) of T,
59: *> (a) Compute T - sigma I = L D L^T, so that L and D
60: *> define all the wanted eigenvalues to high relative accuracy.
61: *> This means that small relative changes in the entries of D and L
62: *> cause only small relative changes in the eigenvalues and
63: *> eigenvectors. The standard (unfactored) representation of the
64: *> tridiagonal matrix T does not have this property in general.
65: *> (b) Compute the eigenvalues to suitable accuracy.
66: *> If the eigenvectors are desired, the algorithm attains full
67: *> accuracy of the computed eigenvalues only right before
68: *> the corresponding vectors have to be computed, see steps c) and d).
69: *> (c) For each cluster of close eigenvalues, select a new
70: *> shift close to the cluster, find a new factorization, and refine
71: *> the shifted eigenvalues to suitable accuracy.
72: *> (d) For each eigenvalue with a large enough relative separation compute
73: *> the corresponding eigenvector by forming a rank revealing twisted
74: *> factorization. Go back to (c) for any clusters that remain.
75: *>
76: *> For more details, see:
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: *> Further Details
89: *> 1.DSTEMR works only on machines which follow IEEE-754
90: *> floating-point standard in their handling of infinities and NaNs.
91: *> This permits the use of efficient inner loops avoiding a check for
92: *> zero divisors.
93: *> \endverbatim
94: *
95: * Arguments:
96: * ==========
97: *
98: *> \param[in] JOBZ
99: *> \verbatim
100: *> JOBZ is CHARACTER*1
101: *> = 'N': Compute eigenvalues only;
102: *> = 'V': Compute eigenvalues and eigenvectors.
103: *> \endverbatim
104: *>
105: *> \param[in] RANGE
106: *> \verbatim
107: *> RANGE is CHARACTER*1
108: *> = 'A': all eigenvalues will be found.
109: *> = 'V': all eigenvalues in the half-open interval (VL,VU]
110: *> will be found.
111: *> = 'I': the IL-th through IU-th eigenvalues will be found.
112: *> \endverbatim
113: *>
114: *> \param[in] N
115: *> \verbatim
116: *> N is INTEGER
117: *> The order of the matrix. N >= 0.
118: *> \endverbatim
119: *>
120: *> \param[in,out] D
121: *> \verbatim
122: *> D is DOUBLE PRECISION array, dimension (N)
123: *> On entry, the N diagonal elements of the tridiagonal matrix
124: *> T. On exit, D is overwritten.
125: *> \endverbatim
126: *>
127: *> \param[in,out] E
128: *> \verbatim
129: *> E is DOUBLE PRECISION array, dimension (N)
130: *> On entry, the (N-1) subdiagonal elements of the tridiagonal
131: *> matrix T in elements 1 to N-1 of E. E(N) need not be set on
132: *> input, but is used internally as workspace.
133: *> On exit, E is overwritten.
134: *> \endverbatim
135: *>
136: *> \param[in] VL
137: *> \verbatim
138: *> VL is DOUBLE PRECISION
139: *> \endverbatim
140: *>
141: *> \param[in] VU
142: *> \verbatim
143: *> VU is DOUBLE PRECISION
144: *>
145: *> If RANGE='V', the lower and upper bounds of the interval to
146: *> be searched for eigenvalues. VL < VU.
147: *> Not referenced if RANGE = 'A' or 'I'.
148: *> \endverbatim
149: *>
150: *> \param[in] IL
151: *> \verbatim
152: *> IL is INTEGER
153: *> \endverbatim
154: *>
155: *> \param[in] IU
156: *> \verbatim
157: *> IU is INTEGER
158: *>
159: *> If RANGE='I', the indices (in ascending order) of the
160: *> smallest and largest eigenvalues to be returned.
161: *> 1 <= IL <= IU <= N, if N > 0.
162: *> Not referenced if RANGE = 'A' or 'V'.
163: *> \endverbatim
164: *>
165: *> \param[out] M
166: *> \verbatim
167: *> M is INTEGER
168: *> The total number of eigenvalues found. 0 <= M <= N.
169: *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
170: *> \endverbatim
171: *>
172: *> \param[out] W
173: *> \verbatim
174: *> W is DOUBLE PRECISION array, dimension (N)
175: *> The first M elements contain the selected eigenvalues in
176: *> ascending order.
177: *> \endverbatim
178: *>
179: *> \param[out] Z
180: *> \verbatim
181: *> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M) )
182: *> If JOBZ = 'V', and if INFO = 0, then the first M columns of Z
183: *> contain the orthonormal eigenvectors of the matrix T
184: *> corresponding to the selected eigenvalues, with the i-th
185: *> column of Z holding the eigenvector associated with W(i).
186: *> If JOBZ = 'N', then Z is not referenced.
187: *> Note: the user must ensure that at least max(1,M) columns are
188: *> supplied in the array Z; if RANGE = 'V', the exact value of M
189: *> is not known in advance and can be computed with a workspace
190: *> query by setting NZC = -1, see below.
191: *> \endverbatim
192: *>
193: *> \param[in] LDZ
194: *> \verbatim
195: *> LDZ is INTEGER
196: *> The leading dimension of the array Z. LDZ >= 1, and if
197: *> JOBZ = 'V', then LDZ >= max(1,N).
198: *> \endverbatim
199: *>
200: *> \param[in] NZC
201: *> \verbatim
202: *> NZC is INTEGER
203: *> The number of eigenvectors to be held in the array Z.
204: *> If RANGE = 'A', then NZC >= max(1,N).
205: *> If RANGE = 'V', then NZC >= the number of eigenvalues in (VL,VU].
206: *> If RANGE = 'I', then NZC >= IU-IL+1.
207: *> If NZC = -1, then a workspace query is assumed; the
208: *> routine calculates the number of columns of the array Z that
209: *> are needed to hold the eigenvectors.
210: *> This value is returned as the first entry of the Z array, and
211: *> no error message related to NZC is issued by XERBLA.
212: *> \endverbatim
213: *>
214: *> \param[out] ISUPPZ
215: *> \verbatim
216: *> ISUPPZ is INTEGER ARRAY, dimension ( 2*max(1,M) )
217: *> The support of the eigenvectors in Z, i.e., the indices
218: *> indicating the nonzero elements in Z. The i-th computed eigenvector
219: *> is nonzero only in elements ISUPPZ( 2*i-1 ) through
220: *> ISUPPZ( 2*i ). This is relevant in the case when the matrix
221: *> is split. ISUPPZ is only accessed when JOBZ is 'V' and N > 0.
222: *> \endverbatim
223: *>
224: *> \param[in,out] TRYRAC
225: *> \verbatim
226: *> TRYRAC is LOGICAL
227: *> If TRYRAC.EQ..TRUE., indicates that the code should check whether
228: *> the tridiagonal matrix defines its eigenvalues to high relative
229: *> accuracy. If so, the code uses relative-accuracy preserving
230: *> algorithms that might be (a bit) slower depending on the matrix.
231: *> If the matrix does not define its eigenvalues to high relative
232: *> accuracy, the code can uses possibly faster algorithms.
233: *> If TRYRAC.EQ..FALSE., the code is not required to guarantee
234: *> relatively accurate eigenvalues and can use the fastest possible
235: *> techniques.
236: *> On exit, a .TRUE. TRYRAC will be set to .FALSE. if the matrix
237: *> does not define its eigenvalues to high relative accuracy.
238: *> \endverbatim
239: *>
240: *> \param[out] WORK
241: *> \verbatim
242: *> WORK is DOUBLE PRECISION array, dimension (LWORK)
243: *> On exit, if INFO = 0, WORK(1) returns the optimal
244: *> (and minimal) LWORK.
245: *> \endverbatim
246: *>
247: *> \param[in] LWORK
248: *> \verbatim
249: *> LWORK is INTEGER
250: *> The dimension of the array WORK. LWORK >= max(1,18*N)
251: *> if JOBZ = 'V', and LWORK >= max(1,12*N) if JOBZ = 'N'.
252: *> If LWORK = -1, then a workspace query is assumed; the routine
253: *> only calculates the optimal size of the WORK array, returns
254: *> this value as the first entry of the WORK array, and no error
255: *> message related to LWORK is issued by XERBLA.
256: *> \endverbatim
257: *>
258: *> \param[out] IWORK
259: *> \verbatim
260: *> IWORK is INTEGER array, dimension (LIWORK)
261: *> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
262: *> \endverbatim
263: *>
264: *> \param[in] LIWORK
265: *> \verbatim
266: *> LIWORK is INTEGER
267: *> The dimension of the array IWORK. LIWORK >= max(1,10*N)
268: *> if the eigenvectors are desired, and LIWORK >= max(1,8*N)
269: *> if only the eigenvalues are to be computed.
270: *> If LIWORK = -1, then a workspace query is assumed; the
271: *> routine only calculates the optimal size of the IWORK array,
272: *> returns this value as the first entry of the IWORK array, and
273: *> no error message related to LIWORK is issued by XERBLA.
274: *> \endverbatim
275: *>
276: *> \param[out] INFO
277: *> \verbatim
278: *> INFO is INTEGER
279: *> On exit, INFO
280: *> = 0: successful exit
281: *> < 0: if INFO = -i, the i-th argument had an illegal value
282: *> > 0: if INFO = 1X, internal error in DLARRE,
283: *> if INFO = 2X, internal error in DLARRV.
284: *> Here, the digit X = ABS( IINFO ) < 10, where IINFO is
285: *> the nonzero error code returned by DLARRE or
286: *> DLARRV, respectively.
287: *> \endverbatim
288: *
289: * Authors:
290: * ========
291: *
292: *> \author Univ. of Tennessee
293: *> \author Univ. of California Berkeley
294: *> \author Univ. of Colorado Denver
295: *> \author NAG Ltd.
296: *
297: *> \date November 2011
298: *
299: *> \ingroup doubleOTHERcomputational
300: *
301: *> \par Contributors:
302: * ==================
303: *>
304: *> Beresford Parlett, University of California, Berkeley, USA \n
305: *> Jim Demmel, University of California, Berkeley, USA \n
306: *> Inderjit Dhillon, University of Texas, Austin, USA \n
307: *> Osni Marques, LBNL/NERSC, USA \n
308: *> Christof Voemel, University of California, Berkeley, USA
309: *
310: * =====================================================================
311: SUBROUTINE DSTEMR( JOBZ, RANGE, N, D, E, VL, VU, IL, IU,
312: $ M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK,
313: $ IWORK, LIWORK, INFO )
314: *
315: * -- LAPACK computational routine (version 3.4.0) --
316: * -- LAPACK is a software package provided by Univ. of Tennessee, --
317: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
318: * November 2011
319: *
320: * .. Scalar Arguments ..
321: CHARACTER JOBZ, RANGE
322: LOGICAL TRYRAC
323: INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
324: DOUBLE PRECISION VL, VU
325: * ..
326: * .. Array Arguments ..
327: INTEGER ISUPPZ( * ), IWORK( * )
328: DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * )
329: DOUBLE PRECISION Z( LDZ, * )
330: * ..
331: *
332: * =====================================================================
333: *
334: * .. Parameters ..
335: DOUBLE PRECISION ZERO, ONE, FOUR, MINRGP
336: PARAMETER ( ZERO = 0.0D0, ONE = 1.0D0,
337: $ FOUR = 4.0D0,
338: $ MINRGP = 1.0D-3 )
339: * ..
340: * .. Local Scalars ..
341: LOGICAL ALLEIG, INDEIG, LQUERY, VALEIG, WANTZ, ZQUERY
342: INTEGER I, IBEGIN, IEND, IFIRST, IIL, IINDBL, IINDW,
343: $ IINDWK, IINFO, IINSPL, IIU, ILAST, IN, INDD,
344: $ INDE2, INDERR, INDGP, INDGRS, INDWRK, ITMP,
345: $ ITMP2, J, JBLK, JJ, LIWMIN, LWMIN, NSPLIT,
346: $ NZCMIN, OFFSET, WBEGIN, WEND
347: DOUBLE PRECISION BIGNUM, CS, EPS, PIVMIN, R1, R2, RMAX, RMIN,
348: $ RTOL1, RTOL2, SAFMIN, SCALE, SMLNUM, SN,
349: $ THRESH, TMP, TNRM, WL, WU
350: * ..
351: * ..
352: * .. External Functions ..
353: LOGICAL LSAME
354: DOUBLE PRECISION DLAMCH, DLANST
355: EXTERNAL LSAME, DLAMCH, DLANST
356: * ..
357: * .. External Subroutines ..
358: EXTERNAL DCOPY, DLAE2, DLAEV2, DLARRC, DLARRE, DLARRJ,
359: $ DLARRR, DLARRV, DLASRT, DSCAL, DSWAP, XERBLA
360: * ..
361: * .. Intrinsic Functions ..
362: INTRINSIC MAX, MIN, SQRT
363:
364:
365: * ..
366: * .. Executable Statements ..
367: *
368: * Test the input parameters.
369: *
370: WANTZ = LSAME( JOBZ, 'V' )
371: ALLEIG = LSAME( RANGE, 'A' )
372: VALEIG = LSAME( RANGE, 'V' )
373: INDEIG = LSAME( RANGE, 'I' )
374: *
375: LQUERY = ( ( LWORK.EQ.-1 ).OR.( LIWORK.EQ.-1 ) )
376: ZQUERY = ( NZC.EQ.-1 )
377:
378: * DSTEMR needs WORK of size 6*N, IWORK of size 3*N.
379: * In addition, DLARRE needs WORK of size 6*N, IWORK of size 5*N.
380: * Furthermore, DLARRV needs WORK of size 12*N, IWORK of size 7*N.
381: IF( WANTZ ) THEN
382: LWMIN = 18*N
383: LIWMIN = 10*N
384: ELSE
385: * need less workspace if only the eigenvalues are wanted
386: LWMIN = 12*N
387: LIWMIN = 8*N
388: ENDIF
389:
390: WL = ZERO
391: WU = ZERO
392: IIL = 0
393: IIU = 0
394:
395: IF( VALEIG ) THEN
396: * We do not reference VL, VU in the cases RANGE = 'I','A'
397: * The interval (WL, WU] contains all the wanted eigenvalues.
398: * It is either given by the user or computed in DLARRE.
399: WL = VL
400: WU = VU
401: ELSEIF( INDEIG ) THEN
402: * We do not reference IL, IU in the cases RANGE = 'V','A'
403: IIL = IL
404: IIU = IU
405: ENDIF
406: *
407: INFO = 0
408: IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
409: INFO = -1
410: ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
411: INFO = -2
412: ELSE IF( N.LT.0 ) THEN
413: INFO = -3
414: ELSE IF( VALEIG .AND. N.GT.0 .AND. WU.LE.WL ) THEN
415: INFO = -7
416: ELSE IF( INDEIG .AND. ( IIL.LT.1 .OR. IIL.GT.N ) ) THEN
417: INFO = -8
418: ELSE IF( INDEIG .AND. ( IIU.LT.IIL .OR. IIU.GT.N ) ) THEN
419: INFO = -9
420: ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
421: INFO = -13
422: ELSE IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
423: INFO = -17
424: ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN
425: INFO = -19
426: END IF
427: *
428: * Get machine constants.
429: *
430: SAFMIN = DLAMCH( 'Safe minimum' )
431: EPS = DLAMCH( 'Precision' )
432: SMLNUM = SAFMIN / EPS
433: BIGNUM = ONE / SMLNUM
434: RMIN = SQRT( SMLNUM )
435: RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) )
436: *
437: IF( INFO.EQ.0 ) THEN
438: WORK( 1 ) = LWMIN
439: IWORK( 1 ) = LIWMIN
440: *
441: IF( WANTZ .AND. ALLEIG ) THEN
442: NZCMIN = N
443: ELSE IF( WANTZ .AND. VALEIG ) THEN
444: CALL DLARRC( 'T', N, VL, VU, D, E, SAFMIN,
445: $ NZCMIN, ITMP, ITMP2, INFO )
446: ELSE IF( WANTZ .AND. INDEIG ) THEN
447: NZCMIN = IIU-IIL+1
448: ELSE
449: * WANTZ .EQ. FALSE.
450: NZCMIN = 0
451: ENDIF
452: IF( ZQUERY .AND. INFO.EQ.0 ) THEN
453: Z( 1,1 ) = NZCMIN
454: ELSE IF( NZC.LT.NZCMIN .AND. .NOT.ZQUERY ) THEN
455: INFO = -14
456: END IF
457: END IF
458:
459: IF( INFO.NE.0 ) THEN
460: *
461: CALL XERBLA( 'DSTEMR', -INFO )
462: *
463: RETURN
464: ELSE IF( LQUERY .OR. ZQUERY ) THEN
465: RETURN
466: END IF
467: *
468: * Handle N = 0, 1, and 2 cases immediately
469: *
470: M = 0
471: IF( N.EQ.0 )
472: $ RETURN
473: *
474: IF( N.EQ.1 ) THEN
475: IF( ALLEIG .OR. INDEIG ) THEN
476: M = 1
477: W( 1 ) = D( 1 )
478: ELSE
479: IF( WL.LT.D( 1 ) .AND. WU.GE.D( 1 ) ) THEN
480: M = 1
481: W( 1 ) = D( 1 )
482: END IF
483: END IF
484: IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN
485: Z( 1, 1 ) = ONE
486: ISUPPZ(1) = 1
487: ISUPPZ(2) = 1
488: END IF
489: RETURN
490: END IF
491: *
492: IF( N.EQ.2 ) THEN
493: IF( .NOT.WANTZ ) THEN
494: CALL DLAE2( D(1), E(1), D(2), R1, R2 )
495: ELSE IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN
496: CALL DLAEV2( D(1), E(1), D(2), R1, R2, CS, SN )
497: END IF
498: IF( ALLEIG.OR.
499: $ (VALEIG.AND.(R2.GT.WL).AND.
500: $ (R2.LE.WU)).OR.
501: $ (INDEIG.AND.(IIL.EQ.1)) ) THEN
502: M = M+1
503: W( M ) = R2
504: IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN
505: Z( 1, M ) = -SN
506: Z( 2, M ) = CS
507: * Note: At most one of SN and CS can be zero.
508: IF (SN.NE.ZERO) THEN
509: IF (CS.NE.ZERO) THEN
510: ISUPPZ(2*M-1) = 1
511: ISUPPZ(2*M) = 2
512: ELSE
513: ISUPPZ(2*M-1) = 1
514: ISUPPZ(2*M) = 1
515: END IF
516: ELSE
517: ISUPPZ(2*M-1) = 2
518: ISUPPZ(2*M) = 2
519: END IF
520: ENDIF
521: ENDIF
522: IF( ALLEIG.OR.
523: $ (VALEIG.AND.(R1.GT.WL).AND.
524: $ (R1.LE.WU)).OR.
525: $ (INDEIG.AND.(IIU.EQ.2)) ) THEN
526: M = M+1
527: W( M ) = R1
528: IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN
529: Z( 1, M ) = CS
530: Z( 2, M ) = SN
531: * Note: At most one of SN and CS can be zero.
532: IF (SN.NE.ZERO) THEN
533: IF (CS.NE.ZERO) THEN
534: ISUPPZ(2*M-1) = 1
535: ISUPPZ(2*M) = 2
536: ELSE
537: ISUPPZ(2*M-1) = 1
538: ISUPPZ(2*M) = 1
539: END IF
540: ELSE
541: ISUPPZ(2*M-1) = 2
542: ISUPPZ(2*M) = 2
543: END IF
544: ENDIF
545: ENDIF
546: RETURN
547: END IF
548:
549: * Continue with general N
550:
551: INDGRS = 1
552: INDERR = 2*N + 1
553: INDGP = 3*N + 1
554: INDD = 4*N + 1
555: INDE2 = 5*N + 1
556: INDWRK = 6*N + 1
557: *
558: IINSPL = 1
559: IINDBL = N + 1
560: IINDW = 2*N + 1
561: IINDWK = 3*N + 1
562: *
563: * Scale matrix to allowable range, if necessary.
564: * The allowable range is related to the PIVMIN parameter; see the
565: * comments in DLARRD. The preference for scaling small values
566: * up is heuristic; we expect users' matrices not to be close to the
567: * RMAX threshold.
568: *
569: SCALE = ONE
570: TNRM = DLANST( 'M', N, D, E )
571: IF( TNRM.GT.ZERO .AND. TNRM.LT.RMIN ) THEN
572: SCALE = RMIN / TNRM
573: ELSE IF( TNRM.GT.RMAX ) THEN
574: SCALE = RMAX / TNRM
575: END IF
576: IF( SCALE.NE.ONE ) THEN
577: CALL DSCAL( N, SCALE, D, 1 )
578: CALL DSCAL( N-1, SCALE, E, 1 )
579: TNRM = TNRM*SCALE
580: IF( VALEIG ) THEN
581: * If eigenvalues in interval have to be found,
582: * scale (WL, WU] accordingly
583: WL = WL*SCALE
584: WU = WU*SCALE
585: ENDIF
586: END IF
587: *
588: * Compute the desired eigenvalues of the tridiagonal after splitting
589: * into smaller subblocks if the corresponding off-diagonal elements
590: * are small
591: * THRESH is the splitting parameter for DLARRE
592: * A negative THRESH forces the old splitting criterion based on the
593: * size of the off-diagonal. A positive THRESH switches to splitting
594: * which preserves relative accuracy.
595: *
596: IF( TRYRAC ) THEN
597: * Test whether the matrix warrants the more expensive relative approach.
598: CALL DLARRR( N, D, E, IINFO )
599: ELSE
600: * The user does not care about relative accurately eigenvalues
601: IINFO = -1
602: ENDIF
603: * Set the splitting criterion
604: IF (IINFO.EQ.0) THEN
605: THRESH = EPS
606: ELSE
607: THRESH = -EPS
608: * relative accuracy is desired but T does not guarantee it
609: TRYRAC = .FALSE.
610: ENDIF
611: *
612: IF( TRYRAC ) THEN
613: * Copy original diagonal, needed to guarantee relative accuracy
614: CALL DCOPY(N,D,1,WORK(INDD),1)
615: ENDIF
616: * Store the squares of the offdiagonal values of T
617: DO 5 J = 1, N-1
618: WORK( INDE2+J-1 ) = E(J)**2
619: 5 CONTINUE
620:
621: * Set the tolerance parameters for bisection
622: IF( .NOT.WANTZ ) THEN
623: * DLARRE computes the eigenvalues to full precision.
624: RTOL1 = FOUR * EPS
625: RTOL2 = FOUR * EPS
626: ELSE
627: * DLARRE computes the eigenvalues to less than full precision.
628: * DLARRV will refine the eigenvalue approximations, and we can
629: * need less accurate initial bisection in DLARRE.
630: * Note: these settings do only affect the subset case and DLARRE
631: RTOL1 = SQRT(EPS)
632: RTOL2 = MAX( SQRT(EPS)*5.0D-3, FOUR * EPS )
633: ENDIF
634: CALL DLARRE( RANGE, N, WL, WU, IIL, IIU, D, E,
635: $ WORK(INDE2), RTOL1, RTOL2, THRESH, NSPLIT,
636: $ IWORK( IINSPL ), M, W, WORK( INDERR ),
637: $ WORK( INDGP ), IWORK( IINDBL ),
638: $ IWORK( IINDW ), WORK( INDGRS ), PIVMIN,
639: $ WORK( INDWRK ), IWORK( IINDWK ), IINFO )
640: IF( IINFO.NE.0 ) THEN
641: INFO = 10 + ABS( IINFO )
642: RETURN
643: END IF
644: * Note that if RANGE .NE. 'V', DLARRE computes bounds on the desired
645: * part of the spectrum. All desired eigenvalues are contained in
646: * (WL,WU]
647:
648:
649: IF( WANTZ ) THEN
650: *
651: * Compute the desired eigenvectors corresponding to the computed
652: * eigenvalues
653: *
654: CALL DLARRV( N, WL, WU, D, E,
655: $ PIVMIN, IWORK( IINSPL ), M,
656: $ 1, M, MINRGP, RTOL1, RTOL2,
657: $ W, WORK( INDERR ), WORK( INDGP ), IWORK( IINDBL ),
658: $ IWORK( IINDW ), WORK( INDGRS ), Z, LDZ,
659: $ ISUPPZ, WORK( INDWRK ), IWORK( IINDWK ), IINFO )
660: IF( IINFO.NE.0 ) THEN
661: INFO = 20 + ABS( IINFO )
662: RETURN
663: END IF
664: ELSE
665: * DLARRE computes eigenvalues of the (shifted) root representation
666: * DLARRV returns the eigenvalues of the unshifted matrix.
667: * However, if the eigenvectors are not desired by the user, we need
668: * to apply the corresponding shifts from DLARRE to obtain the
669: * eigenvalues of the original matrix.
670: DO 20 J = 1, M
671: ITMP = IWORK( IINDBL+J-1 )
672: W( J ) = W( J ) + E( IWORK( IINSPL+ITMP-1 ) )
673: 20 CONTINUE
674: END IF
675: *
676:
677: IF ( TRYRAC ) THEN
678: * Refine computed eigenvalues so that they are relatively accurate
679: * with respect to the original matrix T.
680: IBEGIN = 1
681: WBEGIN = 1
682: DO 39 JBLK = 1, IWORK( IINDBL+M-1 )
683: IEND = IWORK( IINSPL+JBLK-1 )
684: IN = IEND - IBEGIN + 1
685: WEND = WBEGIN - 1
686: * check if any eigenvalues have to be refined in this block
687: 36 CONTINUE
688: IF( WEND.LT.M ) THEN
689: IF( IWORK( IINDBL+WEND ).EQ.JBLK ) THEN
690: WEND = WEND + 1
691: GO TO 36
692: END IF
693: END IF
694: IF( WEND.LT.WBEGIN ) THEN
695: IBEGIN = IEND + 1
696: GO TO 39
697: END IF
698:
699: OFFSET = IWORK(IINDW+WBEGIN-1)-1
700: IFIRST = IWORK(IINDW+WBEGIN-1)
701: ILAST = IWORK(IINDW+WEND-1)
702: RTOL2 = FOUR * EPS
703: CALL DLARRJ( IN,
704: $ WORK(INDD+IBEGIN-1), WORK(INDE2+IBEGIN-1),
705: $ IFIRST, ILAST, RTOL2, OFFSET, W(WBEGIN),
706: $ WORK( INDERR+WBEGIN-1 ),
707: $ WORK( INDWRK ), IWORK( IINDWK ), PIVMIN,
708: $ TNRM, IINFO )
709: IBEGIN = IEND + 1
710: WBEGIN = WEND + 1
711: 39 CONTINUE
712: ENDIF
713: *
714: * If matrix was scaled, then rescale eigenvalues appropriately.
715: *
716: IF( SCALE.NE.ONE ) THEN
717: CALL DSCAL( M, ONE / SCALE, W, 1 )
718: END IF
719: *
720: * If eigenvalues are not in increasing order, then sort them,
721: * possibly along with eigenvectors.
722: *
723: IF( NSPLIT.GT.1 ) THEN
724: IF( .NOT. WANTZ ) THEN
725: CALL DLASRT( 'I', M, W, IINFO )
726: IF( IINFO.NE.0 ) THEN
727: INFO = 3
728: RETURN
729: END IF
730: ELSE
731: DO 60 J = 1, M - 1
732: I = 0
733: TMP = W( J )
734: DO 50 JJ = J + 1, M
735: IF( W( JJ ).LT.TMP ) THEN
736: I = JJ
737: TMP = W( JJ )
738: END IF
739: 50 CONTINUE
740: IF( I.NE.0 ) THEN
741: W( I ) = W( J )
742: W( J ) = TMP
743: IF( WANTZ ) THEN
744: CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
745: ITMP = ISUPPZ( 2*I-1 )
746: ISUPPZ( 2*I-1 ) = ISUPPZ( 2*J-1 )
747: ISUPPZ( 2*J-1 ) = ITMP
748: ITMP = ISUPPZ( 2*I )
749: ISUPPZ( 2*I ) = ISUPPZ( 2*J )
750: ISUPPZ( 2*J ) = ITMP
751: END IF
752: END IF
753: 60 CONTINUE
754: END IF
755: ENDIF
756: *
757: *
758: WORK( 1 ) = LWMIN
759: IWORK( 1 ) = LIWMIN
760: RETURN
761: *
762: * End of DSTEMR
763: *
764: END
CVSweb interface <joel.bertrand@systella.fr>