LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
chegv.f
Go to the documentation of this file.
1*> \brief \b CHEGV
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download CHEGV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chegv.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chegv.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chegv.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE CHEGV( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
22* LWORK, RWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER JOBZ, UPLO
26* INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
27* ..
28* .. Array Arguments ..
29* REAL RWORK( * ), W( * )
30* COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> CHEGV computes all the eigenvalues, and optionally, the eigenvectors
40*> of a complex generalized Hermitian-definite eigenproblem, of the form
41*> A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
42*> Here A and B are assumed to be Hermitian and B is also
43*> positive definite.
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] ITYPE
50*> \verbatim
51*> ITYPE is INTEGER
52*> Specifies the problem type to be solved:
53*> = 1: A*x = (lambda)*B*x
54*> = 2: A*B*x = (lambda)*x
55*> = 3: B*A*x = (lambda)*x
56*> \endverbatim
57*>
58*> \param[in] JOBZ
59*> \verbatim
60*> JOBZ is CHARACTER*1
61*> = 'N': Compute eigenvalues only;
62*> = 'V': Compute eigenvalues and eigenvectors.
63*> \endverbatim
64*>
65*> \param[in] UPLO
66*> \verbatim
67*> UPLO is CHARACTER*1
68*> = 'U': Upper triangles of A and B are stored;
69*> = 'L': Lower triangles of A and B are stored.
70*> \endverbatim
71*>
72*> \param[in] N
73*> \verbatim
74*> N is INTEGER
75*> The order of the matrices A and B. N >= 0.
76*> \endverbatim
77*>
78*> \param[in,out] A
79*> \verbatim
80*> A is COMPLEX array, dimension (LDA, N)
81*> On entry, the Hermitian matrix A. If UPLO = 'U', the
82*> leading N-by-N upper triangular part of A contains the
83*> upper triangular part of the matrix A. If UPLO = 'L',
84*> the leading N-by-N lower triangular part of A contains
85*> the lower triangular part of the matrix A.
86*>
87*> On exit, if JOBZ = 'V', then if INFO = 0, A contains the
88*> matrix Z of eigenvectors. The eigenvectors are normalized
89*> as follows:
90*> if ITYPE = 1 or 2, Z**H*B*Z = I;
91*> if ITYPE = 3, Z**H*inv(B)*Z = I.
92*> If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
93*> or the lower triangle (if UPLO='L') of A, including the
94*> diagonal, is destroyed.
95*> \endverbatim
96*>
97*> \param[in] LDA
98*> \verbatim
99*> LDA is INTEGER
100*> The leading dimension of the array A. LDA >= max(1,N).
101*> \endverbatim
102*>
103*> \param[in,out] B
104*> \verbatim
105*> B is COMPLEX array, dimension (LDB, N)
106*> On entry, the Hermitian positive definite matrix B.
107*> If UPLO = 'U', the leading N-by-N upper triangular part of B
108*> contains the upper triangular part of the matrix B.
109*> If UPLO = 'L', the leading N-by-N lower triangular part of B
110*> contains the lower triangular part of the matrix B.
111*>
112*> On exit, if INFO <= N, the part of B containing the matrix is
113*> overwritten by the triangular factor U or L from the Cholesky
114*> factorization B = U**H*U or B = L*L**H.
115*> \endverbatim
116*>
117*> \param[in] LDB
118*> \verbatim
119*> LDB is INTEGER
120*> The leading dimension of the array B. LDB >= max(1,N).
121*> \endverbatim
122*>
123*> \param[out] W
124*> \verbatim
125*> W is REAL array, dimension (N)
126*> If INFO = 0, the eigenvalues in ascending order.
127*> \endverbatim
128*>
129*> \param[out] WORK
130*> \verbatim
131*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
132*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
133*> \endverbatim
134*>
135*> \param[in] LWORK
136*> \verbatim
137*> LWORK is INTEGER
138*> The length of the array WORK. LWORK >= max(1,2*N-1).
139*> For optimal efficiency, LWORK >= (NB+1)*N,
140*> where NB is the blocksize for CHETRD returned by ILAENV.
141*>
142*> If LWORK = -1, then a workspace query is assumed; the routine
143*> only calculates the optimal size of the WORK array, returns
144*> this value as the first entry of the WORK array, and no error
145*> message related to LWORK is issued by XERBLA.
146*> \endverbatim
147*>
148*> \param[out] RWORK
149*> \verbatim
150*> RWORK is REAL array, dimension (max(1, 3*N-2))
151*> \endverbatim
152*>
153*> \param[out] INFO
154*> \verbatim
155*> INFO is INTEGER
156*> = 0: successful exit
157*> < 0: if INFO = -i, the i-th argument had an illegal value
158*> > 0: CPOTRF or CHEEV returned an error code:
159*> <= N: if INFO = i, CHEEV failed to converge;
160*> i off-diagonal elements of an intermediate
161*> tridiagonal form did not converge to zero;
162*> > N: if INFO = N + i, for 1 <= i <= N, then the leading
163*> principal minor of order i of B is not positive.
164*> The factorization of B could not be completed and
165*> no eigenvalues or eigenvectors were computed.
166*> \endverbatim
167*
168* Authors:
169* ========
170*
171*> \author Univ. of Tennessee
172*> \author Univ. of California Berkeley
173*> \author Univ. of Colorado Denver
174*> \author NAG Ltd.
175*
176*> \ingroup hegv
177*
178* =====================================================================
179 SUBROUTINE chegv( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
180 $ LWORK, RWORK, INFO )
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 REAL SROUNDUP_LWORK
210 EXTERNAL ilaenv, lsame, sroundup_lwork
211* ..
212* .. External Subroutines ..
213 EXTERNAL cheev, chegst, cpotrf, ctrmm, ctrsm, 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, info )
275*
276 IF( wantz ) THEN
277*
278* Backtransform eigenvectors to the original problem.
279*
280 neig = n
281 IF( info.GT.0 )
282 $ neig = info - 1
283 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
284*
285* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
286* backtransform eigenvectors: x = inv(L)**H*y or inv(U)*y
287*
288 IF( upper ) THEN
289 trans = 'N'
290 ELSE
291 trans = 'C'
292 END IF
293*
294 CALL ctrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
295 $ b, ldb, a, lda )
296*
297 ELSE IF( itype.EQ.3 ) THEN
298*
299* For B*A*x=(lambda)*x;
300* backtransform eigenvectors: x = L*y or U**H*y
301*
302 IF( upper ) THEN
303 trans = 'C'
304 ELSE
305 trans = 'N'
306 END IF
307*
308 CALL ctrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
309 $ b, ldb, a, lda )
310 END IF
311 END IF
312*
313 work( 1 ) = sroundup_lwork(lwkopt)
314*
315 RETURN
316*
317* End of CHEGV
318*
319 END
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:140
subroutine chegst(itype, uplo, n, a, lda, b, ldb, info)
CHEGST
Definition chegst.f:128
subroutine chegv(itype, jobz, uplo, n, a, lda, b, ldb, w, work, lwork, rwork, info)
CHEGV
Definition chegv.f:181
subroutine cpotrf(uplo, n, a, lda, info)
CPOTRF
Definition cpotrf.f:107
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