LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhbevd ( character  JOBZ,
character  UPLO,
integer  N,
integer  KD,
complex*16, dimension( ldab, * )  AB,
integer  LDAB,
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 
)

ZHBEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

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

Purpose:
 ZHBEVD computes all the eigenvalues and, optionally, eigenvectors of
 a complex Hermitian band matrix A.  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 triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]KD
          KD is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'.  KD >= 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 KD+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(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).

          On exit, AB is overwritten by values generated during the
          reduction to tridiagonal form.  If UPLO = 'U', the first
          superdiagonal and the diagonal of the tridiagonal matrix T
          are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
          the diagonal and first subdiagonal of T are returned in the
          first two rows of AB.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KD + 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 orthonormal
          eigenvectors of the matrix A, with the i-th column of Z
          holding the eigenvector associated with W(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 >= max(1,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 must be at least 1.
          If JOBZ = 'N' and N > 1, LWORK must be at least N.
          If JOBZ = 'V' and N > 1, LWORK must be at least 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 (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 must be at least 1.
          If JOBZ = 'N' and N > 1, LRWORK must be at least N.
          If JOBZ = 'V' and N > 1, LRWORK must be at least
                        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 must be at least 1.
          If JOBZ = 'V' and N > 1, LIWORK must be at least 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, the algorithm failed to converge; i
                off-diagonal elements of an intermediate tridiagonal
                form did not converge to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 217 of file zhbevd.f.

217 *
218 * -- LAPACK driver routine (version 3.4.0) --
219 * -- LAPACK is a software package provided by Univ. of Tennessee, --
220 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
221 * November 2011
222 *
223 * .. Scalar Arguments ..
224  CHARACTER jobz, uplo
225  INTEGER info, kd, ldab, ldz, liwork, lrwork, lwork, n
226 * ..
227 * .. Array Arguments ..
228  INTEGER iwork( * )
229  DOUBLE PRECISION rwork( * ), w( * )
230  COMPLEX*16 ab( ldab, * ), work( * ), z( ldz, * )
231 * ..
232 *
233 * =====================================================================
234 *
235 * .. Parameters ..
236  DOUBLE PRECISION zero, one
237  parameter ( zero = 0.0d0, one = 1.0d0 )
238  COMPLEX*16 czero, cone
239  parameter ( czero = ( 0.0d0, 0.0d0 ),
240  $ cone = ( 1.0d0, 0.0d0 ) )
241 * ..
242 * .. Local Scalars ..
243  LOGICAL lower, lquery, wantz
244  INTEGER iinfo, imax, inde, indwk2, indwrk, iscale,
245  $ liwmin, llrwk, llwk2, lrwmin, lwmin
246  DOUBLE PRECISION anrm, bignum, eps, rmax, rmin, safmin, sigma,
247  $ smlnum
248 * ..
249 * .. External Functions ..
250  LOGICAL lsame
251  DOUBLE PRECISION dlamch, zlanhb
252  EXTERNAL lsame, dlamch, zlanhb
253 * ..
254 * .. External Subroutines ..
255  EXTERNAL dscal, dsterf, xerbla, zgemm, zhbtrd, zlacpy,
256  $ zlascl, zstedc
257 * ..
258 * .. Intrinsic Functions ..
259  INTRINSIC sqrt
260 * ..
261 * .. Executable Statements ..
262 *
263 * Test the input parameters.
264 *
265  wantz = lsame( jobz, 'V' )
266  lower = lsame( uplo, 'L' )
267  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 .OR. lrwork.EQ.-1 )
268 *
269  info = 0
270  IF( n.LE.1 ) THEN
271  lwmin = 1
272  lrwmin = 1
273  liwmin = 1
274  ELSE
275  IF( wantz ) THEN
276  lwmin = 2*n**2
277  lrwmin = 1 + 5*n + 2*n**2
278  liwmin = 3 + 5*n
279  ELSE
280  lwmin = n
281  lrwmin = n
282  liwmin = 1
283  END IF
284  END IF
285  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
286  info = -1
287  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
288  info = -2
289  ELSE IF( n.LT.0 ) THEN
290  info = -3
291  ELSE IF( kd.LT.0 ) THEN
292  info = -4
293  ELSE IF( ldab.LT.kd+1 ) THEN
294  info = -6
295  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
296  info = -9
297  END IF
298 *
299  IF( info.EQ.0 ) THEN
300  work( 1 ) = lwmin
301  rwork( 1 ) = lrwmin
302  iwork( 1 ) = liwmin
303 *
304  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
305  info = -11
306  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
307  info = -13
308  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
309  info = -15
310  END IF
311  END IF
312 *
313  IF( info.NE.0 ) THEN
314  CALL xerbla( 'ZHBEVD', -info )
315  RETURN
316  ELSE IF( lquery ) THEN
317  RETURN
318  END IF
319 *
320 * Quick return if possible
321 *
322  IF( n.EQ.0 )
323  $ RETURN
324 *
325  IF( n.EQ.1 ) THEN
326  w( 1 ) = ab( 1, 1 )
327  IF( wantz )
328  $ z( 1, 1 ) = cone
329  RETURN
330  END IF
331 *
332 * Get machine constants.
333 *
334  safmin = dlamch( 'Safe minimum' )
335  eps = dlamch( 'Precision' )
336  smlnum = safmin / eps
337  bignum = one / smlnum
338  rmin = sqrt( smlnum )
339  rmax = sqrt( bignum )
340 *
341 * Scale matrix to allowable range, if necessary.
342 *
343  anrm = zlanhb( 'M', uplo, n, kd, ab, ldab, rwork )
344  iscale = 0
345  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
346  iscale = 1
347  sigma = rmin / anrm
348  ELSE IF( anrm.GT.rmax ) THEN
349  iscale = 1
350  sigma = rmax / anrm
351  END IF
352  IF( iscale.EQ.1 ) THEN
353  IF( lower ) THEN
354  CALL zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
355  ELSE
356  CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
357  END IF
358  END IF
359 *
360 * Call ZHBTRD to reduce Hermitian band matrix to tridiagonal form.
361 *
362  inde = 1
363  indwrk = inde + n
364  indwk2 = 1 + n*n
365  llwk2 = lwork - indwk2 + 1
366  llrwk = lrwork - indwrk + 1
367  CALL zhbtrd( jobz, uplo, n, kd, ab, ldab, w, rwork( inde ), z,
368  $ ldz, work, iinfo )
369 *
370 * For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEDC.
371 *
372  IF( .NOT.wantz ) THEN
373  CALL dsterf( n, w, rwork( inde ), info )
374  ELSE
375  CALL zstedc( 'I', n, w, rwork( inde ), work, n, work( indwk2 ),
376  $ llwk2, rwork( indwrk ), llrwk, iwork, liwork,
377  $ info )
378  CALL zgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
379  $ work( indwk2 ), n )
380  CALL zlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
381  END IF
382 *
383 * If matrix was scaled, then rescale eigenvalues appropriately.
384 *
385  IF( iscale.EQ.1 ) THEN
386  IF( info.EQ.0 ) THEN
387  imax = n
388  ELSE
389  imax = info - 1
390  END IF
391  CALL dscal( imax, one / sigma, w, 1 )
392  END IF
393 *
394  work( 1 ) = lwmin
395  rwork( 1 ) = lrwmin
396  iwork( 1 ) = liwmin
397  RETURN
398 *
399 * End of ZHBEVD
400 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
double precision function zlanhb(NORM, UPLO, N, K, AB, LDAB, WORK)
ZLANHB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a Hermitian band matrix.
Definition: zlanhb.f:134
subroutine zstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZSTEDC
Definition: zstedc.f:215
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine zlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: zlascl.f:145
subroutine zhbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
ZHBTRD
Definition: zhbtrd.f:165
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: