LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dsyevx_2stage()

subroutine dsyevx_2stage ( character jobz,
character range,
character uplo,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
double precision vl,
double precision vu,
integer il,
integer iu,
double precision abstol,
integer m,
double precision, dimension( * ) w,
double precision, dimension( ldz, * ) z,
integer ldz,
double precision, dimension( * ) work,
integer lwork,
integer, dimension( * ) iwork,
integer, dimension( * ) ifail,
integer info )

DSYEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices

Download DSYEVX_2STAGE + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DSYEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
!> of a real symmetric matrix A using the 2stage technique for
!> the reduction to tridiagonal.  Eigenvalues and eigenvectors can be
!> selected by specifying either a range of values or a range of indices
!> for the desired eigenvalues.
!> 
Parameters
[in]JOBZ
!>          JOBZ is CHARACTER*1
!>          = 'N':  Compute eigenvalues only;
!>          = 'V':  Compute eigenvalues and eigenvectors.
!>                  Not available in this release.
!> 
[in]RANGE
!>          RANGE is CHARACTER*1
!>          = 'A': all eigenvalues will be found.
!>          = 'V': all eigenvalues in the half-open interval (VL,VU]
!>                 will be found.
!>          = 'I': the IL-th through IU-th eigenvalues will be found.
!> 
[in]UPLO
!>          UPLO is CHARACTER*1
!>          = 'U':  Upper triangle of A is stored;
!>          = 'L':  Lower triangle of A is stored.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array, dimension (LDA, N)
!>          On entry, the symmetric matrix A.  If UPLO = 'U', the
!>          leading N-by-N upper triangular part of A contains the
!>          upper triangular part of the matrix A.  If UPLO = 'L',
!>          the leading N-by-N lower triangular part of A contains
!>          the lower triangular part of the matrix A.
!>          On exit, the lower triangle (if UPLO='L') or the upper
!>          triangle (if UPLO='U') of A, including the diagonal, is
!>          destroyed.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in]VL
!>          VL is DOUBLE PRECISION
!>          If RANGE='V', the lower bound of the interval to
!>          be searched for eigenvalues. VL < VU.
!>          Not referenced if RANGE = 'A' or 'I'.
!> 
[in]VU
!>          VU is DOUBLE PRECISION
!>          If RANGE='V', the upper bound of the interval to
!>          be searched for eigenvalues. VL < VU.
!>          Not referenced if RANGE = 'A' or 'I'.
!> 
[in]IL
!>          IL is INTEGER
!>          If RANGE='I', the index of the
!>          smallest eigenvalue to be returned.
!>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
!>          Not referenced if RANGE = 'A' or 'V'.
!> 
[in]IU
!>          IU is INTEGER
!>          If RANGE='I', the index of the
!>          largest eigenvalue to be returned.
!>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
!>          Not referenced if RANGE = 'A' or 'V'.
!> 
[in]ABSTOL
!>          ABSTOL is DOUBLE PRECISION
!>          The absolute error tolerance for the eigenvalues.
!>          An approximate eigenvalue is accepted as converged
!>          when it is determined to lie in an interval [a,b]
!>          of width less than or equal to
!>
!>                  ABSTOL + EPS *   max( |a|,|b| ) ,
!>
!>          where EPS is the machine precision.  If ABSTOL is less than
!>          or equal to zero, then  EPS*|T|  will be used in its place,
!>          where |T| is the 1-norm of the tridiagonal matrix obtained
!>          by reducing A to tridiagonal form.
!>
!>          Eigenvalues will be computed most accurately when ABSTOL is
!>          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
!>          If this routine returns with INFO>0, indicating that some
!>          eigenvectors did not converge, try setting ABSTOL to
!>          2*DLAMCH('S').
!>
!>          See  by Demmel and
!>          Kahan, LAPACK Working Note #3.
!> 
[out]M
!>          M is INTEGER
!>          The total number of eigenvalues found.  0 <= M <= N.
!>          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
!> 
[out]W
!>          W is DOUBLE PRECISION array, dimension (N)
!>          On normal exit, the first M elements contain the selected
!>          eigenvalues in ascending order.
!> 
[out]Z
!>          Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
!>          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
!>          contain the orthonormal eigenvectors of the matrix A
!>          corresponding to the selected eigenvalues, with the i-th
!>          column of Z holding the eigenvector associated with W(i).
!>          If an eigenvector fails to converge, then that column of Z
!>          contains the latest approximation to the eigenvector, and the
!>          index of the eigenvector is returned in IFAIL.
!>          If JOBZ = 'N', then Z is not referenced.
!>          Note: the user must ensure that at least max(1,M) columns are
!>          supplied in the array Z; if RANGE = 'V', the exact value of M
!>          is not known in advance and an upper bound must be used.
!> 
[in]LDZ
!>          LDZ is INTEGER
!>          The leading dimension of the array Z.  LDZ >= 1, and if
!>          JOBZ = 'V', LDZ >= max(1,N).
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of the array WORK. LWORK >= 1, when N <= 1;
!>          otherwise
!>          If JOBZ = 'N' and N > 1, LWORK must be queried.
!>                                   LWORK = MAX(1, 8*N, dimension) where
!>                                   dimension = max(stage1,stage2) + (KD+1)*N + 3*N
!>                                             = N*KD + N*max(KD+1,FACTOPTNB)
!>                                               + max(2*KD*KD, KD*NTHREADS)
!>                                               + (KD+1)*N + 3*N
!>                                   where KD is the blocking size of the reduction,
!>                                   FACTOPTNB is the blocking used by the QR or LQ
!>                                   algorithm, usually FACTOPTNB=128 is a good choice
!>                                   NTHREADS is the number of threads used when
!>                                   openMP compilation is enabled, otherwise =1.
!>          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (5*N)
!> 
[out]IFAIL
!>          IFAIL is INTEGER array, dimension (N)
!>          If JOBZ = 'V', then if INFO = 0, the first M elements of
!>          IFAIL are zero.  If INFO > 0, then IFAIL contains the
!>          indices of the eigenvectors that failed to converge.
!>          If JOBZ = 'N', then IFAIL is not referenced.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  if INFO = i, then i eigenvectors failed to converge.
!>                Their indices are stored in array IFAIL.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  All details about the 2stage techniques are available in:
!>
!>  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
!>  Parallel reduction to condensed forms for symmetric eigenvalue problems
!>  using aggregated fine-grained and memory-aware kernels. In Proceedings
!>  of 2011 International Conference for High Performance Computing,
!>  Networking, Storage and Analysis (SC '11), New York, NY, USA,
!>  Article 8 , 11 pages.
!>  http://doi.acm.org/10.1145/2063384.2063394
!>
!>  A. Haidar, J. Kurzak, P. Luszczek, 2013.
!>  An improved parallel singular value algorithm and its implementation
!>  for multicore hardware, In Proceedings of 2013 International Conference
!>  for High Performance Computing, Networking, Storage and Analysis (SC '13).
!>  Denver, Colorado, USA, 2013.
!>  Article 90, 12 pages.
!>  http://doi.acm.org/10.1145/2503210.2503292
!>
!>  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
!>  A novel hybrid CPU-GPU generalized eigensolver for electronic structure
!>  calculations based on fine-grained memory aware tasks.
!>  International Journal of High Performance Computing Applications.
!>  Volume 28 Issue 2, Pages 196-209, May 2014.
!>  http://hpc.sagepub.com/content/28/2/196
!>
!> 

Definition at line 295 of file dsyevx_2stage.f.

298*
299 IMPLICIT NONE
300*
301* -- LAPACK driver routine --
302* -- LAPACK is a software package provided by Univ. of Tennessee, --
303* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
304*
305* .. Scalar Arguments ..
306 CHARACTER JOBZ, RANGE, UPLO
307 INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
308 DOUBLE PRECISION ABSTOL, VL, VU
309* ..
310* .. Array Arguments ..
311 INTEGER IFAIL( * ), IWORK( * )
312 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
313* ..
314*
315* =====================================================================
316*
317* .. Parameters ..
318 DOUBLE PRECISION ZERO, ONE
319 parameter( zero = 0.0d+0, one = 1.0d+0 )
320* ..
321* .. Local Scalars ..
322 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
323 $ WANTZ
324 CHARACTER ORDER
325 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
326 $ INDISP, INDIWO, INDTAU, INDWKN, INDWRK, ISCALE,
327 $ ITMP1, J, JJ, LLWORK, LLWRKN,
328 $ NSPLIT, LWMIN, LHTRD, LWTRD, KD, IB, INDHOUS
329 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
330 $ SIGMA, SMLNUM, TMP1, VLL, VUU
331* ..
332* .. External Functions ..
333 LOGICAL LSAME
334 INTEGER ILAENV2STAGE
335 DOUBLE PRECISION DLAMCH, DLANSY
336 EXTERNAL lsame, dlamch, dlansy, ilaenv2stage
337* ..
338* .. External Subroutines ..
339 EXTERNAL dcopy, dlacpy, dorgtr, dormtr, dscal,
340 $ dstebz,
343* ..
344* .. Intrinsic Functions ..
345 INTRINSIC max, min, sqrt
346* ..
347* .. Executable Statements ..
348*
349* Test the input parameters.
350*
351 lower = lsame( uplo, 'L' )
352 wantz = lsame( jobz, 'V' )
353 alleig = lsame( range, 'A' )
354 valeig = lsame( range, 'V' )
355 indeig = lsame( range, 'I' )
356 lquery = ( lwork.EQ.-1 )
357*
358 info = 0
359 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
360 info = -1
361 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
362 info = -2
363 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
364 info = -3
365 ELSE IF( n.LT.0 ) THEN
366 info = -4
367 ELSE IF( lda.LT.max( 1, n ) ) THEN
368 info = -6
369 ELSE
370 IF( valeig ) THEN
371 IF( n.GT.0 .AND. vu.LE.vl )
372 $ info = -8
373 ELSE IF( indeig ) THEN
374 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
375 info = -9
376 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
377 info = -10
378 END IF
379 END IF
380 END IF
381 IF( info.EQ.0 ) THEN
382 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
383 info = -15
384 END IF
385 END IF
386*
387 IF( info.EQ.0 ) THEN
388 IF( n.LE.1 ) THEN
389 lwmin = 1
390 work( 1 ) = lwmin
391 ELSE
392 kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz,
393 $ n, -1, -1, -1 )
394 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz,
395 $ n, kd, -1, -1 )
396 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz,
397 $ n, kd, ib, -1 )
398 lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz,
399 $ n, kd, ib, -1 )
400 lwmin = max( 8*n, 3*n + lhtrd + lwtrd )
401 work( 1 ) = lwmin
402 END IF
403*
404 IF( lwork.LT.lwmin .AND. .NOT.lquery )
405 $ info = -17
406 END IF
407*
408 IF( info.NE.0 ) THEN
409 CALL xerbla( 'DSYEVX_2STAGE', -info )
410 RETURN
411 ELSE IF( lquery ) THEN
412 RETURN
413 END IF
414*
415* Quick return if possible
416*
417 m = 0
418 IF( n.EQ.0 ) THEN
419 RETURN
420 END IF
421*
422 IF( n.EQ.1 ) THEN
423 IF( alleig .OR. indeig ) THEN
424 m = 1
425 w( 1 ) = a( 1, 1 )
426 ELSE
427 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
428 m = 1
429 w( 1 ) = a( 1, 1 )
430 END IF
431 END IF
432 IF( wantz )
433 $ z( 1, 1 ) = one
434 RETURN
435 END IF
436*
437* Get machine constants.
438*
439 safmin = dlamch( 'Safe minimum' )
440 eps = dlamch( 'Precision' )
441 smlnum = safmin / eps
442 bignum = one / smlnum
443 rmin = sqrt( smlnum )
444 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
445*
446* Scale matrix to allowable range, if necessary.
447*
448 iscale = 0
449 abstll = abstol
450 IF( valeig ) THEN
451 vll = vl
452 vuu = vu
453 END IF
454 anrm = dlansy( 'M', uplo, n, a, lda, work )
455 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
456 iscale = 1
457 sigma = rmin / anrm
458 ELSE IF( anrm.GT.rmax ) THEN
459 iscale = 1
460 sigma = rmax / anrm
461 END IF
462 IF( iscale.EQ.1 ) THEN
463 IF( lower ) THEN
464 DO 10 j = 1, n
465 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
466 10 CONTINUE
467 ELSE
468 DO 20 j = 1, n
469 CALL dscal( j, sigma, a( 1, j ), 1 )
470 20 CONTINUE
471 END IF
472 IF( abstol.GT.0 )
473 $ abstll = abstol*sigma
474 IF( valeig ) THEN
475 vll = vl*sigma
476 vuu = vu*sigma
477 END IF
478 END IF
479*
480* Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
481*
482 indtau = 1
483 inde = indtau + n
484 indd = inde + n
485 indhous = indd + n
486 indwrk = indhous + lhtrd
487 llwork = lwork - indwrk + 1
488*
489 CALL dsytrd_2stage( jobz, uplo, n, a, lda, work( indd ),
490 $ work( inde ), work( indtau ), work( indhous ),
491 $ lhtrd, work( indwrk ), llwork, iinfo )
492*
493* If all eigenvalues are desired and ABSTOL is less than or equal to
494* zero, then call DSTERF or DORGTR and SSTEQR. If this fails for
495* some eigenvalue, then try DSTEBZ.
496*
497 test = .false.
498 IF( indeig ) THEN
499 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
500 test = .true.
501 END IF
502 END IF
503 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
504 CALL dcopy( n, work( indd ), 1, w, 1 )
505 indee = indwrk + 2*n
506 IF( .NOT.wantz ) THEN
507 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
508 CALL dsterf( n, w, work( indee ), info )
509 ELSE
510 CALL dlacpy( 'A', n, n, a, lda, z, ldz )
511 CALL dorgtr( uplo, n, z, ldz, work( indtau ),
512 $ work( indwrk ), llwork, iinfo )
513 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
514 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
515 $ work( indwrk ), info )
516 IF( info.EQ.0 ) THEN
517 DO 30 i = 1, n
518 ifail( i ) = 0
519 30 CONTINUE
520 END IF
521 END IF
522 IF( info.EQ.0 ) THEN
523 m = n
524 GO TO 40
525 END IF
526 info = 0
527 END IF
528*
529* Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
530*
531 IF( wantz ) THEN
532 order = 'B'
533 ELSE
534 order = 'E'
535 END IF
536 indibl = 1
537 indisp = indibl + n
538 indiwo = indisp + n
539 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
540 $ work( indd ), work( inde ), m, nsplit, w,
541 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
542 $ iwork( indiwo ), info )
543*
544 IF( wantz ) THEN
545 CALL dstein( n, work( indd ), work( inde ), m, w,
546 $ iwork( indibl ), iwork( indisp ), z, ldz,
547 $ work( indwrk ), iwork( indiwo ), ifail, info )
548*
549* Apply orthogonal matrix used in reduction to tridiagonal
550* form to eigenvectors returned by DSTEIN.
551*
552 indwkn = inde
553 llwrkn = lwork - indwkn + 1
554 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
555 $ 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*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dsytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
DSYTRD_2STAGE
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlansy(norm, uplo, n, a, lda, work)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansy.f:120
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
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 dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:129
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84
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:121
subroutine dormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
DORMTR
Definition dormtr.f:170
Here is the call graph for this function:
Here is the caller graph for this function: