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

◆ zhbgvd()

subroutine zhbgvd ( character  JOBZ,
character  UPLO,
integer  N,
integer  KA,
integer  KB,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
complex*16, dimension( ldbb, * )  BB,
integer  LDBB,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

ZHBGVD

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

Purpose:
 ZHBGVD computes all the eigenvalues, and optionally, the eigenvectors
 of a complex generalized Hermitian-definite banded eigenproblem, of
 the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
 and banded, and B is also positive definite.  If eigenvectors are
 desired, it uses a divide and conquer algorithm.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
Parameters
[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]KA
          KA is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'. KA >= 0.
[in]KB
          KB is INTEGER
          The number of superdiagonals of the matrix B if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'. KB >= 0.
[in,out]AB
          AB is COMPLEX*16 array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix A, stored in the first ka+1 rows of the array.  The
          j-th column of A is stored in the j-th column of the array AB
          as follows:
          if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+ka).

          On exit, the contents of AB are destroyed.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KA+1.
[in,out]BB
          BB is COMPLEX*16 array, dimension (LDBB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix B, stored in the first kb+1 rows of the array.  The
          j-th column of B is stored in the j-th column of the array BB
          as follows:
          if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
          if UPLO = 'L', BB(1+i-j,j)    = B(i,j) for j<=i<=min(n,j+kb).

          On exit, the factor S from the split Cholesky factorization
          B = S**H*S, as returned by ZPBSTF.
[in]LDBB
          LDBB is INTEGER
          The leading dimension of the array BB.  LDBB >= KB+1.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
          eigenvectors, with the i-th column of Z holding the
          eigenvector associated with W(i). The eigenvectors are
          normalized so that Z**H*B*Z = I.
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= N.
[out]WORK
          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
          On exit, if INFO=0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= N.
          If JOBZ = 'V' and N > 1, LWORK >= 2*N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK, RWORK and
          IWORK arrays, returns these values as the first entries of
          the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
          On exit, if INFO=0, RWORK(1) returns the optimal LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The dimension of array RWORK.
          If N <= 1,               LRWORK >= 1.
          If JOBZ = 'N' and N > 1, LRWORK >= N.
          If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.

          If LRWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO=0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of array IWORK.
          If JOBZ = 'N' or N <= 1, LIWORK >= 1.
          If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, and i is:
             <= N:  the algorithm 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 ZPBSTF
                    returned INFO = i: 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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 249 of file zhbgvd.f.

252*
253* -- LAPACK driver routine --
254* -- LAPACK is a software package provided by Univ. of Tennessee, --
255* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
256*
257* .. Scalar Arguments ..
258 CHARACTER JOBZ, UPLO
259 INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LRWORK,
260 $ LWORK, N
261* ..
262* .. Array Arguments ..
263 INTEGER IWORK( * )
264 DOUBLE PRECISION RWORK( * ), W( * )
265 COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
266 $ Z( LDZ, * )
267* ..
268*
269* =====================================================================
270*
271* .. Parameters ..
272 COMPLEX*16 CONE, CZERO
273 parameter( cone = ( 1.0d+0, 0.0d+0 ),
274 $ czero = ( 0.0d+0, 0.0d+0 ) )
275* ..
276* .. Local Scalars ..
277 LOGICAL LQUERY, UPPER, WANTZ
278 CHARACTER VECT
279 INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLRWK,
280 $ LLWK2, LRWMIN, LWMIN
281* ..
282* .. External Functions ..
283 LOGICAL LSAME
284 EXTERNAL lsame
285* ..
286* .. External Subroutines ..
287 EXTERNAL dsterf, xerbla, zgemm, zhbgst, zhbtrd, zlacpy,
288 $ zpbstf, zstedc
289* ..
290* .. Executable Statements ..
291*
292* Test the input parameters.
293*
294 wantz = lsame( jobz, 'V' )
295 upper = lsame( uplo, 'U' )
296 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
297*
298 info = 0
299 IF( n.LE.1 ) THEN
300 lwmin = 1+n
301 lrwmin = 1+n
302 liwmin = 1
303 ELSE IF( wantz ) THEN
304 lwmin = 2*n**2
305 lrwmin = 1 + 5*n + 2*n**2
306 liwmin = 3 + 5*n
307 ELSE
308 lwmin = n
309 lrwmin = n
310 liwmin = 1
311 END IF
312 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
313 info = -1
314 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
315 info = -2
316 ELSE IF( n.LT.0 ) THEN
317 info = -3
318 ELSE IF( ka.LT.0 ) THEN
319 info = -4
320 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
321 info = -5
322 ELSE IF( ldab.LT.ka+1 ) THEN
323 info = -7
324 ELSE IF( ldbb.LT.kb+1 ) THEN
325 info = -9
326 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
327 info = -12
328 END IF
329*
330 IF( info.EQ.0 ) THEN
331 work( 1 ) = lwmin
332 rwork( 1 ) = lrwmin
333 iwork( 1 ) = liwmin
334*
335 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
336 info = -14
337 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
338 info = -16
339 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
340 info = -18
341 END IF
342 END IF
343*
344 IF( info.NE.0 ) THEN
345 CALL xerbla( 'ZHBGVD', -info )
346 RETURN
347 ELSE IF( lquery ) THEN
348 RETURN
349 END IF
350*
351* Quick return if possible
352*
353 IF( n.EQ.0 )
354 $ RETURN
355*
356* Form a split Cholesky factorization of B.
357*
358 CALL zpbstf( uplo, n, kb, bb, ldbb, info )
359 IF( info.NE.0 ) THEN
360 info = n + info
361 RETURN
362 END IF
363*
364* Transform problem to standard eigenvalue problem.
365*
366 inde = 1
367 indwrk = inde + n
368 indwk2 = 1 + n*n
369 llwk2 = lwork - indwk2 + 2
370 llrwk = lrwork - indwrk + 2
371 CALL zhbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
372 $ work, rwork, iinfo )
373*
374* Reduce Hermitian band matrix to tridiagonal form.
375*
376 IF( wantz ) THEN
377 vect = 'U'
378 ELSE
379 vect = 'N'
380 END IF
381 CALL zhbtrd( vect, uplo, n, ka, ab, ldab, w, rwork( inde ), z,
382 $ ldz, work, iinfo )
383*
384* For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEDC.
385*
386 IF( .NOT.wantz ) THEN
387 CALL dsterf( n, w, rwork( inde ), info )
388 ELSE
389 CALL zstedc( 'I', n, w, rwork( inde ), work, n, work( indwk2 ),
390 $ llwk2, rwork( indwrk ), llrwk, iwork, liwork,
391 $ info )
392 CALL zgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
393 $ work( indwk2 ), n )
394 CALL zlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
395 END IF
396*
397 work( 1 ) = lwmin
398 rwork( 1 ) = lrwmin
399 iwork( 1 ) = liwmin
400 RETURN
401*
402* End of ZHBGVD
403*
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:86
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
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 zpbstf(UPLO, N, KD, AB, LDAB, INFO)
ZPBSTF
Definition: zpbstf.f:153
subroutine zstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZSTEDC
Definition: zstedc.f:212
subroutine zhbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
ZHBTRD
Definition: zhbtrd.f:163
Here is the call graph for this function:
Here is the caller graph for this function: