LAPACK 3.12.1
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
!>                    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 177 of file chegv.f.

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