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

◆ dsygv()

subroutine dsygv ( integer  itype,
character  jobz,
character  uplo,
integer  n,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( ldb, * )  b,
integer  ldb,
double precision, dimension( * )  w,
double precision, dimension( * )  work,
integer  lwork,
integer  info 
)

DSYGV

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

Purpose:
 DSYGV 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.
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]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and 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, 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.
[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 positive definite 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).
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[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,3*N-1).
          For optimal efficiency, LWORK >= (NB+2)*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]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  DPOTRF or DSYEV returned an error code:
             <= N:  if INFO = i, DSYEV failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not converge to zero;
             > 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.

Definition at line 173 of file dsygv.f.

175*
176* -- LAPACK driver routine --
177* -- LAPACK is a software package provided by Univ. of Tennessee, --
178* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
179*
180* .. Scalar Arguments ..
181 CHARACTER JOBZ, UPLO
182 INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
183* ..
184* .. Array Arguments ..
185 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
186* ..
187*
188* =====================================================================
189*
190* .. Parameters ..
191 DOUBLE PRECISION ONE
192 parameter( one = 1.0d+0 )
193* ..
194* .. Local Scalars ..
195 LOGICAL LQUERY, UPPER, WANTZ
196 CHARACTER TRANS
197 INTEGER LWKMIN, LWKOPT, NB, NEIG
198* ..
199* .. External Functions ..
200 LOGICAL LSAME
201 INTEGER ILAENV
202 EXTERNAL lsame, ilaenv
203* ..
204* .. External Subroutines ..
205 EXTERNAL dpotrf, dsyev, dsygst, dtrmm, dtrsm, xerbla
206* ..
207* .. Intrinsic Functions ..
208 INTRINSIC max
209* ..
210* .. Executable Statements ..
211*
212* Test the input parameters.
213*
214 wantz = lsame( jobz, 'V' )
215 upper = lsame( uplo, 'U' )
216 lquery = ( lwork.EQ.-1 )
217*
218 info = 0
219 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
220 info = -1
221 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
222 info = -2
223 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
224 info = -3
225 ELSE IF( n.LT.0 ) THEN
226 info = -4
227 ELSE IF( lda.LT.max( 1, n ) ) THEN
228 info = -6
229 ELSE IF( ldb.LT.max( 1, n ) ) THEN
230 info = -8
231 END IF
232*
233 IF( info.EQ.0 ) THEN
234 lwkmin = max( 1, 3*n - 1 )
235 nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
236 lwkopt = max( lwkmin, ( nb + 2 )*n )
237 work( 1 ) = lwkopt
238*
239 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
240 info = -11
241 END IF
242 END IF
243*
244 IF( info.NE.0 ) THEN
245 CALL xerbla( 'DSYGV ', -info )
246 RETURN
247 ELSE IF( lquery ) THEN
248 RETURN
249 END IF
250*
251* Quick return if possible
252*
253 IF( n.EQ.0 )
254 $ RETURN
255*
256* Form a Cholesky factorization of B.
257*
258 CALL dpotrf( uplo, n, b, ldb, info )
259 IF( info.NE.0 ) THEN
260 info = n + info
261 RETURN
262 END IF
263*
264* Transform problem to standard eigenvalue problem and solve.
265*
266 CALL dsygst( itype, uplo, n, a, lda, b, ldb, info )
267 CALL dsyev( jobz, uplo, n, a, lda, w, work, lwork, info )
268*
269 IF( wantz ) THEN
270*
271* Backtransform eigenvectors to the original problem.
272*
273 neig = n
274 IF( info.GT.0 )
275 $ neig = info - 1
276 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
277*
278* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
279* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
280*
281 IF( upper ) THEN
282 trans = 'N'
283 ELSE
284 trans = 'T'
285 END IF
286*
287 CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
288 $ b, ldb, a, lda )
289*
290 ELSE IF( itype.EQ.3 ) THEN
291*
292* For B*A*x=(lambda)*x;
293* backtransform eigenvectors: x = L*y or U**T*y
294*
295 IF( upper ) THEN
296 trans = 'T'
297 ELSE
298 trans = 'N'
299 END IF
300*
301 CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
302 $ b, ldb, a, lda )
303 END IF
304 END IF
305*
306 work( 1 ) = lwkopt
307 RETURN
308*
309* End of DSYGV
310*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsyev(jobz, uplo, n, a, lda, w, work, lwork, info)
DSYEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Definition dsyev.f:132
subroutine dsygst(itype, uplo, n, a, lda, b, ldb, info)
DSYGST
Definition dsygst.f:127
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:107
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: