LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsyevr_2stage.f
Go to the documentation of this file.
1*> \brief <b> DSYEVR_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*> Download DSYEVR_2STAGE + dependencies
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyevr_2stage.f">
12*> [TGZ]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyevr_2stage.f">
14*> [ZIP]</a>
15*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyevr_2stage.f">
16*> [TXT]</a>
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DSYEVR_2STAGE( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
22* IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK,
23* LWORK, IWORK, LIWORK, INFO )
24*
25* IMPLICIT NONE
26*
27* .. Scalar Arguments ..
28* CHARACTER JOBZ, RANGE, UPLO
29* INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
30* DOUBLE PRECISION ABSTOL, VL, VU
31* ..
32* .. Array Arguments ..
33* INTEGER ISUPPZ( * ), IWORK( * )
34* DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> DSYEVR_2STAGE computes selected eigenvalues and, optionally, eigenvectors
44*> of a real symmetric matrix A using the 2stage technique for
45*> the reduction to tridiagonal. Eigenvalues and eigenvectors can be
46*> selected by specifying either a range of values or a range of
47*> indices for the desired eigenvalues.
48*>
49*> DSYEVR_2STAGE first reduces the matrix A to tridiagonal form T with a call
50*> to DSYTRD. Then, whenever possible, DSYEVR_2STAGE calls DSTEMR to compute
51*> the eigenspectrum using Relatively Robust Representations. DSTEMR
52*> computes eigenvalues by the dqds algorithm, while orthogonal
53*> eigenvectors are computed from various "good" L D L^T representations
54*> (also known as Relatively Robust Representations). Gram-Schmidt
55*> orthogonalization is avoided as far as possible. More specifically,
56*> the various steps of the algorithm are as 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*> The desired accuracy of the output can be specified by the input
77*> parameter ABSTOL.
78*>
79*> For more details, see DSTEMR's documentation and:
80*> - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
81*> to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
82*> Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
83*> - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
84*> Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
85*> 2004. Also LAPACK Working Note 154.
86*> - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
87*> tridiagonal eigenvalue/eigenvector problem",
88*> Computer Science Division Technical Report No. UCB/CSD-97-971,
89*> UC Berkeley, May 1997.
90*>
91*>
92*> Note 1 : DSYEVR_2STAGE calls DSTEMR when the full spectrum is requested
93*> on machines which conform to the ieee-754 floating point standard.
94*> DSYEVR_2STAGE calls DSTEBZ and SSTEIN on non-ieee machines and
95*> when partial spectrum requests are made.
96*>
97*> Normal execution of DSTEMR may create NaNs and infinities and
98*> hence may abort due to a floating point exception in environments
99*> which do not handle NaNs and infinities in the ieee standard default
100*> manner.
101*> \endverbatim
102*
103* Arguments:
104* ==========
105*
106*> \param[in] JOBZ
107*> \verbatim
108*> JOBZ is CHARACTER*1
109*> = 'N': Compute eigenvalues only;
110*> = 'V': Compute eigenvalues and eigenvectors.
111*> Not available in this release.
112*> \endverbatim
113*>
114*> \param[in] RANGE
115*> \verbatim
116*> RANGE is CHARACTER*1
117*> = 'A': all eigenvalues will be found.
118*> = 'V': all eigenvalues in the half-open interval (VL,VU]
119*> will be found.
120*> = 'I': the IL-th through IU-th eigenvalues will be found.
121*> For RANGE = 'V' or 'I' and IU - IL < N - 1, DSTEBZ and
122*> DSTEIN are called
123*> \endverbatim
124*>
125*> \param[in] UPLO
126*> \verbatim
127*> UPLO is CHARACTER*1
128*> = 'U': Upper triangle of A is stored;
129*> = 'L': Lower triangle of A is stored.
130*> \endverbatim
131*>
132*> \param[in] N
133*> \verbatim
134*> N is INTEGER
135*> The order of the matrix A. N >= 0.
136*> \endverbatim
137*>
138*> \param[in,out] A
139*> \verbatim
140*> A is DOUBLE PRECISION array, dimension (LDA, N)
141*> On entry, the symmetric matrix A. If UPLO = 'U', the
142*> leading N-by-N upper triangular part of A contains the
143*> upper triangular part of the matrix A. If UPLO = 'L',
144*> the leading N-by-N lower triangular part of A contains
145*> the lower triangular part of the matrix A.
146*> On exit, the lower triangle (if UPLO='L') or the upper
147*> triangle (if UPLO='U') of A, including the diagonal, is
148*> destroyed.
149*> \endverbatim
150*>
151*> \param[in] LDA
152*> \verbatim
153*> LDA is INTEGER
154*> The leading dimension of the array A. LDA >= max(1,N).
155*> \endverbatim
156*>
157*> \param[in] VL
158*> \verbatim
159*> VL is DOUBLE PRECISION
160*> If RANGE='V', the lower bound of the interval to
161*> be searched for eigenvalues. VL < VU.
162*> Not referenced if RANGE = 'A' or 'I'.
163*> \endverbatim
164*>
165*> \param[in] VU
166*> \verbatim
167*> VU is DOUBLE PRECISION
168*> If RANGE='V', the upper bound of the interval to
169*> be searched for eigenvalues. VL < VU.
170*> Not referenced if RANGE = 'A' or 'I'.
171*> \endverbatim
172*>
173*> \param[in] IL
174*> \verbatim
175*> IL is INTEGER
176*> If RANGE='I', the index of the
177*> smallest eigenvalue to be returned.
178*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
179*> Not referenced if RANGE = 'A' or 'V'.
180*> \endverbatim
181*>
182*> \param[in] IU
183*> \verbatim
184*> IU is INTEGER
185*> If RANGE='I', the index of the
186*> largest eigenvalue to be returned.
187*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
188*> Not referenced if RANGE = 'A' or 'V'.
189*> \endverbatim
190*>
191*> \param[in] ABSTOL
192*> \verbatim
193*> ABSTOL is DOUBLE PRECISION
194*> The absolute error tolerance for the eigenvalues.
195*> An approximate eigenvalue is accepted as converged
196*> when it is determined to lie in an interval [a,b]
197*> of width less than or equal to
198*>
199*> ABSTOL + EPS * max( |a|,|b| ) ,
200*>
201*> where EPS is the machine precision. If ABSTOL is less than
202*> or equal to zero, then EPS*|T| will be used in its place,
203*> where |T| is the 1-norm of the tridiagonal matrix obtained
204*> by reducing A to tridiagonal form.
205*>
206*> See "Computing Small Singular Values of Bidiagonal Matrices
207*> with Guaranteed High Relative Accuracy," by Demmel and
208*> Kahan, LAPACK Working Note #3.
209*>
210*> If high relative accuracy is important, set ABSTOL to
211*> DLAMCH( 'Safe minimum' ). Doing so will guarantee that
212*> eigenvalues are computed to high relative accuracy when
213*> possible in future releases. The current code does not
214*> make any guarantees about high relative accuracy, but
215*> future releases will. See J. Barlow and J. Demmel,
216*> "Computing Accurate Eigensystems of Scaled Diagonally
217*> Dominant Matrices", LAPACK Working Note #7, for a discussion
218*> of which matrices define their eigenvalues to high relative
219*> accuracy.
220*> \endverbatim
221*>
222*> \param[out] M
223*> \verbatim
224*> M is INTEGER
225*> The total number of eigenvalues found. 0 <= M <= N.
226*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
227*> \endverbatim
228*>
229*> \param[out] W
230*> \verbatim
231*> W is DOUBLE PRECISION array, dimension (N)
232*> The first M elements contain the selected eigenvalues in
233*> ascending order.
234*> \endverbatim
235*>
236*> \param[out] Z
237*> \verbatim
238*> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
239*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
240*> contain the orthonormal eigenvectors of the matrix A
241*> corresponding to the selected eigenvalues, with the i-th
242*> column of Z holding the eigenvector associated with W(i).
243*> If JOBZ = 'N', then Z is not referenced.
244*> Note: the user must ensure that at least max(1,M) columns are
245*> supplied in the array Z; if RANGE = 'V', the exact value of M
246*> is not known in advance and an upper bound must be used.
247*> Supplying N columns is always safe.
248*> \endverbatim
249*>
250*> \param[in] LDZ
251*> \verbatim
252*> LDZ is INTEGER
253*> The leading dimension of the array Z. LDZ >= 1, and if
254*> JOBZ = 'V', LDZ >= max(1,N).
255*> \endverbatim
256*>
257*> \param[out] ISUPPZ
258*> \verbatim
259*> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
260*> The support of the eigenvectors in Z, i.e., the indices
261*> indicating the nonzero elements in Z. The i-th eigenvector
262*> is nonzero only in elements ISUPPZ( 2*i-1 ) through
263*> ISUPPZ( 2*i ). This is an output of DSTEMR (tridiagonal
264*> matrix). The support of the eigenvectors of A is typically
265*> 1:N because of the orthogonal transformations applied by DORMTR.
266*> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
267*> \endverbatim
268*>
269*> \param[out] WORK
270*> \verbatim
271*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
272*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
273*> \endverbatim
274*>
275*> \param[in] LWORK
276*> \verbatim
277*> LWORK is INTEGER
278*> The dimension of the array WORK.
279*> If N <= 1, LWORK must be at least 1.
280*> If JOBZ = 'N' and N > 1, LWORK must be queried.
281*> LWORK = MAX(1, 26*N, dimension) where
282*> dimension = max(stage1,stage2) + (KD+1)*N + 5*N
283*> = N*KD + N*max(KD+1,FACTOPTNB)
284*> + max(2*KD*KD, KD*NTHREADS)
285*> + (KD+1)*N + 5*N
286*> where KD is the blocking size of the reduction,
287*> FACTOPTNB is the blocking used by the QR or LQ
288*> algorithm, usually FACTOPTNB=128 is a good choice
289*> NTHREADS is the number of threads used when
290*> openMP compilation is enabled, otherwise =1.
291*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
292*>
293*> If LWORK = -1, then a workspace query is assumed; the routine
294*> only calculates the optimal size of the WORK array, returns
295*> this value as the first entry of the WORK array, and no error
296*> message related to LWORK is issued by XERBLA.
297*> \endverbatim
298*>
299*> \param[out] IWORK
300*> \verbatim
301*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
302*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
303*> \endverbatim
304*>
305*> \param[in] LIWORK
306*> \verbatim
307*> LIWORK is INTEGER
308*> The dimension of the array IWORK.
309*> If N <= 1, LIWORK >= 1, else LIWORK >= 10*N.
310*>
311*> If LIWORK = -1, then a workspace query is assumed; the
312*> routine only calculates the optimal size of the IWORK array,
313*> returns this value as the first entry of the IWORK array, and
314*> no error message related to LIWORK is issued by XERBLA.
315*> \endverbatim
316*>
317*> \param[out] INFO
318*> \verbatim
319*> INFO is INTEGER
320*> = 0: successful exit
321*> < 0: if INFO = -i, the i-th argument had an illegal value
322*> > 0: Internal error
323*> \endverbatim
324*
325* Authors:
326* ========
327*
328*> \author Univ. of Tennessee
329*> \author Univ. of California Berkeley
330*> \author Univ. of Colorado Denver
331*> \author NAG Ltd.
332*
333*> \ingroup heevr_2stage
334*
335*> \par Contributors:
336* ==================
337*>
338*> Inderjit Dhillon, IBM Almaden, USA \n
339*> Osni Marques, LBNL/NERSC, USA \n
340*> Ken Stanley, Computer Science Division, University of
341*> California at Berkeley, USA \n
342*> Jason Riedy, Computer Science Division, University of
343*> California at Berkeley, USA \n
344*>
345*> \par Further Details:
346* =====================
347*>
348*> \verbatim
349*>
350*> All details about the 2stage techniques are available in:
351*>
352*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
353*> Parallel reduction to condensed forms for symmetric eigenvalue problems
354*> using aggregated fine-grained and memory-aware kernels. In Proceedings
355*> of 2011 International Conference for High Performance Computing,
356*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
357*> Article 8 , 11 pages.
358*> http://doi.acm.org/10.1145/2063384.2063394
359*>
360*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
361*> An improved parallel singular value algorithm and its implementation
362*> for multicore hardware, In Proceedings of 2013 International Conference
363*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
364*> Denver, Colorado, USA, 2013.
365*> Article 90, 12 pages.
366*> http://doi.acm.org/10.1145/2503210.2503292
367*>
368*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
369*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
370*> calculations based on fine-grained memory aware tasks.
371*> International Journal of High Performance Computing Applications.
372*> Volume 28 Issue 2, Pages 196-209, May 2014.
373*> http://hpc.sagepub.com/content/28/2/196
374*>
375*> \endverbatim
376*
377* =====================================================================
378 SUBROUTINE dsyevr_2stage( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
379 $ IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK,
380 $ LWORK, IWORK, LIWORK, INFO )
381*
382 IMPLICIT NONE
383*
384* -- LAPACK driver routine --
385* -- LAPACK is a software package provided by Univ. of Tennessee, --
386* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
387*
388* .. Scalar Arguments ..
389 CHARACTER JOBZ, RANGE, UPLO
390 INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
391 DOUBLE PRECISION ABSTOL, VL, VU
392* ..
393* .. Array Arguments ..
394 INTEGER ISUPPZ( * ), IWORK( * )
395 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
396* ..
397*
398* =====================================================================
399*
400* .. Parameters ..
401 DOUBLE PRECISION ZERO, ONE, TWO
402 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
403* ..
404* .. Local Scalars ..
405 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, VALEIG, WANTZ,
406 $ TRYRAC
407 CHARACTER ORDER
408 INTEGER I, IEEEOK, IINFO, IMAX, INDD, INDDD, INDE,
409 $ indee, indibl, indifl, indisp, indiwo, indtau,
410 $ indwk, indwkn, iscale, j, jj, liwmin,
411 $ llwork, llwrkn, lwmin, nsplit,
412 $ lhtrd, lwtrd, kd, ib, indhous
413 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
414 $ SIGMA, SMLNUM, TMP1, VLL, VUU
415* ..
416* .. External Functions ..
417 LOGICAL LSAME
418 INTEGER ILAENV, ILAENV2STAGE
419 DOUBLE PRECISION DLAMCH, DLANSY
420 EXTERNAL lsame, dlamch, dlansy, ilaenv,
421 $ ilaenv2stage
422* ..
423* .. External Subroutines ..
424 EXTERNAL dcopy, dormtr, dscal, dstebz, dstemr,
425 $ dstein,
427* ..
428* .. Intrinsic Functions ..
429 INTRINSIC max, min, sqrt
430* ..
431* .. Executable Statements ..
432*
433* Test the input parameters.
434*
435 ieeeok = ilaenv( 10, 'DSYEVR', 'N', 1, 2, 3, 4 )
436*
437 lower = lsame( uplo, 'L' )
438 wantz = lsame( jobz, 'V' )
439 alleig = lsame( range, 'A' )
440 valeig = lsame( range, 'V' )
441 indeig = lsame( range, 'I' )
442*
443 lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
444*
445 kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz, n, -1, -1,
446 $ -1 )
447 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1,
448 $ -1 )
449 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib,
450 $ -1 )
451 lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz, n, kd, ib,
452 $ -1 )
453*
454 IF( n.LE.1 ) THEN
455 lwmin = 1
456 liwmin = 1
457 ELSE
458 lwmin = max( 26*n, 5*n + lhtrd + lwtrd )
459 liwmin = 10*n
460 END IF
461*
462 info = 0
463 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
464 info = -1
465 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
466 info = -2
467 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
468 info = -3
469 ELSE IF( n.LT.0 ) THEN
470 info = -4
471 ELSE IF( lda.LT.max( 1, n ) ) THEN
472 info = -6
473 ELSE
474 IF( valeig ) THEN
475 IF( n.GT.0 .AND. vu.LE.vl )
476 $ info = -8
477 ELSE IF( indeig ) THEN
478 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
479 info = -9
480 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
481 info = -10
482 END IF
483 END IF
484 END IF
485 IF( info.EQ.0 ) THEN
486 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
487 info = -15
488 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
489 info = -18
490 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
491 info = -20
492 END IF
493 END IF
494*
495 IF( info.EQ.0 ) THEN
496* NB = ILAENV( 1, 'DSYTRD', UPLO, N, -1, -1, -1 )
497* NB = MAX( NB, ILAENV( 1, 'DORMTR', UPLO, N, -1, -1, -1 ) )
498* LWKOPT = MAX( ( NB+1 )*N, LWMIN )
499 work( 1 ) = lwmin
500 iwork( 1 ) = liwmin
501 END IF
502*
503 IF( info.NE.0 ) THEN
504 CALL xerbla( 'DSYEVR_2STAGE', -info )
505 RETURN
506 ELSE IF( lquery ) THEN
507 RETURN
508 END IF
509*
510* Quick return if possible
511*
512 m = 0
513 IF( n.EQ.0 ) THEN
514 work( 1 ) = 1
515 RETURN
516 END IF
517*
518 IF( n.EQ.1 ) THEN
519 work( 1 ) = 1
520 IF( alleig .OR. indeig ) THEN
521 m = 1
522 w( 1 ) = a( 1, 1 )
523 ELSE
524 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
525 m = 1
526 w( 1 ) = a( 1, 1 )
527 END IF
528 END IF
529 IF( wantz ) THEN
530 z( 1, 1 ) = one
531 isuppz( 1 ) = 1
532 isuppz( 2 ) = 1
533 END IF
534 RETURN
535 END IF
536*
537* Get machine constants.
538*
539 safmin = dlamch( 'Safe minimum' )
540 eps = dlamch( 'Precision' )
541 smlnum = safmin / eps
542 bignum = one / smlnum
543 rmin = sqrt( smlnum )
544 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
545*
546* Scale matrix to allowable range, if necessary.
547*
548 iscale = 0
549 abstll = abstol
550 IF (valeig) THEN
551 vll = vl
552 vuu = vu
553 END IF
554 anrm = dlansy( 'M', uplo, n, a, lda, work )
555 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
556 iscale = 1
557 sigma = rmin / anrm
558 ELSE IF( anrm.GT.rmax ) THEN
559 iscale = 1
560 sigma = rmax / anrm
561 END IF
562 IF( iscale.EQ.1 ) THEN
563 IF( lower ) THEN
564 DO 10 j = 1, n
565 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
566 10 CONTINUE
567 ELSE
568 DO 20 j = 1, n
569 CALL dscal( j, sigma, a( 1, j ), 1 )
570 20 CONTINUE
571 END IF
572 IF( abstol.GT.0 )
573 $ abstll = abstol*sigma
574 IF( valeig ) THEN
575 vll = vl*sigma
576 vuu = vu*sigma
577 END IF
578 END IF
579
580* Initialize indices into workspaces. Note: The IWORK indices are
581* used only if DSTERF or DSTEMR fail.
582
583* WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
584* elementary reflectors used in DSYTRD.
585 indtau = 1
586* WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
587 indd = indtau + n
588* WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
589* tridiagonal matrix from DSYTRD.
590 inde = indd + n
591* WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
592* -written by DSTEMR (the DSTERF path copies the diagonal to W).
593 inddd = inde + n
594* WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
595* -written while computing the eigenvalues in DSTERF and DSTEMR.
596 indee = inddd + n
597* INDHOUS is the starting offset Householder storage of stage 2
598 indhous = indee + n
599* INDWK is the starting offset of the left-over workspace, and
600* LLWORK is the remaining workspace size.
601 indwk = indhous + lhtrd
602 llwork = lwork - indwk + 1
603
604
605* IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
606* stores the block indices of each of the M<=N eigenvalues.
607 indibl = 1
608* IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
609* stores the starting and finishing indices of each block.
610 indisp = indibl + n
611* IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
612* that corresponding to eigenvectors that fail to converge in
613* DSTEIN. This information is discarded; if any fail, the driver
614* returns INFO > 0.
615 indifl = indisp + n
616* INDIWO is the offset of the remaining integer workspace.
617 indiwo = indifl + n
618
619*
620* Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
621*
622*
623 CALL dsytrd_2stage( jobz, uplo, n, a, lda, work( indd ),
624 $ work( inde ), work( indtau ), work( indhous ),
625 $ lhtrd, work( indwk ), llwork, iinfo )
626*
627* If all eigenvalues are desired
628* then call DSTERF or DSTEMR and DORMTR.
629*
630 IF( ( alleig .OR. ( indeig .AND. il.EQ.1 .AND. iu.EQ.n ) ) .AND.
631 $ ieeeok.EQ.1 ) THEN
632 IF( .NOT.wantz ) THEN
633 CALL dcopy( n, work( indd ), 1, w, 1 )
634 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
635 CALL dsterf( n, w, work( indee ), info )
636 ELSE
637 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
638 CALL dcopy( n, work( indd ), 1, work( inddd ), 1 )
639*
640 IF (abstol .LE. two*n*eps) THEN
641 tryrac = .true.
642 ELSE
643 tryrac = .false.
644 END IF
645 CALL dstemr( jobz, 'A', n, work( inddd ), work( indee ),
646 $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
647 $ tryrac, work( indwk ), lwork, iwork, liwork,
648 $ info )
649*
650*
651*
652* Apply orthogonal matrix used in reduction to tridiagonal
653* form to eigenvectors returned by DSTEMR.
654*
655 IF( wantz .AND. info.EQ.0 ) THEN
656 indwkn = inde
657 llwrkn = lwork - indwkn + 1
658 CALL dormtr( 'L', uplo, 'N', n, m, a, lda,
659 $ work( indtau ), z, ldz, work( indwkn ),
660 $ llwrkn, iinfo )
661 END IF
662 END IF
663*
664*
665 IF( info.EQ.0 ) THEN
666* Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
667* undefined.
668 m = n
669 GO TO 30
670 END IF
671 info = 0
672 END IF
673*
674* Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
675* Also call DSTEBZ and DSTEIN if DSTEMR fails.
676*
677 IF( wantz ) THEN
678 order = 'B'
679 ELSE
680 order = 'E'
681 END IF
682
683 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
684 $ work( indd ), work( inde ), m, nsplit, w,
685 $ iwork( indibl ), iwork( indisp ), work( indwk ),
686 $ iwork( indiwo ), info )
687*
688 IF( wantz ) THEN
689 CALL dstein( n, work( indd ), work( inde ), m, w,
690 $ iwork( indibl ), iwork( indisp ), z, ldz,
691 $ work( indwk ), iwork( indiwo ), iwork( indifl ),
692 $ info )
693*
694* Apply orthogonal matrix used in reduction to tridiagonal
695* form to eigenvectors returned by DSTEIN.
696*
697 indwkn = inde
698 llwrkn = lwork - indwkn + 1
699 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
700 $ z,
701 $ ldz, work( indwkn ), llwrkn, iinfo )
702 END IF
703*
704* If matrix was scaled, then rescale eigenvalues appropriately.
705*
706* Jump here if DSTEMR/DSTEIN succeeded.
707 30 CONTINUE
708 IF( iscale.EQ.1 ) THEN
709 IF( info.EQ.0 ) THEN
710 imax = m
711 ELSE
712 imax = info - 1
713 END IF
714 CALL dscal( imax, one / sigma, w, 1 )
715 END IF
716*
717* If eigenvalues are not in order, then sort them, along with
718* eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
719* It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
720* not return this detailed information to the user.
721*
722 IF( wantz ) THEN
723 DO 50 j = 1, m - 1
724 i = 0
725 tmp1 = w( j )
726 DO 40 jj = j + 1, m
727 IF( w( jj ).LT.tmp1 ) THEN
728 i = jj
729 tmp1 = w( jj )
730 END IF
731 40 CONTINUE
732*
733 IF( i.NE.0 ) THEN
734 w( i ) = w( j )
735 w( j ) = tmp1
736 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
737 END IF
738 50 CONTINUE
739 END IF
740*
741* Set WORK(1) to optimal workspace size.
742*
743 work( 1 ) = lwmin
744 iwork( 1 ) = liwmin
745*
746 RETURN
747*
748* End of DSYEVR_2STAGE
749*
750 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dsyevr_2stage(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, iwork, liwork, info)
DSYEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...
subroutine dsytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
DSYTRD_2STAGE
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
DSTEBZ
Definition dstebz.f:272
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:172
subroutine dstemr(jobz, range, n, d, e, vl, vu, il, iu, m, w, z, ldz, nzc, isuppz, tryrac, work, lwork, iwork, liwork, info)
DSTEMR
Definition dstemr.f:320
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
DORMTR
Definition dormtr.f:170