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

◆ dsygvx()

subroutine dsygvx ( integer itype,
character jobz,
character range,
character uplo,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldb, * ) b,
integer ldb,
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 )

DSYGVX

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

Purpose:
!>
!> DSYGVX computes selected eigenvalues, and optionally, eigenvectors
!> of a real generalized symmetric-definite eigenproblem, of the form
!> A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.  Here A
!> and B are assumed to be symmetric and B is also positive definite.
!> Eigenvalues and eigenvectors can be selected by specifying either a
!> range of values or a range of indices for the desired eigenvalues.
!> 
Parameters
[in]ITYPE
!>          ITYPE is INTEGER
!>          Specifies the problem type to be solved:
!>          = 1:  A*x = (lambda)*B*x
!>          = 2:  A*B*x = (lambda)*x
!>          = 3:  B*A*x = (lambda)*x
!> 
[in]JOBZ
!>          JOBZ is CHARACTER*1
!>          = 'N':  Compute eigenvalues only;
!>          = 'V':  Compute eigenvalues and eigenvectors.
!> 
[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 and B are stored;
!>          = 'L':  Lower triangle of A and B are stored.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix pencil (A,B).  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,out]B
!>          B is DOUBLE PRECISION array, dimension (LDB, N)
!>          On entry, the symmetric matrix B.  If UPLO = 'U', the
!>          leading N-by-N upper triangular part of B contains the
!>          upper triangular part of the matrix B.  If UPLO = 'L',
!>          the leading N-by-N lower triangular part of B contains
!>          the lower triangular part of the matrix B.
!>
!>          On exit, if INFO <= N, the part of B containing the matrix is
!>          overwritten by the triangular factor U or L from the Cholesky
!>          factorization B = U**T*U or B = L*L**T.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= 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 C to tridiagonal form, where C is the symmetric
!>          matrix of the standard symmetric problem to which the
!>          generalized problem is transformed.
!>
!>          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').
!> 
[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 = 'N', then Z is not referenced.
!>          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).
!>          The eigenvectors are normalized as follows:
!>          if ITYPE = 1 or 2, Z**T*B*Z = I;
!>          if ITYPE = 3, Z**T*inv(B)*Z = 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.
!>          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 >= max(1,8*N).
!>          For optimal efficiency, LWORK >= (NB+3)*N,
!>          where NB is the blocksize for DSYTRD returned by ILAENV.
!>
!>          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:  DPOTRF or DSYEVX returned an error code:
!>             <= N:  if INFO = i, DSYEVX failed to converge;
!>                    i eigenvectors failed to converge.  Their indices
!>                    are stored in array IFAIL.
!>             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
!>                    principal minor of order i of B is not positive.
!>                    The factorization of B could not be completed and
!>                    no eigenvalues or eigenvectors were computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 292 of file dsygvx.f.

295*
296* -- LAPACK driver routine --
297* -- LAPACK is a software package provided by Univ. of Tennessee, --
298* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
299*
300* .. Scalar Arguments ..
301 CHARACTER JOBZ, RANGE, UPLO
302 INTEGER IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N
303 DOUBLE PRECISION ABSTOL, VL, VU
304* ..
305* .. Array Arguments ..
306 INTEGER IFAIL( * ), IWORK( * )
307 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * ),
308 $ Z( LDZ, * )
309* ..
310*
311* =====================================================================
312*
313* .. Parameters ..
314 DOUBLE PRECISION ONE
315 parameter( one = 1.0d+0 )
316* ..
317* .. Local Scalars ..
318 LOGICAL ALLEIG, INDEIG, LQUERY, UPPER, VALEIG, WANTZ
319 CHARACTER TRANS
320 INTEGER LWKMIN, LWKOPT, NB
321* ..
322* .. External Functions ..
323 LOGICAL LSAME
324 INTEGER ILAENV
325 EXTERNAL lsame, ilaenv
326* ..
327* .. External Subroutines ..
328 EXTERNAL dpotrf, dsyevx, dsygst, dtrmm, dtrsm,
329 $ xerbla
330* ..
331* .. Intrinsic Functions ..
332 INTRINSIC max, min
333* ..
334* .. Executable Statements ..
335*
336* Test the input parameters.
337*
338 upper = lsame( uplo, 'U' )
339 wantz = lsame( jobz, 'V' )
340 alleig = lsame( range, 'A' )
341 valeig = lsame( range, 'V' )
342 indeig = lsame( range, 'I' )
343 lquery = ( lwork.EQ.-1 )
344*
345 info = 0
346 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
347 info = -1
348 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
349 info = -2
350 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
351 info = -3
352 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
353 info = -4
354 ELSE IF( n.LT.0 ) THEN
355 info = -5
356 ELSE IF( lda.LT.max( 1, n ) ) THEN
357 info = -7
358 ELSE IF( ldb.LT.max( 1, n ) ) THEN
359 info = -9
360 ELSE
361 IF( valeig ) THEN
362 IF( n.GT.0 .AND. vu.LE.vl )
363 $ info = -11
364 ELSE IF( indeig ) THEN
365 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
366 info = -12
367 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
368 info = -13
369 END IF
370 END IF
371 END IF
372 IF (info.EQ.0) THEN
373 IF (ldz.LT.1 .OR. (wantz .AND. ldz.LT.n)) THEN
374 info = -18
375 END IF
376 END IF
377*
378 IF( info.EQ.0 ) THEN
379 lwkmin = max( 1, 8*n )
380 nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
381 lwkopt = max( lwkmin, ( nb + 3 )*n )
382 work( 1 ) = lwkopt
383*
384 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
385 info = -20
386 END IF
387 END IF
388*
389 IF( info.NE.0 ) THEN
390 CALL xerbla( 'DSYGVX', -info )
391 RETURN
392 ELSE IF( lquery ) THEN
393 RETURN
394 END IF
395*
396* Quick return if possible
397*
398 m = 0
399 IF( n.EQ.0 ) THEN
400 RETURN
401 END IF
402*
403* Form a Cholesky factorization of B.
404*
405 CALL dpotrf( uplo, n, b, ldb, info )
406 IF( info.NE.0 ) THEN
407 info = n + info
408 RETURN
409 END IF
410*
411* Transform problem to standard eigenvalue problem and solve.
412*
413 CALL dsygst( itype, uplo, n, a, lda, b, ldb, info )
414 CALL dsyevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu,
415 $ abstol,
416 $ m, w, z, ldz, work, lwork, iwork, ifail, info )
417*
418 IF( wantz ) THEN
419*
420* Backtransform eigenvectors to the original problem.
421*
422 IF( info.GT.0 )
423 $ m = info - 1
424 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
425*
426* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
427* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
428*
429 IF( upper ) THEN
430 trans = 'N'
431 ELSE
432 trans = 'T'
433 END IF
434*
435 CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, m, one,
436 $ b,
437 $ ldb, z, ldz )
438*
439 ELSE IF( itype.EQ.3 ) THEN
440*
441* For B*A*x=(lambda)*x;
442* backtransform eigenvectors: x = L*y or U**T*y
443*
444 IF( upper ) THEN
445 trans = 'T'
446 ELSE
447 trans = 'N'
448 END IF
449*
450 CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, m, one,
451 $ b,
452 $ ldb, z, ldz )
453 END IF
454 END IF
455*
456* Set WORK(1) to optimal workspace size.
457*
458 work( 1 ) = lwkopt
459*
460 RETURN
461*
462* End of DSYGVX
463*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsyevx(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, iwork, ifail, info)
DSYEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Definition dsyevx.f:252
subroutine dsygst(itype, uplo, n, a, lda, b, ldb, info)
DSYGST
Definition dsygst.f:125
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:105
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: