1*> \brief <b> DSTEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER 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 DSTEVX + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dstevx.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dstevx.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dstevx.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18*  Definition:
19*  ===========
20*
21*       SUBROUTINE DSTEVX( JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL,
22*                          M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO )
23*
24*       .. Scalar Arguments ..
25*       CHARACTER          JOBZ, RANGE
26*       INTEGER            IL, INFO, IU, LDZ, M, N
27*       DOUBLE PRECISION   ABSTOL, VL, VU
28*       ..
29*       .. Array Arguments ..
30*       INTEGER            IFAIL( * ), IWORK( * )
31*       DOUBLE PRECISION   D( * ), E( * ), W( * ), WORK( * ), Z( LDZ, * )
32*       ..
33*
34*
35*> \par Purpose:
36*  =============
37*>
38*> \verbatim
39*>
40*> DSTEVX computes selected eigenvalues and, optionally, eigenvectors
41*> of a real symmetric tridiagonal matrix A.  Eigenvalues and
42*> eigenvectors can be selected by specifying either a range of values
43*> or a range of indices for the desired eigenvalues.
44*> \endverbatim
45*
46*  Arguments:
47*  ==========
48*
49*> \param[in] JOBZ
50*> \verbatim
51*>          JOBZ is CHARACTER*1
52*>          = 'N':  Compute eigenvalues only;
53*>          = 'V':  Compute eigenvalues and eigenvectors.
54*> \endverbatim
55*>
56*> \param[in] RANGE
57*> \verbatim
58*>          RANGE is CHARACTER*1
59*>          = 'A': all eigenvalues will be found.
60*>          = 'V': all eigenvalues in the half-open interval (VL,VU]
61*>                 will be found.
62*>          = 'I': the IL-th through IU-th eigenvalues will be found.
63*> \endverbatim
64*>
65*> \param[in] N
66*> \verbatim
67*>          N is INTEGER
68*>          The order of the matrix.  N >= 0.
69*> \endverbatim
70*>
71*> \param[in,out] D
72*> \verbatim
73*>          D is DOUBLE PRECISION array, dimension (N)
74*>          On entry, the n diagonal elements of the tridiagonal matrix
75*>          A.
76*>          On exit, D may be multiplied by a constant factor chosen
77*>          to avoid over/underflow in computing the eigenvalues.
78*> \endverbatim
79*>
80*> \param[in,out] E
81*> \verbatim
82*>          E is DOUBLE PRECISION array, dimension (max(1,N-1))
83*>          On entry, the (n-1) subdiagonal elements of the tridiagonal
84*>          matrix A in elements 1 to N-1 of E.
85*>          On exit, E may be multiplied by a constant factor chosen
86*>          to avoid over/underflow in computing the eigenvalues.
87*> \endverbatim
88*>
89*> \param[in] VL
90*> \verbatim
91*>          VL is DOUBLE PRECISION
92*> \endverbatim
93*>
94*> \param[in] VU
95*> \verbatim
96*>          VU is DOUBLE PRECISION
97*>          If RANGE='V', the lower and upper bounds of the interval to
98*>          be searched for eigenvalues. VL < VU.
99*>          Not referenced if RANGE = 'A' or 'I'.
100*> \endverbatim
101*>
102*> \param[in] IL
103*> \verbatim
104*>          IL is INTEGER
105*> \endverbatim
106*>
107*> \param[in] IU
108*> \verbatim
109*>          IU is INTEGER
110*>          If RANGE='I', the indices (in ascending order) of the
111*>          smallest and largest eigenvalues to be returned.
112*>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
113*>          Not referenced if RANGE = 'A' or 'V'.
114*> \endverbatim
115*>
116*> \param[in] ABSTOL
117*> \verbatim
118*>          ABSTOL is DOUBLE PRECISION
119*>          The absolute error tolerance for the eigenvalues.
120*>          An approximate eigenvalue is accepted as converged
121*>          when it is determined to lie in an interval [a,b]
122*>          of width less than or equal to
123*>
124*>                  ABSTOL + EPS *   max( |a|,|b| ) ,
125*>
126*>          where EPS is the machine precision.  If ABSTOL is less
127*>          than or equal to zero, then  EPS*|T|  will be used in
128*>          its place, where |T| is the 1-norm of the tridiagonal
129*>          matrix.
130*>
131*>          Eigenvalues will be computed most accurately when ABSTOL is
132*>          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
133*>          If this routine returns with INFO>0, indicating that some
134*>          eigenvectors did not converge, try setting ABSTOL to
135*>          2*DLAMCH('S').
136*>
137*>          See "Computing Small Singular Values of Bidiagonal Matrices
138*>          with Guaranteed High Relative Accuracy," by Demmel and
139*>          Kahan, LAPACK Working Note #3.
140*> \endverbatim
141*>
142*> \param[out] M
143*> \verbatim
144*>          M is INTEGER
145*>          The total number of eigenvalues found.  0 <= M <= N.
146*>          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
147*> \endverbatim
148*>
149*> \param[out] W
150*> \verbatim
151*>          W is DOUBLE PRECISION array, dimension (N)
152*>          The first M elements contain the selected eigenvalues in
153*>          ascending order.
154*> \endverbatim
155*>
156*> \param[out] Z
157*> \verbatim
158*>          Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M) )
159*>          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
160*>          contain the orthonormal eigenvectors of the matrix A
161*>          corresponding to the selected eigenvalues, with the i-th
162*>          column of Z holding the eigenvector associated with W(i).
163*>          If an eigenvector fails to converge (INFO > 0), then that
164*>          column of Z contains the latest approximation to the
165*>          eigenvector, and the index of the eigenvector is returned
166*>          in IFAIL.  If JOBZ = 'N', then Z is not referenced.
167*>          Note: the user must ensure that at least max(1,M) columns are
168*>          supplied in the array Z; if RANGE = 'V', the exact value of M
169*>          is not known in advance and an upper bound must be used.
170*> \endverbatim
171*>
172*> \param[in] LDZ
173*> \verbatim
174*>          LDZ is INTEGER
175*>          The leading dimension of the array Z.  LDZ >= 1, and if
176*>          JOBZ = 'V', LDZ >= max(1,N).
177*> \endverbatim
178*>
179*> \param[out] WORK
180*> \verbatim
181*>          WORK is DOUBLE PRECISION array, dimension (5*N)
182*> \endverbatim
183*>
184*> \param[out] IWORK
185*> \verbatim
186*>          IWORK is INTEGER array, dimension (5*N)
187*> \endverbatim
188*>
189*> \param[out] IFAIL
190*> \verbatim
191*>          IFAIL is INTEGER array, dimension (N)
192*>          If JOBZ = 'V', then if INFO = 0, the first M elements of
193*>          IFAIL are zero.  If INFO > 0, then IFAIL contains the
194*>          indices of the eigenvectors that failed to converge.
195*>          If JOBZ = 'N', then IFAIL is not referenced.
196*> \endverbatim
197*>
198*> \param[out] INFO
199*> \verbatim
200*>          INFO is INTEGER
201*>          = 0:  successful exit
202*>          < 0:  if INFO = -i, the i-th argument had an illegal value
203*>          > 0:  if INFO = i, then i eigenvectors failed to converge.
204*>                Their indices are stored in array IFAIL.
205*> \endverbatim
206*
207*  Authors:
208*  ========
209*
210*> \author Univ. of Tennessee
211*> \author Univ. of California Berkeley
212*> \author Univ. of Colorado Denver
213*> \author NAG Ltd.
214*
215*> \date November 2011
216*
217*> \ingroup doubleOTHEReigen
218*
219*  =====================================================================
220      SUBROUTINE DSTEVX( JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL,
221     $                   M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO )
222*
223*  -- LAPACK driver routine (version 3.4.0) --
224*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
225*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
226*     November 2011
227*
228*     .. Scalar Arguments ..
229      CHARACTER          JOBZ, RANGE
230      INTEGER            IL, INFO, IU, LDZ, M, N
231      DOUBLE PRECISION   ABSTOL, VL, VU
232*     ..
233*     .. Array Arguments ..
234      INTEGER            IFAIL( * ), IWORK( * )
235      DOUBLE PRECISION   D( * ), E( * ), W( * ), WORK( * ), Z( LDZ, * )
236*     ..
237*
238*  =====================================================================
239*
240*     .. Parameters ..
241      DOUBLE PRECISION   ZERO, ONE
242      PARAMETER          ( ZERO = 0.0D0, ONE = 1.0D0 )
243*     ..
244*     .. Local Scalars ..
245      LOGICAL            ALLEIG, INDEIG, TEST, VALEIG, WANTZ
246      CHARACTER          ORDER
247      INTEGER            I, IMAX, INDIBL, INDISP, INDIWO, INDWRK,
248     $                   ISCALE, ITMP1, J, JJ, NSPLIT
249      DOUBLE PRECISION   BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA, SMLNUM,
250     $                   TMP1, TNRM, VLL, VUU
251*     ..
252*     .. External Functions ..
253      LOGICAL            LSAME
254      DOUBLE PRECISION   DLAMCH, DLANST
255      EXTERNAL           LSAME, DLAMCH, DLANST
256*     ..
257*     .. External Subroutines ..
258      EXTERNAL           DCOPY, DSCAL, DSTEBZ, DSTEIN, DSTEQR, DSTERF,
259     $                   DSWAP, XERBLA
260*     ..
261*     .. Intrinsic Functions ..
262      INTRINSIC          MAX, MIN, SQRT
263*     ..
264*     .. Executable Statements ..
265*
266*     Test the input parameters.
267*
268      WANTZ = LSAME( JOBZ, 'V' )
269      ALLEIG = LSAME( RANGE, 'A' )
270      VALEIG = LSAME( RANGE, 'V' )
271      INDEIG = LSAME( RANGE, 'I' )
272*
273      INFO = 0
274      IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
275         INFO = -1
276      ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
277         INFO = -2
278      ELSE IF( N.LT.0 ) THEN
279         INFO = -3
280      ELSE
281         IF( VALEIG ) THEN
282            IF( N.GT.0 .AND. VU.LE.VL )
283     $         INFO = -7
284         ELSE IF( INDEIG ) THEN
285            IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
286               INFO = -8
287            ELSE IF( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
288               INFO = -9
289            END IF
290         END IF
291      END IF
292      IF( INFO.EQ.0 ) THEN
293         IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) )
294     $      INFO = -14
295      END IF
296*
297      IF( INFO.NE.0 ) THEN
298         CALL XERBLA( 'DSTEVX', -INFO )
299         RETURN
300      END IF
301*
302*     Quick return if possible
303*
304      M = 0
305      IF( N.EQ.0 )
306     $   RETURN
307*
308      IF( N.EQ.1 ) THEN
309         IF( ALLEIG .OR. INDEIG ) THEN
310            M = 1
311            W( 1 ) = D( 1 )
312         ELSE
313            IF( VL.LT.D( 1 ) .AND. VU.GE.D( 1 ) ) THEN
314               M = 1
315               W( 1 ) = D( 1 )
316            END IF
317         END IF
318         IF( WANTZ )
319     $      Z( 1, 1 ) = ONE
320         RETURN
321      END IF
322*
323*     Get machine constants.
324*
325      SAFMIN = DLAMCH( 'Safe minimum' )
326      EPS = DLAMCH( 'Precision' )
327      SMLNUM = SAFMIN / EPS
328      BIGNUM = ONE / SMLNUM
329      RMIN = SQRT( SMLNUM )
330      RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) )
331*
332*     Scale matrix to allowable range, if necessary.
333*
334      ISCALE = 0
335      IF( VALEIG ) THEN
336         VLL = VL
337         VUU = VU
338      ELSE
339         VLL = ZERO
340         VUU = ZERO
341      END IF
342      TNRM = DLANST( 'M', N, D, E )
343      IF( TNRM.GT.ZERO .AND. TNRM.LT.RMIN ) THEN
344         ISCALE = 1
345         SIGMA = RMIN / TNRM
346      ELSE IF( TNRM.GT.RMAX ) THEN
347         ISCALE = 1
348         SIGMA = RMAX / TNRM
349      END IF
350      IF( ISCALE.EQ.1 ) THEN
351         CALL DSCAL( N, SIGMA, D, 1 )
352         CALL DSCAL( N-1, SIGMA, E( 1 ), 1 )
353         IF( VALEIG ) THEN
354            VLL = VL*SIGMA
355            VUU = VU*SIGMA
356         END IF
357      END IF
358*
359*     If all eigenvalues are desired and ABSTOL is less than zero, then
360*     call DSTERF or SSTEQR.  If this fails for some eigenvalue, then
361*     try DSTEBZ.
362*
363      TEST = .FALSE.
364      IF( INDEIG ) THEN
365         IF( IL.EQ.1 .AND. IU.EQ.N ) THEN
366            TEST = .TRUE.
367         END IF
368      END IF
369      IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN
370         CALL DCOPY( N, D, 1, W, 1 )
371         CALL DCOPY( N-1, E( 1 ), 1, WORK( 1 ), 1 )
372         INDWRK = N + 1
373         IF( .NOT.WANTZ ) THEN
374            CALL DSTERF( N, W, WORK, INFO )
375         ELSE
376            CALL DSTEQR( 'I', N, W, WORK, Z, LDZ, WORK( INDWRK ), INFO )
377            IF( INFO.EQ.0 ) THEN
378               DO 10 I = 1, N
379                  IFAIL( I ) = 0
380   10          CONTINUE
381            END IF
382         END IF
383         IF( INFO.EQ.0 ) THEN
384            M = N
385            GO TO 20
386         END IF
387         INFO = 0
388      END IF
389*
390*     Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
391*
392      IF( WANTZ ) THEN
393         ORDER = 'B'
394      ELSE
395         ORDER = 'E'
396      END IF
397      INDWRK = 1
398      INDIBL = 1
399      INDISP = INDIBL + N
400      INDIWO = INDISP + N
401      CALL DSTEBZ( RANGE, ORDER, N, VLL, VUU, IL, IU, ABSTOL, D, E, M,
402     $             NSPLIT, W, IWORK( INDIBL ), IWORK( INDISP ),
403     $             WORK( INDWRK ), IWORK( INDIWO ), INFO )
404*
405      IF( WANTZ ) THEN
406         CALL DSTEIN( N, D, E, M, W, IWORK( INDIBL ), IWORK( INDISP ),
407     $                Z, LDZ, WORK( INDWRK ), IWORK( INDIWO ), IFAIL,
408     $                INFO )
409      END IF
410*
411*     If matrix was scaled, then rescale eigenvalues appropriately.
412*
413   20 CONTINUE
414      IF( ISCALE.EQ.1 ) THEN
415         IF( INFO.EQ.0 ) THEN
416            IMAX = M
417         ELSE
418            IMAX = INFO - 1
419         END IF
420         CALL DSCAL( IMAX, ONE / SIGMA, W, 1 )
421      END IF
422*
423*     If eigenvalues are not in order, then sort them, along with
424*     eigenvectors.
425*
426      IF( WANTZ ) THEN
427         DO 40 J = 1, M - 1
428            I = 0
429            TMP1 = W( J )
430            DO 30 JJ = J + 1, M
431               IF( W( JJ ).LT.TMP1 ) THEN
432                  I = JJ
433                  TMP1 = W( JJ )
434               END IF
435   30       CONTINUE
436*
437            IF( I.NE.0 ) THEN
438               ITMP1 = IWORK( INDIBL+I-1 )
439               W( I ) = W( J )
440               IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
441               W( J ) = TMP1
442               IWORK( INDIBL+J-1 ) = ITMP1
443               CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
444               IF( INFO.NE.0 ) THEN
445                  ITMP1 = IFAIL( I )
446                  IFAIL( I ) = IFAIL( J )
447                  IFAIL( J ) = ITMP1
448               END IF
449            END IF
450   40    CONTINUE
451      END IF
452*
453      RETURN
454*
455*     End of DSTEVX
456*
457      END
458