Page MenuHomec4science

dsygvd.f
No OneTemporary

File Metadata

Created
Sun, Jun 9, 14:10

dsygvd.f

*> \brief \b DSYGST
*
* =========== DOCUMENTATION ===========
*
* Online html documentation available at
* http://www.netlib.org/lapack/explore-html/
*
*> \htmlonly
*> Download DSYGVD + dependencies
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsygvd.f">
*> [TGZ]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsygvd.f">
*> [ZIP]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsygvd.f">
*> [TXT]</a>
*> \endhtmlonly
*
* Definition:
* ===========
*
* SUBROUTINE DSYGVD( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
* LWORK, IWORK, LIWORK, INFO )
*
* .. Scalar Arguments ..
* CHARACTER JOBZ, UPLO
* INTEGER INFO, ITYPE, LDA, LDB, LIWORK, LWORK, N
* ..
* .. Array Arguments ..
* INTEGER IWORK( * )
* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
* ..
*
*
*> \par Purpose:
* =============
*>
*> \verbatim
*>
*> DSYGVD computes all the eigenvalues, and optionally, the 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.
*> If eigenvectors are desired, it uses a divide and conquer algorithm.
*>
*> The divide and conquer algorithm makes very mild assumptions about
*> floating point arithmetic. It will work on machines with a guard
*> digit in add/subtract, or on those binary machines without guard
*> digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
*> Cray-2. It could conceivably fail on hexadecimal or decimal machines
*> without guard digits, but we know of none.
*> \endverbatim
*
* Arguments:
* ==========
*
*> \param[in] ITYPE
*> \verbatim
*> 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
*> \endverbatim
*>
*> \param[in] JOBZ
*> \verbatim
*> JOBZ is CHARACTER*1
*> = 'N': Compute eigenvalues only;
*> = 'V': Compute eigenvalues and eigenvectors.
*> \endverbatim
*>
*> \param[in] UPLO
*> \verbatim
*> UPLO is CHARACTER*1
*> = 'U': Upper triangles of A and B are stored;
*> = 'L': Lower triangles of A and B are stored.
*> \endverbatim
*>
*> \param[in] N
*> \verbatim
*> N is INTEGER
*> The order of the matrices A and B. N >= 0.
*> \endverbatim
*>
*> \param[in,out] A
*> \verbatim
*> 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, if JOBZ = 'V', then if INFO = 0, A contains the
*> matrix Z of eigenvectors. 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 JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
*> or the lower triangle (if UPLO='L') of A, including the
*> diagonal, is destroyed.
*> \endverbatim
*>
*> \param[in] LDA
*> \verbatim
*> LDA is INTEGER
*> The leading dimension of the array A. LDA >= max(1,N).
*> \endverbatim
*>
*> \param[in,out] B
*> \verbatim
*> 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.
*> \endverbatim
*>
*> \param[in] LDB
*> \verbatim
*> LDB is INTEGER
*> The leading dimension of the array B. LDB >= max(1,N).
*> \endverbatim
*>
*> \param[out] W
*> \verbatim
*> W is DOUBLE PRECISION array, dimension (N)
*> If INFO = 0, the eigenvalues in ascending order.
*> \endverbatim
*>
*> \param[out] WORK
*> \verbatim
*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
*> \endverbatim
*>
*> \param[in] LWORK
*> \verbatim
*> LWORK is INTEGER
*> The dimension of the array WORK.
*> If N <= 1, LWORK >= 1.
*> If JOBZ = 'N' and N > 1, LWORK >= 2*N+1.
*> If JOBZ = 'V' and N > 1, LWORK >= 1 + 6*N + 2*N**2.
*>
*> If LWORK = -1, then a workspace query is assumed; the routine
*> only calculates the optimal sizes of the WORK and IWORK
*> arrays, returns these values as the first entries of the WORK
*> and IWORK arrays, and no error message related to LWORK or
*> LIWORK is issued by XERBLA.
*> \endverbatim
*>
*> \param[out] IWORK
*> \verbatim
*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
*> \endverbatim
*>
*> \param[in] LIWORK
*> \verbatim
*> LIWORK is INTEGER
*> The dimension of the array IWORK.
*> If N <= 1, LIWORK >= 1.
*> If JOBZ = 'N' and N > 1, LIWORK >= 1.
*> If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.
*>
*> If LIWORK = -1, then a workspace query is assumed; the
*> routine only calculates the optimal sizes of the WORK and
*> IWORK arrays, returns these values as the first entries of
*> the WORK and IWORK arrays, and no error message related to
*> LWORK or LIWORK is issued by XERBLA.
*> \endverbatim
*>
*> \param[out] INFO
*> \verbatim
*> INFO is INTEGER
*> = 0: successful exit
*> < 0: if INFO = -i, the i-th argument had an illegal value
*> > 0: DPOTRF or DSYEVD returned an error code:
*> <= N: if INFO = i and JOBZ = 'N', then the algorithm
*> failed to converge; i off-diagonal elements of an
*> intermediate tridiagonal form did not converge to
*> zero;
*> if INFO = i and JOBZ = 'V', then the algorithm
*> failed to compute an eigenvalue while working on
*> the submatrix lying in rows and columns INFO/(N+1)
*> through mod(INFO,N+1);
*> > N: if INFO = N + i, for 1 <= i <= N, then the leading
*> minor of order i of B is not positive definite.
*> The factorization of B could not be completed and
*> no eigenvalues or eigenvectors were computed.
*> \endverbatim
*
* Authors:
* ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \date November 2011
*
*> \ingroup doubleSYeigen
*
*> \par Further Details:
* =====================
*>
*> \verbatim
*>
*> Modified so that no backsubstitution is performed if DSYEVD fails to
*> converge (NEIG in old code could be greater than N causing out of
*> bounds reference to A - reported by Ralf Meyer). Also corrected the
*> description of INFO and the test on ITYPE. Sven, 16 Feb 05.
*> \endverbatim
*
*> \par Contributors:
* ==================
*>
*> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
*>
* =====================================================================
SUBROUTINE DSYGVD( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
$ LWORK, IWORK, LIWORK, INFO )
*
* -- LAPACK driver routine (version 3.4.0) --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
* November 2011
*
* .. Scalar Arguments ..
CHARACTER JOBZ, UPLO
INTEGER INFO, ITYPE, LDA, LDB, LIWORK, LWORK, N
* ..
* .. Array Arguments ..
INTEGER IWORK( * )
DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
* ..
*
* =====================================================================
*
* .. Parameters ..
DOUBLE PRECISION ONE
PARAMETER ( ONE = 1.0D+0 )
* ..
* .. Local Scalars ..
LOGICAL LQUERY, UPPER, WANTZ
CHARACTER TRANS
INTEGER LIOPT, LIWMIN, LOPT, LWMIN
* ..
* .. External Functions ..
LOGICAL LSAME
EXTERNAL LSAME
* ..
* .. External Subroutines ..
EXTERNAL DPOTRF, DSYEVD, DSYGST, DTRMM, DTRSM, XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC DBLE, MAX
* ..
* .. Executable Statements ..
*
* Test the input parameters.
*
WANTZ = LSAME( JOBZ, 'V' )
UPPER = LSAME( UPLO, 'U' )
LQUERY = ( LWORK.EQ.-1 .OR. LIWORK.EQ.-1 )
*
INFO = 0
IF( N.LE.1 ) THEN
LIWMIN = 1
LWMIN = 1
ELSE IF( WANTZ ) THEN
LIWMIN = 3 + 5*N
LWMIN = 1 + 6*N + 2*N**2
ELSE
LIWMIN = 1
LWMIN = 2*N + 1
END IF
LOPT = LWMIN
LIOPT = LIWMIN
IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
INFO = -1
ELSE IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
INFO = -2
ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
INFO = -3
ELSE IF( N.LT.0 ) THEN
INFO = -4
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
INFO = -6
ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
INFO = -8
END IF
*
IF( INFO.EQ.0 ) THEN
WORK( 1 ) = LOPT
IWORK( 1 ) = LIOPT
*
IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
INFO = -11
ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN
INFO = -13
END IF
END IF
*
IF( INFO.NE.0 ) THEN
CALL XERBLA( 'DSYGVD', -INFO )
RETURN
ELSE IF( LQUERY ) THEN
RETURN
END IF
*
* Quick return if possible
*
IF( N.EQ.0 )
$ RETURN
*
* Form a Cholesky factorization of B.
*
CALL DPOTRF( UPLO, N, B, LDB, INFO )
IF( INFO.NE.0 ) THEN
INFO = N + INFO
RETURN
END IF
*
* Transform problem to standard eigenvalue problem and solve.
*
CALL DSYGST( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
CALL DSYEVD( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK,
$ INFO )
LOPT = MAX( DBLE( LOPT ), DBLE( WORK( 1 ) ) )
LIOPT = MAX( DBLE( LIOPT ), DBLE( IWORK( 1 ) ) )
*
IF( WANTZ .AND. INFO.EQ.0 ) THEN
*
* Backtransform eigenvectors to the original problem.
*
IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
*
* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
*
IF( UPPER ) THEN
TRANS = 'N'
ELSE
TRANS = 'T'
END IF
*
CALL DTRSM( 'Left', UPLO, TRANS, 'Non-unit', N, N, ONE,
$ B, LDB, A, LDA )
*
ELSE IF( ITYPE.EQ.3 ) THEN
*
* For B*A*x=(lambda)*x;
* backtransform eigenvectors: x = L*y or U**T*y
*
IF( UPPER ) THEN
TRANS = 'T'
ELSE
TRANS = 'N'
END IF
*
CALL DTRMM( 'Left', UPLO, TRANS, 'Non-unit', N, N, ONE,
$ B, LDB, A, LDA )
END IF
END IF
*
WORK( 1 ) = LOPT
IWORK( 1 ) = LIOPT
*
RETURN
*
* End of DSYGVD
*
END

Event Timeline