LAPACK 3.12.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.
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 243 of file zhbgvd.f.

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