OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
zhbgvx.f
Go to the documentation of this file.
1*> \brief \b ZHBGVX
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZHBGVX + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhbgvx.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhbgvx.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhbgvx.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZHBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
22* LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
23* LDZ, WORK, RWORK, IWORK, IFAIL, INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER JOBZ, RANGE, UPLO
27* INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
28* $ N
29* DOUBLE PRECISION ABSTOL, VL, VU
30* ..
31* .. Array Arguments ..
32* INTEGER IFAIL( * ), IWORK( * )
33* DOUBLE PRECISION RWORK( * ), W( * )
34* COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
35* $ WORK( * ), Z( LDZ, * )
36* ..
37*
38*
39*> \par Purpose:
40* =============
41*>
42*> \verbatim
43*>
44*> ZHBGVX computes all the eigenvalues, and optionally, the eigenvectors
45*> of a complex generalized Hermitian-definite banded eigenproblem, of
46*> the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
47*> and banded, and B is also positive definite. Eigenvalues and
48*> eigenvectors can be selected by specifying either all eigenvalues,
49*> a range of values or a range of indices for the desired eigenvalues.
50*> \endverbatim
51*
52* Arguments:
53* ==========
54*
55*> \param[in] JOBZ
56*> \verbatim
57*> JOBZ is CHARACTER*1
58*> = 'N': Compute eigenvalues only;
59*> = 'V': Compute eigenvalues and eigenvectors.
60*> \endverbatim
61*>
62*> \param[in] RANGE
63*> \verbatim
64*> RANGE is CHARACTER*1
65*> = 'A': all eigenvalues will be found;
66*> = 'V': all eigenvalues in the half-open interval (VL,VU]
67*> will be found;
68*> = 'I': the IL-th through IU-th eigenvalues will be found.
69*> \endverbatim
70*>
71*> \param[in] UPLO
72*> \verbatim
73*> UPLO is CHARACTER*1
74*> = 'U': Upper triangles of A and B are stored;
75*> = 'L': Lower triangles of A and B are stored.
76*> \endverbatim
77*>
78*> \param[in] N
79*> \verbatim
80*> N is INTEGER
81*> The order of the matrices A and B. N >= 0.
82*> \endverbatim
83*>
84*> \param[in] KA
85*> \verbatim
86*> KA is INTEGER
87*> The number of superdiagonals of the matrix A if UPLO = 'U',
88*> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
89*> \endverbatim
90*>
91*> \param[in] KB
92*> \verbatim
93*> KB is INTEGER
94*> The number of superdiagonals of the matrix B if UPLO = 'U',
95*> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
96*> \endverbatim
97*>
98*> \param[in,out] AB
99*> \verbatim
100*> AB is COMPLEX*16 array, dimension (LDAB, N)
101*> On entry, the upper or lower triangle of the Hermitian band
102*> matrix A, stored in the first ka+1 rows of the array. The
103*> j-th column of A is stored in the j-th column of the array AB
104*> as follows:
105*> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
106*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
107*>
108*> On exit, the contents of AB are destroyed.
109*> \endverbatim
110*>
111*> \param[in] LDAB
112*> \verbatim
113*> LDAB is INTEGER
114*> The leading dimension of the array AB. LDAB >= KA+1.
115*> \endverbatim
116*>
117*> \param[in,out] BB
118*> \verbatim
119*> BB is COMPLEX*16 array, dimension (LDBB, N)
120*> On entry, the upper or lower triangle of the Hermitian band
121*> matrix B, stored in the first kb+1 rows of the array. The
122*> j-th column of B is stored in the j-th column of the array BB
123*> as follows:
124*> if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
125*> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
126*>
127*> On exit, the factor S from the split Cholesky factorization
128*> B = S**H*S, as returned by ZPBSTF.
129*> \endverbatim
130*>
131*> \param[in] LDBB
132*> \verbatim
133*> LDBB is INTEGER
134*> The leading dimension of the array BB. LDBB >= KB+1.
135*> \endverbatim
136*>
137*> \param[out] Q
138*> \verbatim
139*> Q is COMPLEX*16 array, dimension (LDQ, N)
140*> If JOBZ = 'V', the n-by-n matrix used in the reduction of
141*> A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
142*> and consequently C to tridiagonal form.
143*> If JOBZ = 'N', the array Q is not referenced.
144*> \endverbatim
145*>
146*> \param[in] LDQ
147*> \verbatim
148*> LDQ is INTEGER
149*> The leading dimension of the array Q. If JOBZ = 'N',
150*> LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).
151*> \endverbatim
152*>
153*> \param[in] VL
154*> \verbatim
155*> VL is DOUBLE PRECISION
156*>
157*> If RANGE='V', the lower bound of the interval to
158*> be searched for eigenvalues. VL < VU.
159*> Not referenced if RANGE = 'A' or 'I'.
160*> \endverbatim
161*>
162*> \param[in] VU
163*> \verbatim
164*> VU is DOUBLE PRECISION
165*>
166*> If RANGE='V', the upper bound of the interval to
167*> be searched for eigenvalues. VL < VU.
168*> Not referenced if RANGE = 'A' or 'I'.
169*> \endverbatim
170*>
171*> \param[in] IL
172*> \verbatim
173*> IL is INTEGER
174*>
175*> If RANGE='I', the index of the
176*> smallest eigenvalue to be returned.
177*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
178*> Not referenced if RANGE = 'A' or 'V'.
179*> \endverbatim
180*>
181*> \param[in] IU
182*> \verbatim
183*> IU is INTEGER
184*>
185*> If RANGE='I', the index of the
186*> largest eigenvalue to be returned.
187*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
188*> Not referenced if RANGE = 'A' or 'V'.
189*> \endverbatim
190*>
191*> \param[in] ABSTOL
192*> \verbatim
193*> ABSTOL is DOUBLE PRECISION
194*> The absolute error tolerance for the eigenvalues.
195*> An approximate eigenvalue is accepted as converged
196*> when it is determined to lie in an interval [a,b]
197*> of width less than or equal to
198*>
199*> ABSTOL + EPS * max( |a|,|b| ) ,
200*>
201*> where EPS is the machine precision. If ABSTOL is less than
202*> or equal to zero, then EPS*|T| will be used in its place,
203*> where |T| is the 1-norm of the tridiagonal matrix obtained
204*> by reducing AP to tridiagonal form.
205*>
206*> Eigenvalues will be computed most accurately when ABSTOL is
207*> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
208*> If this routine returns with INFO>0, indicating that some
209*> eigenvectors did not converge, try setting ABSTOL to
210*> 2*DLAMCH('S').
211*> \endverbatim
212*>
213*> \param[out] M
214*> \verbatim
215*> M is INTEGER
216*> The total number of eigenvalues found. 0 <= M <= N.
217*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
218*> \endverbatim
219*>
220*> \param[out] W
221*> \verbatim
222*> W is DOUBLE PRECISION array, dimension (N)
223*> If INFO = 0, the eigenvalues in ascending order.
224*> \endverbatim
225*>
226*> \param[out] Z
227*> \verbatim
228*> Z is COMPLEX*16 array, dimension (LDZ, N)
229*> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
230*> eigenvectors, with the i-th column of Z holding the
231*> eigenvector associated with W(i). The eigenvectors are
232*> normalized so that Z**H*B*Z = I.
233*> If JOBZ = 'N', then Z is not referenced.
234*> \endverbatim
235*>
236*> \param[in] LDZ
237*> \verbatim
238*> LDZ is INTEGER
239*> The leading dimension of the array Z. LDZ >= 1, and if
240*> JOBZ = 'V', LDZ >= N.
241*> \endverbatim
242*>
243*> \param[out] WORK
244*> \verbatim
245*> WORK is COMPLEX*16 array, dimension (N)
246*> \endverbatim
247*>
248*> \param[out] RWORK
249*> \verbatim
250*> RWORK is DOUBLE PRECISION array, dimension (7*N)
251*> \endverbatim
252*>
253*> \param[out] IWORK
254*> \verbatim
255*> IWORK is INTEGER array, dimension (5*N)
256*> \endverbatim
257*>
258*> \param[out] IFAIL
259*> \verbatim
260*> IFAIL is INTEGER array, dimension (N)
261*> If JOBZ = 'V', then if INFO = 0, the first M elements of
262*> IFAIL are zero. If INFO > 0, then IFAIL contains the
263*> indices of the eigenvectors that failed to converge.
264*> If JOBZ = 'N', then IFAIL is not referenced.
265*> \endverbatim
266*>
267*> \param[out] INFO
268*> \verbatim
269*> INFO is INTEGER
270*> = 0: successful exit
271*> < 0: if INFO = -i, the i-th argument had an illegal value
272*> > 0: if INFO = i, and i is:
273*> <= N: then i eigenvectors failed to converge. Their
274*> indices are stored in array IFAIL.
275*> > N: if INFO = N + i, for 1 <= i <= N, then ZPBSTF
276*> returned INFO = i: B is not positive definite.
277*> The factorization of B could not be completed and
278*> no eigenvalues or eigenvectors were computed.
279*> \endverbatim
280*
281* Authors:
282* ========
283*
284*> \author Univ. of Tennessee
285*> \author Univ. of California Berkeley
286*> \author Univ. of Colorado Denver
287*> \author NAG Ltd.
288*
289*> \ingroup complex16OTHEReigen
290*
291*> \par Contributors:
292* ==================
293*>
294*> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
295*
296* =====================================================================
297 SUBROUTINE zhbgvx( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
298 $ LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
299 $ LDZ, WORK, RWORK, IWORK, IFAIL, INFO )
300*
301* -- LAPACK driver routine --
302* -- LAPACK is a software package provided by Univ. of Tennessee, --
303* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
304*
305* .. Scalar Arguments ..
306 CHARACTER JOBZ, RANGE, UPLO
307 INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
308 $ n
309 DOUBLE PRECISION ABSTOL, VL, VU
310* ..
311* .. Array Arguments ..
312 INTEGER IFAIL( * ), IWORK( * )
313 DOUBLE PRECISION RWORK( * ), W( * )
314 COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
315 $ work( * ), z( ldz, * )
316* ..
317*
318* =====================================================================
319*
320* .. Parameters ..
321 DOUBLE PRECISION ZERO
322 PARAMETER ( ZERO = 0.0d+0 )
323 COMPLEX*16 CZERO, CONE
324 parameter( czero = ( 0.0d+0, 0.0d+0 ),
325 $ cone = ( 1.0d+0, 0.0d+0 ) )
326* ..
327* .. Local Scalars ..
328 LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
329 CHARACTER ORDER, VECT
330 INTEGER I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP,
331 $ indiwk, indrwk, indwrk, itmp1, j, jj, nsplit
332 DOUBLE PRECISION TMP1
333* ..
334* .. External Functions ..
335 LOGICAL LSAME
336 EXTERNAL LSAME
337* ..
338* .. External Subroutines ..
339 EXTERNAL dcopy, dstebz, dsterf, xerbla, zcopy, zgemv,
341 $ zswap
342* ..
343* .. Intrinsic Functions ..
344 INTRINSIC min
345* ..
346* .. Executable Statements ..
347*
348* Test the input parameters.
349*
350 wantz = lsame( jobz, 'V' )
351 upper = lsame( uplo, 'U' )
352 alleig = lsame( range, 'A' )
353 valeig = lsame( range, 'V' )
354 indeig = lsame( range, 'I' )
355*
356 info = 0
357 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
358 info = -1
359 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
360 info = -2
361 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
362 info = -3
363 ELSE IF( n.LT.0 ) THEN
364 info = -4
365 ELSE IF( ka.LT.0 ) THEN
366 info = -5
367 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
368 info = -6
369 ELSE IF( ldab.LT.ka+1 ) THEN
370 info = -8
371 ELSE IF( ldbb.LT.kb+1 ) THEN
372 info = -10
373 ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
374 info = -12
375 ELSE
376 IF( valeig ) THEN
377 IF( n.GT.0 .AND. vu.LE.vl )
378 $ info = -14
379 ELSE IF( indeig ) THEN
380 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
381 info = -15
382 ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
383 info = -16
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 = -21
390 END IF
391 END IF
392*
393 IF( info.NE.0 ) THEN
394 CALL xerbla( 'ZHBGVX', -info )
395 RETURN
396 END IF
397*
398* Quick return if possible
399*
400 m = 0
401 IF( n.EQ.0 )
402 $ RETURN
403*
404* Form a split Cholesky factorization of B.
405*
406 CALL zpbstf( uplo, n, kb, bb, ldbb, info )
407 IF( info.NE.0 ) THEN
408 info = n + info
409 RETURN
410 END IF
411*
412* Transform problem to standard eigenvalue problem.
413*
414 CALL zhbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
415 $ work, rwork, iinfo )
416*
417* Solve the standard eigenvalue problem.
418* Reduce Hermitian band matrix to tridiagonal form.
419*
420 indd = 1
421 inde = indd + n
422 indrwk = inde + n
423 indwrk = 1
424 IF( wantz ) THEN
425 vect = 'u'
426 ELSE
427 VECT = 'n'
428 END IF
429 CALL ZHBTRD( VECT, UPLO, N, KA, AB, LDAB, RWORK( INDD ),
430 $ RWORK( INDE ), Q, LDQ, WORK( INDWRK ), IINFO )
431*
432* If all eigenvalues are desired and ABSTOL is less than or equal
433* to zero, then call DSTERF or ZSTEQR. If this fails for some
434* eigenvalue, then try DSTEBZ.
435*
436 TEST = .FALSE.
437 IF( INDEIG ) THEN
438.EQ..AND..EQ. IF( IL1 IUN ) THEN
439 TEST = .TRUE.
440 END IF
441 END IF
442.OR..AND..LE. IF( ( ALLEIG TEST ) ( ABSTOLZERO ) ) THEN
443 CALL DCOPY( N, RWORK( INDD ), 1, W, 1 )
444 INDEE = INDRWK + 2*N
445 CALL DCOPY( N-1, RWORK( INDE ), 1, RWORK( INDEE ), 1 )
446.NOT. IF( WANTZ ) THEN
447 CALL DSTERF( N, W, RWORK( INDEE ), INFO )
448 ELSE
449 CALL ZLACPY( 'a', N, N, Q, LDQ, Z, LDZ )
450 CALL ZSTEQR( JOBZ, N, W, RWORK( INDEE ), Z, LDZ,
451 $ RWORK( INDRWK ), INFO )
452.EQ. IF( INFO0 ) THEN
453 DO 10 I = 1, N
454 IFAIL( I ) = 0
455 10 CONTINUE
456 END IF
457 END IF
458.EQ. IF( INFO0 ) THEN
459 M = N
460 GO TO 30
461 END IF
462 INFO = 0
463 END IF
464*
465* Otherwise, call DSTEBZ and, if eigenvectors are desired,
466* call ZSTEIN.
467*
468 IF( WANTZ ) THEN
469 ORDER = 'b'
470 ELSE
471 ORDER = 'e'
472 END IF
473 INDIBL = 1
474 INDISP = INDIBL + N
475 INDIWK = INDISP + N
476 CALL DSTEBZ( RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL,
477 $ RWORK( INDD ), RWORK( INDE ), M, NSPLIT, W,
478 $ IWORK( INDIBL ), IWORK( INDISP ), RWORK( INDRWK ),
479 $ IWORK( INDIWK ), INFO )
480*
481 IF( WANTZ ) THEN
482 CALL ZSTEIN( N, RWORK( INDD ), RWORK( INDE ), M, W,
483 $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
484 $ RWORK( INDRWK ), IWORK( INDIWK ), IFAIL, INFO )
485*
486* Apply unitary matrix used in reduction to tridiagonal
487* form to eigenvectors returned by ZSTEIN.
488*
489 DO 20 J = 1, M
490 CALL ZCOPY( N, Z( 1, J ), 1, WORK( 1 ), 1 )
491 CALL ZGEMV( 'n', N, N, CONE, Q, LDQ, WORK, 1, CZERO,
492 $ Z( 1, J ), 1 )
493 20 CONTINUE
494 END IF
495*
496 30 CONTINUE
497*
498* If eigenvalues are not in order, then sort them, along with
499* eigenvectors.
500*
501 IF( WANTZ ) THEN
502 DO 50 J = 1, M - 1
503 I = 0
504 TMP1 = W( J )
505 DO 40 JJ = J + 1, M
506.LT. IF( W( JJ )TMP1 ) THEN
507 I = JJ
508 TMP1 = W( JJ )
509 END IF
510 40 CONTINUE
511*
512.NE. IF( I0 ) THEN
513 ITMP1 = IWORK( INDIBL+I-1 )
514 W( I ) = W( J )
515 IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
516 W( J ) = TMP1
517 IWORK( INDIBL+J-1 ) = ITMP1
518 CALL ZSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
519.NE. IF( INFO0 ) THEN
520 ITMP1 = IFAIL( I )
521 IFAIL( I ) = IFAIL( J )
522 IFAIL( J ) = ITMP1
523 END IF
524 END IF
525 50 CONTINUE
526 END IF
527*
528 RETURN
529*
530* End of ZHBGVX
531*
532 END
subroutine dstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
DSTEBZ
Definition dstebz.f:273
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zhbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
ZHBGST
Definition zhbgst.f:165
subroutine zstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
ZSTEIN
Definition zstein.f:182
subroutine zpbstf(uplo, n, kd, ab, ldab, info)
ZPBSTF
Definition zpbstf.f:153
subroutine zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
Definition zsteqr.f:132
subroutine zhbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
ZHBTRD
Definition zhbtrd.f:163
subroutine zhbgvx(jobz, range, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
ZHBGVX
Definition zhbgvx.f:300
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:158
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21