LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsyevx_2stage.f
Go to the documentation of this file.
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*> \ingroup heevx_2stage
263*
264*> \par Further Details:
265* =====================
266*>
267*> \verbatim
268*>
269*> All details about the 2stage techniques are available in:
270*>
271*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
272*> Parallel reduction to condensed forms for symmetric eigenvalue problems
273*> using aggregated fine-grained and memory-aware kernels. In Proceedings
274*> of 2011 International Conference for High Performance Computing,
275*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
276*> Article 8 , 11 pages.
277*> http://doi.acm.org/10.1145/2063384.2063394
278*>
279*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
280*> An improved parallel singular value algorithm and its implementation
281*> for multicore hardware, In Proceedings of 2013 International Conference
282*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
283*> Denver, Colorado, USA, 2013.
284*> Article 90, 12 pages.
285*> http://doi.acm.org/10.1145/2503210.2503292
286*>
287*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
288*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
289*> calculations based on fine-grained memory aware tasks.
290*> International Journal of High Performance Computing Applications.
291*> Volume 28 Issue 2, Pages 196-209, May 2014.
292*> http://hpc.sagepub.com/content/28/2/196
293*>
294*> \endverbatim
295*
296* =====================================================================
297 SUBROUTINE dsyevx_2stage( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
298 $ IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
299 $ LWORK, IWORK, IFAIL, INFO )
300*
301 IMPLICIT NONE
302*
303* -- LAPACK driver routine --
304* -- LAPACK is a software package provided by Univ. of Tennessee, --
305* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
306*
307* .. Scalar Arguments ..
308 CHARACTER JOBZ, RANGE, UPLO
309 INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
310 DOUBLE PRECISION ABSTOL, VL, VU
311* ..
312* .. Array Arguments ..
313 INTEGER IFAIL( * ), IWORK( * )
314 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
315* ..
316*
317* =====================================================================
318*
319* .. Parameters ..
320 DOUBLE PRECISION ZERO, ONE
321 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
322* ..
323* .. Local Scalars ..
324 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
325 $ WANTZ
326 CHARACTER ORDER
327 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
328 $ indisp, indiwo, indtau, indwkn, indwrk, iscale,
329 $ itmp1, j, jj, llwork, llwrkn,
330 $ nsplit, lwmin, lhtrd, lwtrd, kd, ib, indhous
331 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
332 $ SIGMA, SMLNUM, TMP1, VLL, VUU
333* ..
334* .. External Functions ..
335 LOGICAL LSAME
336 INTEGER ILAENV2STAGE
337 DOUBLE PRECISION DLAMCH, DLANSY
338 EXTERNAL lsame, dlamch, dlansy, ilaenv2stage
339* ..
340* .. External Subroutines ..
341 EXTERNAL dcopy, dlacpy, dorgtr, dormtr, dscal, dstebz,
344* ..
345* .. Intrinsic Functions ..
346 INTRINSIC max, min, sqrt
347* ..
348* .. Executable Statements ..
349*
350* Test the input parameters.
351*
352 lower = lsame( uplo, 'L' )
353 wantz = lsame( jobz, 'V' )
354 alleig = lsame( range, 'A' )
355 valeig = lsame( range, 'V' )
356 indeig = lsame( range, 'I' )
357 lquery = ( lwork.EQ.-1 )
358*
359 info = 0
360 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
361 info = -1
362 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
363 info = -2
364 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
365 info = -3
366 ELSE IF( n.LT.0 ) THEN
367 info = -4
368 ELSE IF( lda.LT.max( 1, n ) ) THEN
369 info = -6
370 ELSE
371 IF( valeig ) THEN
372 IF( n.GT.0 .AND. vu.LE.vl )
373 $ info = -8
374 ELSE IF( indeig ) THEN
375 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
376 info = -9
377 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
378 info = -10
379 END IF
380 END IF
381 END IF
382 IF( info.EQ.0 ) THEN
383 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
384 info = -15
385 END IF
386 END IF
387*
388 IF( info.EQ.0 ) THEN
389 IF( n.LE.1 ) THEN
390 lwmin = 1
391 work( 1 ) = lwmin
392 ELSE
393 kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz,
394 $ n, -1, -1, -1 )
395 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz,
396 $ n, kd, -1, -1 )
397 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz,
398 $ n, kd, ib, -1 )
399 lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz,
400 $ n, kd, ib, -1 )
401 lwmin = max( 8*n, 3*n + lhtrd + lwtrd )
402 work( 1 ) = lwmin
403 END IF
404*
405 IF( lwork.LT.lwmin .AND. .NOT.lquery )
406 $ info = -17
407 END IF
408*
409 IF( info.NE.0 ) THEN
410 CALL xerbla( 'DSYEVX_2STAGE', -info )
411 RETURN
412 ELSE IF( lquery ) THEN
413 RETURN
414 END IF
415*
416* Quick return if possible
417*
418 m = 0
419 IF( n.EQ.0 ) THEN
420 RETURN
421 END IF
422*
423 IF( n.EQ.1 ) THEN
424 IF( alleig .OR. indeig ) THEN
425 m = 1
426 w( 1 ) = a( 1, 1 )
427 ELSE
428 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
429 m = 1
430 w( 1 ) = a( 1, 1 )
431 END IF
432 END IF
433 IF( wantz )
434 $ z( 1, 1 ) = one
435 RETURN
436 END IF
437*
438* Get machine constants.
439*
440 safmin = dlamch( 'Safe minimum' )
441 eps = dlamch( 'Precision' )
442 smlnum = safmin / eps
443 bignum = one / smlnum
444 rmin = sqrt( smlnum )
445 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
446*
447* Scale matrix to allowable range, if necessary.
448*
449 iscale = 0
450 abstll = abstol
451 IF( valeig ) THEN
452 vll = vl
453 vuu = vu
454 END IF
455 anrm = dlansy( 'M', uplo, n, a, lda, work )
456 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
457 iscale = 1
458 sigma = rmin / anrm
459 ELSE IF( anrm.GT.rmax ) THEN
460 iscale = 1
461 sigma = rmax / anrm
462 END IF
463 IF( iscale.EQ.1 ) THEN
464 IF( lower ) THEN
465 DO 10 j = 1, n
466 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
467 10 CONTINUE
468 ELSE
469 DO 20 j = 1, n
470 CALL dscal( j, sigma, a( 1, j ), 1 )
471 20 CONTINUE
472 END IF
473 IF( abstol.GT.0 )
474 $ abstll = abstol*sigma
475 IF( valeig ) THEN
476 vll = vl*sigma
477 vuu = vu*sigma
478 END IF
479 END IF
480*
481* Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
482*
483 indtau = 1
484 inde = indtau + n
485 indd = inde + n
486 indhous = indd + n
487 indwrk = indhous + lhtrd
488 llwork = lwork - indwrk + 1
489*
490 CALL dsytrd_2stage( jobz, uplo, n, a, lda, work( indd ),
491 $ work( inde ), work( indtau ), work( indhous ),
492 $ lhtrd, work( indwrk ), llwork, iinfo )
493*
494* If all eigenvalues are desired and ABSTOL is less than or equal to
495* zero, then call DSTERF or DORGTR and SSTEQR. If this fails for
496* some eigenvalue, then try DSTEBZ.
497*
498 test = .false.
499 IF( indeig ) THEN
500 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
501 test = .true.
502 END IF
503 END IF
504 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
505 CALL dcopy( n, work( indd ), 1, w, 1 )
506 indee = indwrk + 2*n
507 IF( .NOT.wantz ) THEN
508 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
509 CALL dsterf( n, w, work( indee ), info )
510 ELSE
511 CALL dlacpy( 'A', n, n, a, lda, z, ldz )
512 CALL dorgtr( uplo, n, z, ldz, work( indtau ),
513 $ work( indwrk ), llwork, iinfo )
514 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
515 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
516 $ work( indwrk ), info )
517 IF( info.EQ.0 ) THEN
518 DO 30 i = 1, n
519 ifail( i ) = 0
520 30 CONTINUE
521 END IF
522 END IF
523 IF( info.EQ.0 ) THEN
524 m = n
525 GO TO 40
526 END IF
527 info = 0
528 END IF
529*
530* Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
531*
532 IF( wantz ) THEN
533 order = 'B'
534 ELSE
535 order = 'E'
536 END IF
537 indibl = 1
538 indisp = indibl + n
539 indiwo = indisp + n
540 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
541 $ work( indd ), work( inde ), m, nsplit, w,
542 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
543 $ iwork( indiwo ), info )
544*
545 IF( wantz ) THEN
546 CALL dstein( n, work( indd ), work( inde ), m, w,
547 $ iwork( indibl ), iwork( indisp ), z, ldz,
548 $ work( indwrk ), iwork( indiwo ), ifail, info )
549*
550* Apply orthogonal matrix used in reduction to tridiagonal
551* form to eigenvectors returned by DSTEIN.
552*
553 indwkn = inde
554 llwrkn = lwork - indwkn + 1
555 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
556 $ ldz, work( indwkn ), llwrkn, iinfo )
557 END IF
558*
559* If matrix was scaled, then rescale eigenvalues appropriately.
560*
561 40 CONTINUE
562 IF( iscale.EQ.1 ) THEN
563 IF( info.EQ.0 ) THEN
564 imax = m
565 ELSE
566 imax = info - 1
567 END IF
568 CALL dscal( imax, one / sigma, w, 1 )
569 END IF
570*
571* If eigenvalues are not in order, then sort them, along with
572* eigenvectors.
573*
574 IF( wantz ) THEN
575 DO 60 j = 1, m - 1
576 i = 0
577 tmp1 = w( j )
578 DO 50 jj = j + 1, m
579 IF( w( jj ).LT.tmp1 ) THEN
580 i = jj
581 tmp1 = w( jj )
582 END IF
583 50 CONTINUE
584*
585 IF( i.NE.0 ) THEN
586 itmp1 = iwork( indibl+i-1 )
587 w( i ) = w( j )
588 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
589 w( j ) = tmp1
590 iwork( indibl+j-1 ) = itmp1
591 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
592 IF( info.NE.0 ) THEN
593 itmp1 = ifail( i )
594 ifail( i ) = ifail( j )
595 ifail( j ) = itmp1
596 END IF
597 END IF
598 60 CONTINUE
599 END IF
600*
601* Set WORK(1) to optimal workspace size.
602*
603 work( 1 ) = lwmin
604*
605 RETURN
606*
607* End of DSYEVX_2STAGE
608*
609 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dsyevx_2stage(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, iwork, ifail, info)
DSYEVX_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 dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
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 dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:131
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dorgtr(uplo, n, a, lda, tau, work, lwork, info)
DORGTR
Definition dorgtr.f:123
subroutine dormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
DORMTR
Definition dormtr.f:171