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

◆ chegv()

subroutine chegv ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  W,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CHEGV

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

Purpose:
 CHEGV computes all the eigenvalues, and optionally, the eigenvectors
 of a complex generalized Hermitian-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 Hermitian 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 COMPLEX array, dimension (LDA, N)
          On entry, the Hermitian 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**H*B*Z = I;
          if ITYPE = 3, Z**H*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 COMPLEX array, dimension (LDB, N)
          On entry, the Hermitian 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**H*U or B = L*L**H.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is COMPLEX 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,2*N-1).
          For optimal efficiency, LWORK >= (NB+1)*N,
          where NB is the blocksize for CHETRD 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]RWORK
          RWORK is REAL array, dimension (max(1, 3*N-2))
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  CPOTRF or CHEEV returned an error code:
             <= N:  if INFO = i, CHEEV 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
                    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.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 179 of file chegv.f.

181*
182* -- LAPACK driver routine --
183* -- LAPACK is a software package provided by Univ. of Tennessee, --
184* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
185*
186* .. Scalar Arguments ..
187 CHARACTER JOBZ, UPLO
188 INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
189* ..
190* .. Array Arguments ..
191 REAL RWORK( * ), W( * )
192 COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
193* ..
194*
195* =====================================================================
196*
197* .. Parameters ..
198 COMPLEX ONE
199 parameter( one = ( 1.0e+0, 0.0e+0 ) )
200* ..
201* .. Local Scalars ..
202 LOGICAL LQUERY, UPPER, WANTZ
203 CHARACTER TRANS
204 INTEGER LWKOPT, NB, NEIG
205* ..
206* .. External Functions ..
207 LOGICAL LSAME
208 INTEGER ILAENV
209 EXTERNAL ilaenv, lsame
210* ..
211* .. External Subroutines ..
212 EXTERNAL cheev, chegst, cpotrf, ctrmm, ctrsm, xerbla
213* ..
214* .. Intrinsic Functions ..
215 INTRINSIC max
216* ..
217* .. Executable Statements ..
218*
219* Test the input parameters.
220*
221 wantz = lsame( jobz, 'V' )
222 upper = lsame( uplo, 'U' )
223 lquery = ( lwork.EQ. -1 )
224*
225 info = 0
226 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
227 info = -1
228 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
229 info = -2
230 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
231 info = -3
232 ELSE IF( n.LT.0 ) THEN
233 info = -4
234 ELSE IF( lda.LT.max( 1, n ) ) THEN
235 info = -6
236 ELSE IF( ldb.LT.max( 1, n ) ) THEN
237 info = -8
238 END IF
239*
240 IF( info.EQ.0 ) THEN
241 nb = ilaenv( 1, 'CHETRD', uplo, n, -1, -1, -1 )
242 lwkopt = max( 1, ( nb + 1 )*n )
243 work( 1 ) = lwkopt
244*
245 IF( lwork.LT.max( 1, 2*n-1 ) .AND. .NOT.lquery ) THEN
246 info = -11
247 END IF
248 END IF
249*
250 IF( info.NE.0 ) THEN
251 CALL xerbla( 'CHEGV ', -info )
252 RETURN
253 ELSE IF( lquery ) THEN
254 RETURN
255 END IF
256*
257* Quick return if possible
258*
259 IF( n.EQ.0 )
260 $ RETURN
261*
262* Form a Cholesky factorization of B.
263*
264 CALL cpotrf( uplo, n, b, ldb, info )
265 IF( info.NE.0 ) THEN
266 info = n + info
267 RETURN
268 END IF
269*
270* Transform problem to standard eigenvalue problem and solve.
271*
272 CALL chegst( itype, uplo, n, a, lda, b, ldb, info )
273 CALL cheev( jobz, uplo, n, a, lda, w, work, lwork, rwork, info )
274*
275 IF( wantz ) THEN
276*
277* Backtransform eigenvectors to the original problem.
278*
279 neig = n
280 IF( info.GT.0 )
281 $ neig = info - 1
282 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
283*
284* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
285* backtransform eigenvectors: x = inv(L)**H*y or inv(U)*y
286*
287 IF( upper ) THEN
288 trans = 'N'
289 ELSE
290 trans = 'C'
291 END IF
292*
293 CALL ctrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
294 $ b, ldb, a, lda )
295*
296 ELSE IF( itype.EQ.3 ) THEN
297*
298* For B*A*x=(lambda)*x;
299* backtransform eigenvectors: x = L*y or U**H*y
300*
301 IF( upper ) THEN
302 trans = 'C'
303 ELSE
304 trans = 'N'
305 END IF
306*
307 CALL ctrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
308 $ b, ldb, a, lda )
309 END IF
310 END IF
311*
312 work( 1 ) = lwkopt
313*
314 RETURN
315*
316* End of CHEGV
317*
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine ctrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRMM
Definition: ctrmm.f:177
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:180
subroutine chegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
CHEGST
Definition: chegst.f:128
subroutine cheev(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, INFO)
CHEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition: cheev.f:140
subroutine cpotrf(UPLO, N, A, LDA, INFO)
CPOTRF
Definition: cpotrf.f:107
Here is the call graph for this function:
Here is the caller graph for this function: