LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zhegvx.f
Go to the documentation of this file.
1 *> \brief \b ZHEGVX
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZHEGVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhegvx.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhegvx.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhegvx.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZHEGVX( ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB,
22 * VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
23 * LWORK, RWORK, IWORK, IFAIL, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER JOBZ, RANGE, UPLO
27 * INTEGER IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N
28 * DOUBLE PRECISION ABSTOL, VL, VU
29 * ..
30 * .. Array Arguments ..
31 * INTEGER IFAIL( * ), IWORK( * )
32 * DOUBLE PRECISION RWORK( * ), W( * )
33 * COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * ),
34 * $ Z( LDZ, * )
35 * ..
36 *
37 *
38 *> \par Purpose:
39 * =============
40 *>
41 *> \verbatim
42 *>
43 *> ZHEGVX computes selected eigenvalues, and optionally, eigenvectors
44 *> of a complex generalized Hermitian-definite eigenproblem, of the form
45 *> A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. Here A and
46 *> B are assumed to be Hermitian and B is also positive definite.
47 *> Eigenvalues and eigenvectors can be selected by specifying either a
48 *> range of values or a range of indices for the desired eigenvalues.
49 *> \endverbatim
50 *
51 * Arguments:
52 * ==========
53 *
54 *> \param[in] ITYPE
55 *> \verbatim
56 *> ITYPE is INTEGER
57 *> Specifies the problem type to be solved:
58 *> = 1: A*x = (lambda)*B*x
59 *> = 2: A*B*x = (lambda)*x
60 *> = 3: B*A*x = (lambda)*x
61 *> \endverbatim
62 *>
63 *> \param[in] JOBZ
64 *> \verbatim
65 *> JOBZ is CHARACTER*1
66 *> = 'N': Compute eigenvalues only;
67 *> = 'V': Compute eigenvalues and eigenvectors.
68 *> \endverbatim
69 *>
70 *> \param[in] RANGE
71 *> \verbatim
72 *> RANGE is CHARACTER*1
73 *> = 'A': all eigenvalues will be found.
74 *> = 'V': all eigenvalues in the half-open interval (VL,VU]
75 *> will be found.
76 *> = 'I': the IL-th through IU-th eigenvalues will be found.
77 *> \endverbatim
78 *>
79 *> \param[in] UPLO
80 *> \verbatim
81 *> UPLO is CHARACTER*1
82 *> = 'U': Upper triangles of A and B are stored;
83 *> = 'L': Lower triangles of A and B are stored.
84 *> \endverbatim
85 *>
86 *> \param[in] N
87 *> \verbatim
88 *> N is INTEGER
89 *> The order of the matrices A and B. N >= 0.
90 *> \endverbatim
91 *>
92 *> \param[in,out] A
93 *> \verbatim
94 *> A is COMPLEX*16 array, dimension (LDA, N)
95 *> On entry, the Hermitian matrix A. If UPLO = 'U', the
96 *> leading N-by-N upper triangular part of A contains the
97 *> upper triangular part of the matrix A. If UPLO = 'L',
98 *> the leading N-by-N lower triangular part of A contains
99 *> the lower triangular part of the matrix A.
100 *>
101 *> On exit, the lower triangle (if UPLO='L') or the upper
102 *> triangle (if UPLO='U') of A, including the diagonal, is
103 *> destroyed.
104 *> \endverbatim
105 *>
106 *> \param[in] LDA
107 *> \verbatim
108 *> LDA is INTEGER
109 *> The leading dimension of the array A. LDA >= max(1,N).
110 *> \endverbatim
111 *>
112 *> \param[in,out] B
113 *> \verbatim
114 *> B is COMPLEX*16 array, dimension (LDB, N)
115 *> On entry, the Hermitian matrix B. If UPLO = 'U', the
116 *> leading N-by-N upper triangular part of B contains the
117 *> upper triangular part of the matrix B. If UPLO = 'L',
118 *> the leading N-by-N lower triangular part of B contains
119 *> the lower triangular part of the matrix B.
120 *>
121 *> On exit, if INFO <= N, the part of B containing the matrix is
122 *> overwritten by the triangular factor U or L from the Cholesky
123 *> factorization B = U**H*U or B = L*L**H.
124 *> \endverbatim
125 *>
126 *> \param[in] LDB
127 *> \verbatim
128 *> LDB is INTEGER
129 *> The leading dimension of the array B. LDB >= max(1,N).
130 *> \endverbatim
131 *>
132 *> \param[in] VL
133 *> \verbatim
134 *> VL is DOUBLE PRECISION
135 *>
136 *> If RANGE='V', the lower bound of the interval to
137 *> be searched for eigenvalues. VL < VU.
138 *> Not referenced if RANGE = 'A' or 'I'.
139 *> \endverbatim
140 *>
141 *> \param[in] VU
142 *> \verbatim
143 *> VU is DOUBLE PRECISION
144 *>
145 *> If RANGE='V', the upper bound of the interval to
146 *> be searched for eigenvalues. VL < VU.
147 *> Not referenced if RANGE = 'A' or 'I'.
148 *> \endverbatim
149 *>
150 *> \param[in] IL
151 *> \verbatim
152 *> IL is INTEGER
153 *>
154 *> If RANGE='I', the index of the
155 *> smallest eigenvalue to be returned.
156 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
157 *> Not referenced if RANGE = 'A' or 'V'.
158 *> \endverbatim
159 *>
160 *> \param[in] IU
161 *> \verbatim
162 *> IU is INTEGER
163 *>
164 *> If RANGE='I', the index of the
165 *> largest eigenvalue to be returned.
166 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
167 *> Not referenced if RANGE = 'A' or 'V'.
168 *> \endverbatim
169 *>
170 *> \param[in] ABSTOL
171 *> \verbatim
172 *> ABSTOL is DOUBLE PRECISION
173 *> The absolute error tolerance for the eigenvalues.
174 *> An approximate eigenvalue is accepted as converged
175 *> when it is determined to lie in an interval [a,b]
176 *> of width less than or equal to
177 *>
178 *> ABSTOL + EPS * max( |a|,|b| ) ,
179 *>
180 *> where EPS is the machine precision. If ABSTOL is less than
181 *> or equal to zero, then EPS*|T| will be used in its place,
182 *> where |T| is the 1-norm of the tridiagonal matrix obtained
183 *> by reducing C to tridiagonal form, where C is the symmetric
184 *> matrix of the standard symmetric problem to which the
185 *> generalized problem is transformed.
186 *>
187 *> Eigenvalues will be computed most accurately when ABSTOL is
188 *> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
189 *> If this routine returns with INFO>0, indicating that some
190 *> eigenvectors did not converge, try setting ABSTOL to
191 *> 2*DLAMCH('S').
192 *> \endverbatim
193 *>
194 *> \param[out] M
195 *> \verbatim
196 *> M is INTEGER
197 *> The total number of eigenvalues found. 0 <= M <= N.
198 *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
199 *> \endverbatim
200 *>
201 *> \param[out] W
202 *> \verbatim
203 *> W is DOUBLE PRECISION array, dimension (N)
204 *> The first M elements contain the selected
205 *> eigenvalues in ascending order.
206 *> \endverbatim
207 *>
208 *> \param[out] Z
209 *> \verbatim
210 *> Z is COMPLEX*16 array, dimension (LDZ, max(1,M))
211 *> If JOBZ = 'N', then Z is not referenced.
212 *> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
213 *> contain the orthonormal eigenvectors of the matrix A
214 *> corresponding to the selected eigenvalues, with the i-th
215 *> column of Z holding the eigenvector associated with W(i).
216 *> The eigenvectors are normalized as follows:
217 *> if ITYPE = 1 or 2, Z**T*B*Z = I;
218 *> if ITYPE = 3, Z**T*inv(B)*Z = I.
219 *>
220 *> If an eigenvector fails to converge, then that column of Z
221 *> contains the latest approximation to the eigenvector, and the
222 *> index of the eigenvector is returned in IFAIL.
223 *> Note: the user must ensure that at least max(1,M) columns are
224 *> supplied in the array Z; if RANGE = 'V', the exact value of M
225 *> is not known in advance and an upper bound must be used.
226 *> \endverbatim
227 *>
228 *> \param[in] LDZ
229 *> \verbatim
230 *> LDZ is INTEGER
231 *> The leading dimension of the array Z. LDZ >= 1, and if
232 *> JOBZ = 'V', LDZ >= max(1,N).
233 *> \endverbatim
234 *>
235 *> \param[out] WORK
236 *> \verbatim
237 *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
238 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
239 *> \endverbatim
240 *>
241 *> \param[in] LWORK
242 *> \verbatim
243 *> LWORK is INTEGER
244 *> The length of the array WORK. LWORK >= max(1,2*N).
245 *> For optimal efficiency, LWORK >= (NB+1)*N,
246 *> where NB is the blocksize for ZHETRD returned by ILAENV.
247 *>
248 *> If LWORK = -1, then a workspace query is assumed; the routine
249 *> only calculates the optimal size of the WORK array, returns
250 *> this value as the first entry of the WORK array, and no error
251 *> message related to LWORK is issued by XERBLA.
252 *> \endverbatim
253 *>
254 *> \param[out] RWORK
255 *> \verbatim
256 *> RWORK is DOUBLE PRECISION array, dimension (7*N)
257 *> \endverbatim
258 *>
259 *> \param[out] IWORK
260 *> \verbatim
261 *> IWORK is INTEGER array, dimension (5*N)
262 *> \endverbatim
263 *>
264 *> \param[out] IFAIL
265 *> \verbatim
266 *> IFAIL is INTEGER array, dimension (N)
267 *> If JOBZ = 'V', then if INFO = 0, the first M elements of
268 *> IFAIL are zero. If INFO > 0, then IFAIL contains the
269 *> indices of the eigenvectors that failed to converge.
270 *> If JOBZ = 'N', then IFAIL is not referenced.
271 *> \endverbatim
272 *>
273 *> \param[out] INFO
274 *> \verbatim
275 *> INFO is INTEGER
276 *> = 0: successful exit
277 *> < 0: if INFO = -i, the i-th argument had an illegal value
278 *> > 0: ZPOTRF or ZHEEVX returned an error code:
279 *> <= N: if INFO = i, ZHEEVX failed to converge;
280 *> i eigenvectors failed to converge. Their indices
281 *> are stored in array IFAIL.
282 *> > N: if INFO = N + i, for 1 <= i <= N, then the leading
283 *> minor of order i of B is not positive definite.
284 *> The factorization of B could not be completed and
285 *> no eigenvalues or eigenvectors were computed.
286 *> \endverbatim
287 *
288 * Authors:
289 * ========
290 *
291 *> \author Univ. of Tennessee
292 *> \author Univ. of California Berkeley
293 *> \author Univ. of Colorado Denver
294 *> \author NAG Ltd.
295 *
296 *> \date June 2016
297 *
298 *> \ingroup complex16HEeigen
299 *
300 *> \par Contributors:
301 * ==================
302 *>
303 *> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
304 *
305 * =====================================================================
306  SUBROUTINE zhegvx( ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB,
307  $ vl, vu, il, iu, abstol, m, w, z, ldz, work,
308  $ lwork, rwork, iwork, ifail, info )
309 *
310 * -- LAPACK driver routine (version 3.6.1) --
311 * -- LAPACK is a software package provided by Univ. of Tennessee, --
312 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
313 * June 2016
314 *
315 * .. Scalar Arguments ..
316  CHARACTER JOBZ, RANGE, UPLO
317  INTEGER IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N
318  DOUBLE PRECISION ABSTOL, VL, VU
319 * ..
320 * .. Array Arguments ..
321  INTEGER IFAIL( * ), IWORK( * )
322  DOUBLE PRECISION RWORK( * ), W( * )
323  COMPLEX*16 A( lda, * ), B( ldb, * ), WORK( * ),
324  $ z( ldz, * )
325 * ..
326 *
327 * =====================================================================
328 *
329 * .. Parameters ..
330  COMPLEX*16 CONE
331  parameter ( cone = ( 1.0d+0, 0.0d+0 ) )
332 * ..
333 * .. Local Scalars ..
334  LOGICAL ALLEIG, INDEIG, LQUERY, UPPER, VALEIG, WANTZ
335  CHARACTER TRANS
336  INTEGER LWKOPT, NB
337 * ..
338 * .. External Functions ..
339  LOGICAL LSAME
340  INTEGER ILAENV
341  EXTERNAL lsame, ilaenv
342 * ..
343 * .. External Subroutines ..
344  EXTERNAL xerbla, zheevx, zhegst, zpotrf, ztrmm, ztrsm
345 * ..
346 * .. Intrinsic Functions ..
347  INTRINSIC max, min
348 * ..
349 * .. Executable Statements ..
350 *
351 * Test the input parameters.
352 *
353  wantz = lsame( jobz, 'V' )
354  upper = lsame( uplo, 'U' )
355  alleig = lsame( range, 'A' )
356  valeig = lsame( range, 'V' )
357  indeig = lsame( range, 'I' )
358  lquery = ( lwork.EQ.-1 )
359 *
360  info = 0
361  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
362  info = -1
363  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
364  info = -2
365  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
366  info = -3
367  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
368  info = -4
369  ELSE IF( n.LT.0 ) THEN
370  info = -5
371  ELSE IF( lda.LT.max( 1, n ) ) THEN
372  info = -7
373  ELSE IF( ldb.LT.max( 1, n ) ) THEN
374  info = -9
375  ELSE
376  IF( valeig ) THEN
377  IF( n.GT.0 .AND. vu.LE.vl )
378  $ info = -11
379  ELSE IF( indeig ) THEN
380  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
381  info = -12
382  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
383  info = -13
384  END IF
385  END IF
386  END IF
387  IF (info.EQ.0) THEN
388  IF (ldz.LT.1 .OR. (wantz .AND. ldz.LT.n)) THEN
389  info = -18
390  END IF
391  END IF
392 *
393  IF( info.EQ.0 ) THEN
394  nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
395  lwkopt = max( 1, ( nb + 1 )*n )
396  work( 1 ) = lwkopt
397 *
398  IF( lwork.LT.max( 1, 2*n ) .AND. .NOT.lquery ) THEN
399  info = -20
400  END IF
401  END IF
402 *
403  IF( info.NE.0 ) THEN
404  CALL xerbla( 'ZHEGVX', -info )
405  RETURN
406  ELSE IF( lquery ) THEN
407  RETURN
408  END IF
409 *
410 * Quick return if possible
411 *
412  m = 0
413  IF( n.EQ.0 ) THEN
414  RETURN
415  END IF
416 *
417 * Form a Cholesky factorization of B.
418 *
419  CALL zpotrf( uplo, n, b, ldb, info )
420  IF( info.NE.0 ) THEN
421  info = n + info
422  RETURN
423  END IF
424 *
425 * Transform problem to standard eigenvalue problem and solve.
426 *
427  CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
428  CALL zheevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol,
429  $ m, w, z, ldz, work, lwork, rwork, iwork, ifail,
430  $ info )
431 *
432  IF( wantz ) THEN
433 *
434 * Backtransform eigenvectors to the original problem.
435 *
436  IF( info.GT.0 )
437  $ m = info - 1
438  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
439 *
440 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
441 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
442 *
443  IF( upper ) THEN
444  trans = 'N'
445  ELSE
446  trans = 'C'
447  END IF
448 *
449  CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, m, cone, b,
450  $ ldb, z, ldz )
451 *
452  ELSE IF( itype.EQ.3 ) THEN
453 *
454 * For B*A*x=(lambda)*x;
455 * backtransform eigenvectors: x = L*y or U**H *y
456 *
457  IF( upper ) THEN
458  trans = 'C'
459  ELSE
460  trans = 'N'
461  END IF
462 *
463  CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, m, cone, b,
464  $ ldb, z, ldz )
465  END IF
466  END IF
467 *
468 * Set WORK(1) to optimal complex workspace size.
469 *
470  work( 1 ) = lwkopt
471 *
472  RETURN
473 *
474 * End of ZHEGVX
475 *
476  END
subroutine zpotrf(UPLO, N, A, LDA, INFO)
ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
Definition: zpotrf.f:102
subroutine zhegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
ZHEGST
Definition: zhegst.f:129
subroutine zhegvx(ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
ZHEGVX
Definition: zhegvx.f:309
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ztrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRMM
Definition: ztrmm.f:179
subroutine zheevx(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
ZHEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: zheevx.f:261
subroutine ztrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRSM
Definition: ztrsm.f:182