LAPACK 3.12.0
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*> \htmlonly
11*> Download DSYEVR_2STAGE + dependencies
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyevr_2stage.f">
13*> [TGZ]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyevr_2stage.f">
15*> [ZIP]</a>
16*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyevr_2stage.f">
17*> [TXT]</a>
18*> \endhtmlonly
19*
20* Definition:
21* ===========
22*
23* SUBROUTINE DSYEVR_2STAGE( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
24* IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK,
25* LWORK, IWORK, LIWORK, INFO )
26*
27* IMPLICIT NONE
28*
29* .. Scalar Arguments ..
30* CHARACTER JOBZ, RANGE, UPLO
31* INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
32* DOUBLE PRECISION ABSTOL, VL, VU
33* ..
34* .. Array Arguments ..
35* INTEGER ISUPPZ( * ), IWORK( * )
36* DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
37* ..
38*
39*
40*> \par Purpose:
41* =============
42*>
43*> \verbatim
44*>
45*> DSYEVR_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
49*> indices for the desired eigenvalues.
50*>
51*> DSYEVR_2STAGE first reduces the matrix A to tridiagonal form T with a call
52*> to DSYTRD. Then, whenever possible, DSYEVR_2STAGE calls DSTEMR to compute
53*> the eigenspectrum using Relatively Robust Representations. DSTEMR
54*> computes eigenvalues by the dqds algorithm, while orthogonal
55*> eigenvectors are computed from various "good" L D L^T representations
56*> (also known as Relatively Robust Representations). Gram-Schmidt
57*> orthogonalization is avoided as far as possible. More specifically,
58*> the various steps of the algorithm are as follows.
59*>
60*> For each unreduced block (submatrix) of T,
61*> (a) Compute T - sigma I = L D L^T, so that L and D
62*> define all the wanted eigenvalues to high relative accuracy.
63*> This means that small relative changes in the entries of D and L
64*> cause only small relative changes in the eigenvalues and
65*> eigenvectors. The standard (unfactored) representation of the
66*> tridiagonal matrix T does not have this property in general.
67*> (b) Compute the eigenvalues to suitable accuracy.
68*> If the eigenvectors are desired, the algorithm attains full
69*> accuracy of the computed eigenvalues only right before
70*> the corresponding vectors have to be computed, see steps c) and d).
71*> (c) For each cluster of close eigenvalues, select a new
72*> shift close to the cluster, find a new factorization, and refine
73*> the shifted eigenvalues to suitable accuracy.
74*> (d) For each eigenvalue with a large enough relative separation compute
75*> the corresponding eigenvector by forming a rank revealing twisted
76*> factorization. Go back to (c) for any clusters that remain.
77*>
78*> The desired accuracy of the output can be specified by the input
79*> parameter ABSTOL.
80*>
81*> For more details, see DSTEMR's documentation and:
82*> - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
83*> to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
84*> Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
85*> - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
86*> Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
87*> 2004. Also LAPACK Working Note 154.
88*> - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
89*> tridiagonal eigenvalue/eigenvector problem",
90*> Computer Science Division Technical Report No. UCB/CSD-97-971,
91*> UC Berkeley, May 1997.
92*>
93*>
94*> Note 1 : DSYEVR_2STAGE calls DSTEMR when the full spectrum is requested
95*> on machines which conform to the ieee-754 floating point standard.
96*> DSYEVR_2STAGE calls DSTEBZ and SSTEIN on non-ieee machines and
97*> when partial spectrum requests are made.
98*>
99*> Normal execution of DSTEMR may create NaNs and infinities and
100*> hence may abort due to a floating point exception in environments
101*> which do not handle NaNs and infinities in the ieee standard default
102*> manner.
103*> \endverbatim
104*
105* Arguments:
106* ==========
107*
108*> \param[in] JOBZ
109*> \verbatim
110*> JOBZ is CHARACTER*1
111*> = 'N': Compute eigenvalues only;
112*> = 'V': Compute eigenvalues and eigenvectors.
113*> Not available in this release.
114*> \endverbatim
115*>
116*> \param[in] RANGE
117*> \verbatim
118*> RANGE is CHARACTER*1
119*> = 'A': all eigenvalues will be found.
120*> = 'V': all eigenvalues in the half-open interval (VL,VU]
121*> will be found.
122*> = 'I': the IL-th through IU-th eigenvalues will be found.
123*> For RANGE = 'V' or 'I' and IU - IL < N - 1, DSTEBZ and
124*> DSTEIN are called
125*> \endverbatim
126*>
127*> \param[in] UPLO
128*> \verbatim
129*> UPLO is CHARACTER*1
130*> = 'U': Upper triangle of A is stored;
131*> = 'L': Lower triangle of A is stored.
132*> \endverbatim
133*>
134*> \param[in] N
135*> \verbatim
136*> N is INTEGER
137*> The order of the matrix A. N >= 0.
138*> \endverbatim
139*>
140*> \param[in,out] A
141*> \verbatim
142*> A is DOUBLE PRECISION array, dimension (LDA, N)
143*> On entry, the symmetric matrix A. If UPLO = 'U', the
144*> leading N-by-N upper triangular part of A contains the
145*> upper triangular part of the matrix A. If UPLO = 'L',
146*> the leading N-by-N lower triangular part of A contains
147*> the lower triangular part of the matrix A.
148*> On exit, the lower triangle (if UPLO='L') or the upper
149*> triangle (if UPLO='U') of A, including the diagonal, is
150*> destroyed.
151*> \endverbatim
152*>
153*> \param[in] LDA
154*> \verbatim
155*> LDA is INTEGER
156*> The leading dimension of the array A. LDA >= max(1,N).
157*> \endverbatim
158*>
159*> \param[in] VL
160*> \verbatim
161*> VL is DOUBLE PRECISION
162*> If RANGE='V', the lower bound of the interval to
163*> be searched for eigenvalues. VL < VU.
164*> Not referenced if RANGE = 'A' or 'I'.
165*> \endverbatim
166*>
167*> \param[in] VU
168*> \verbatim
169*> VU is DOUBLE PRECISION
170*> If RANGE='V', the upper bound of the interval to
171*> be searched for eigenvalues. VL < VU.
172*> Not referenced if RANGE = 'A' or 'I'.
173*> \endverbatim
174*>
175*> \param[in] IL
176*> \verbatim
177*> IL is INTEGER
178*> If RANGE='I', the index of the
179*> smallest eigenvalue to be returned.
180*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
181*> Not referenced if RANGE = 'A' or 'V'.
182*> \endverbatim
183*>
184*> \param[in] IU
185*> \verbatim
186*> IU is INTEGER
187*> If RANGE='I', the index of the
188*> largest eigenvalue to be returned.
189*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
190*> Not referenced if RANGE = 'A' or 'V'.
191*> \endverbatim
192*>
193*> \param[in] ABSTOL
194*> \verbatim
195*> ABSTOL is DOUBLE PRECISION
196*> The absolute error tolerance for the eigenvalues.
197*> An approximate eigenvalue is accepted as converged
198*> when it is determined to lie in an interval [a,b]
199*> of width less than or equal to
200*>
201*> ABSTOL + EPS * max( |a|,|b| ) ,
202*>
203*> where EPS is the machine precision. If ABSTOL is less than
204*> or equal to zero, then EPS*|T| will be used in its place,
205*> where |T| is the 1-norm of the tridiagonal matrix obtained
206*> by reducing A to tridiagonal form.
207*>
208*> See "Computing Small Singular Values of Bidiagonal Matrices
209*> with Guaranteed High Relative Accuracy," by Demmel and
210*> Kahan, LAPACK Working Note #3.
211*>
212*> If high relative accuracy is important, set ABSTOL to
213*> DLAMCH( 'Safe minimum' ). Doing so will guarantee that
214*> eigenvalues are computed to high relative accuracy when
215*> possible in future releases. The current code does not
216*> make any guarantees about high relative accuracy, but
217*> future releases will. See J. Barlow and J. Demmel,
218*> "Computing Accurate Eigensystems of Scaled Diagonally
219*> Dominant Matrices", LAPACK Working Note #7, for a discussion
220*> of which matrices define their eigenvalues to high relative
221*> accuracy.
222*> \endverbatim
223*>
224*> \param[out] M
225*> \verbatim
226*> M is INTEGER
227*> The total number of eigenvalues found. 0 <= M <= N.
228*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
229*> \endverbatim
230*>
231*> \param[out] W
232*> \verbatim
233*> W is DOUBLE PRECISION array, dimension (N)
234*> The first M elements contain the selected eigenvalues in
235*> ascending order.
236*> \endverbatim
237*>
238*> \param[out] Z
239*> \verbatim
240*> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
241*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
242*> contain the orthonormal eigenvectors of the matrix A
243*> corresponding to the selected eigenvalues, with the i-th
244*> column of Z holding the eigenvector associated with W(i).
245*> If JOBZ = 'N', then Z is not referenced.
246*> Note: the user must ensure that at least max(1,M) columns are
247*> supplied in the array Z; if RANGE = 'V', the exact value of M
248*> is not known in advance and an upper bound must be used.
249*> Supplying N columns is always safe.
250*> \endverbatim
251*>
252*> \param[in] LDZ
253*> \verbatim
254*> LDZ is INTEGER
255*> The leading dimension of the array Z. LDZ >= 1, and if
256*> JOBZ = 'V', LDZ >= max(1,N).
257*> \endverbatim
258*>
259*> \param[out] ISUPPZ
260*> \verbatim
261*> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
262*> The support of the eigenvectors in Z, i.e., the indices
263*> indicating the nonzero elements in Z. The i-th eigenvector
264*> is nonzero only in elements ISUPPZ( 2*i-1 ) through
265*> ISUPPZ( 2*i ). This is an output of DSTEMR (tridiagonal
266*> matrix). The support of the eigenvectors of A is typically
267*> 1:N because of the orthogonal transformations applied by DORMTR.
268*> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
269*> \endverbatim
270*>
271*> \param[out] WORK
272*> \verbatim
273*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
274*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
275*> \endverbatim
276*>
277*> \param[in] LWORK
278*> \verbatim
279*> LWORK is INTEGER
280*> The dimension of the array WORK.
281*> If JOBZ = 'N' and N > 1, LWORK must be queried.
282*> LWORK = MAX(1, 26*N, dimension) where
283*> dimension = max(stage1,stage2) + (KD+1)*N + 5*N
284*> = N*KD + N*max(KD+1,FACTOPTNB)
285*> + max(2*KD*KD, KD*NTHREADS)
286*> + (KD+1)*N + 5*N
287*> where KD is the blocking size of the reduction,
288*> FACTOPTNB is the blocking used by the QR or LQ
289*> algorithm, usually FACTOPTNB=128 is a good choice
290*> NTHREADS is the number of threads used when
291*> openMP compilation is enabled, otherwise =1.
292*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
293*>
294*> If LWORK = -1, then a workspace query is assumed; the routine
295*> only calculates the optimal size of the WORK array, returns
296*> this value as the first entry of the WORK array, and no error
297*> message related to LWORK is issued by XERBLA.
298*> \endverbatim
299*>
300*> \param[out] IWORK
301*> \verbatim
302*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
303*> On exit, if INFO = 0, IWORK(1) returns the optimal LWORK.
304*> \endverbatim
305*>
306*> \param[in] LIWORK
307*> \verbatim
308*> LIWORK is INTEGER
309*> The dimension of the array IWORK. LIWORK >= max(1,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, ilaenv2stage
421* ..
422* .. External Subroutines ..
423 EXTERNAL dcopy, dormtr, dscal, dstebz, dstemr, dstein,
425* ..
426* .. Intrinsic Functions ..
427 INTRINSIC max, min, sqrt
428* ..
429* .. Executable Statements ..
430*
431* Test the input parameters.
432*
433 ieeeok = ilaenv( 10, 'DSYEVR', 'N', 1, 2, 3, 4 )
434*
435 lower = lsame( uplo, 'L' )
436 wantz = lsame( jobz, 'V' )
437 alleig = lsame( range, 'A' )
438 valeig = lsame( range, 'V' )
439 indeig = lsame( range, 'I' )
440*
441 lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
442*
443 kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz, n, -1, -1, -1 )
444 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1, -1 )
445 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
446 lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
447 lwmin = max( 26*n, 5*n + lhtrd + lwtrd )
448 liwmin = max( 1, 10*n )
449*
450 info = 0
451 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
452 info = -1
453 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
454 info = -2
455 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
456 info = -3
457 ELSE IF( n.LT.0 ) THEN
458 info = -4
459 ELSE IF( lda.LT.max( 1, n ) ) THEN
460 info = -6
461 ELSE
462 IF( valeig ) THEN
463 IF( n.GT.0 .AND. vu.LE.vl )
464 $ info = -8
465 ELSE IF( indeig ) THEN
466 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
467 info = -9
468 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
469 info = -10
470 END IF
471 END IF
472 END IF
473 IF( info.EQ.0 ) THEN
474 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
475 info = -15
476 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
477 info = -18
478 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
479 info = -20
480 END IF
481 END IF
482*
483 IF( info.EQ.0 ) THEN
484* NB = ILAENV( 1, 'DSYTRD', UPLO, N, -1, -1, -1 )
485* NB = MAX( NB, ILAENV( 1, 'DORMTR', UPLO, N, -1, -1, -1 ) )
486* LWKOPT = MAX( ( NB+1 )*N, LWMIN )
487 work( 1 ) = lwmin
488 iwork( 1 ) = liwmin
489 END IF
490*
491 IF( info.NE.0 ) THEN
492 CALL xerbla( 'DSYEVR_2STAGE', -info )
493 RETURN
494 ELSE IF( lquery ) THEN
495 RETURN
496 END IF
497*
498* Quick return if possible
499*
500 m = 0
501 IF( n.EQ.0 ) THEN
502 work( 1 ) = 1
503 RETURN
504 END IF
505*
506 IF( n.EQ.1 ) THEN
507 work( 1 ) = 7
508 IF( alleig .OR. indeig ) THEN
509 m = 1
510 w( 1 ) = a( 1, 1 )
511 ELSE
512 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
513 m = 1
514 w( 1 ) = a( 1, 1 )
515 END IF
516 END IF
517 IF( wantz ) THEN
518 z( 1, 1 ) = one
519 isuppz( 1 ) = 1
520 isuppz( 2 ) = 1
521 END IF
522 RETURN
523 END IF
524*
525* Get machine constants.
526*
527 safmin = dlamch( 'Safe minimum' )
528 eps = dlamch( 'Precision' )
529 smlnum = safmin / eps
530 bignum = one / smlnum
531 rmin = sqrt( smlnum )
532 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
533*
534* Scale matrix to allowable range, if necessary.
535*
536 iscale = 0
537 abstll = abstol
538 IF (valeig) THEN
539 vll = vl
540 vuu = vu
541 END IF
542 anrm = dlansy( 'M', uplo, n, a, lda, work )
543 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
544 iscale = 1
545 sigma = rmin / anrm
546 ELSE IF( anrm.GT.rmax ) THEN
547 iscale = 1
548 sigma = rmax / anrm
549 END IF
550 IF( iscale.EQ.1 ) THEN
551 IF( lower ) THEN
552 DO 10 j = 1, n
553 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
554 10 CONTINUE
555 ELSE
556 DO 20 j = 1, n
557 CALL dscal( j, sigma, a( 1, j ), 1 )
558 20 CONTINUE
559 END IF
560 IF( abstol.GT.0 )
561 $ abstll = abstol*sigma
562 IF( valeig ) THEN
563 vll = vl*sigma
564 vuu = vu*sigma
565 END IF
566 END IF
567
568* Initialize indices into workspaces. Note: The IWORK indices are
569* used only if DSTERF or DSTEMR fail.
570
571* WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
572* elementary reflectors used in DSYTRD.
573 indtau = 1
574* WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
575 indd = indtau + n
576* WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
577* tridiagonal matrix from DSYTRD.
578 inde = indd + n
579* WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
580* -written by DSTEMR (the DSTERF path copies the diagonal to W).
581 inddd = inde + n
582* WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
583* -written while computing the eigenvalues in DSTERF and DSTEMR.
584 indee = inddd + n
585* INDHOUS is the starting offset Householder storage of stage 2
586 indhous = indee + n
587* INDWK is the starting offset of the left-over workspace, and
588* LLWORK is the remaining workspace size.
589 indwk = indhous + lhtrd
590 llwork = lwork - indwk + 1
591
592
593* IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
594* stores the block indices of each of the M<=N eigenvalues.
595 indibl = 1
596* IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
597* stores the starting and finishing indices of each block.
598 indisp = indibl + n
599* IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
600* that corresponding to eigenvectors that fail to converge in
601* DSTEIN. This information is discarded; if any fail, the driver
602* returns INFO > 0.
603 indifl = indisp + n
604* INDIWO is the offset of the remaining integer workspace.
605 indiwo = indifl + n
606
607*
608* Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
609*
610*
611 CALL dsytrd_2stage( jobz, uplo, n, a, lda, work( indd ),
612 $ work( inde ), work( indtau ), work( indhous ),
613 $ lhtrd, work( indwk ), llwork, iinfo )
614*
615* If all eigenvalues are desired
616* then call DSTERF or DSTEMR and DORMTR.
617*
618 IF( ( alleig .OR. ( indeig .AND. il.EQ.1 .AND. iu.EQ.n ) ) .AND.
619 $ ieeeok.EQ.1 ) THEN
620 IF( .NOT.wantz ) THEN
621 CALL dcopy( n, work( indd ), 1, w, 1 )
622 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
623 CALL dsterf( n, w, work( indee ), info )
624 ELSE
625 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
626 CALL dcopy( n, work( indd ), 1, work( inddd ), 1 )
627*
628 IF (abstol .LE. two*n*eps) THEN
629 tryrac = .true.
630 ELSE
631 tryrac = .false.
632 END IF
633 CALL dstemr( jobz, 'A', n, work( inddd ), work( indee ),
634 $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
635 $ tryrac, work( indwk ), lwork, iwork, liwork,
636 $ info )
637*
638*
639*
640* Apply orthogonal matrix used in reduction to tridiagonal
641* form to eigenvectors returned by DSTEMR.
642*
643 IF( wantz .AND. info.EQ.0 ) THEN
644 indwkn = inde
645 llwrkn = lwork - indwkn + 1
646 CALL dormtr( 'L', uplo, 'N', n, m, a, lda,
647 $ work( indtau ), z, ldz, work( indwkn ),
648 $ llwrkn, iinfo )
649 END IF
650 END IF
651*
652*
653 IF( info.EQ.0 ) THEN
654* Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
655* undefined.
656 m = n
657 GO TO 30
658 END IF
659 info = 0
660 END IF
661*
662* Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
663* Also call DSTEBZ and DSTEIN if DSTEMR fails.
664*
665 IF( wantz ) THEN
666 order = 'B'
667 ELSE
668 order = 'E'
669 END IF
670
671 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
672 $ work( indd ), work( inde ), m, nsplit, w,
673 $ iwork( indibl ), iwork( indisp ), work( indwk ),
674 $ iwork( indiwo ), info )
675*
676 IF( wantz ) THEN
677 CALL dstein( n, work( indd ), work( inde ), m, w,
678 $ iwork( indibl ), iwork( indisp ), z, ldz,
679 $ work( indwk ), iwork( indiwo ), iwork( indifl ),
680 $ info )
681*
682* Apply orthogonal matrix used in reduction to tridiagonal
683* form to eigenvectors returned by DSTEIN.
684*
685 indwkn = inde
686 llwrkn = lwork - indwkn + 1
687 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
688 $ ldz, work( indwkn ), llwrkn, iinfo )
689 END IF
690*
691* If matrix was scaled, then rescale eigenvalues appropriately.
692*
693* Jump here if DSTEMR/DSTEIN succeeded.
694 30 CONTINUE
695 IF( iscale.EQ.1 ) THEN
696 IF( info.EQ.0 ) THEN
697 imax = m
698 ELSE
699 imax = info - 1
700 END IF
701 CALL dscal( imax, one / sigma, w, 1 )
702 END IF
703*
704* If eigenvalues are not in order, then sort them, along with
705* eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
706* It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
707* not return this detailed information to the user.
708*
709 IF( wantz ) THEN
710 DO 50 j = 1, m - 1
711 i = 0
712 tmp1 = w( j )
713 DO 40 jj = j + 1, m
714 IF( w( jj ).LT.tmp1 ) THEN
715 i = jj
716 tmp1 = w( jj )
717 END IF
718 40 CONTINUE
719*
720 IF( i.NE.0 ) THEN
721 w( i ) = w( j )
722 w( j ) = tmp1
723 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
724 END IF
725 50 CONTINUE
726 END IF
727*
728* Set WORK(1) to optimal workspace size.
729*
730 work( 1 ) = lwmin
731 iwork( 1 ) = liwmin
732*
733 RETURN
734*
735* End of DSYEVR_2STAGE
736*
737 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:273
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:174
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:322
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
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:171