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

◆ zhbevd()

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.
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.

Definition at line 207 of file zhbevd.f.

209*
210* -- LAPACK driver routine --
211* -- LAPACK is a software package provided by Univ. of Tennessee, --
212* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
213*
214* .. Scalar Arguments ..
215 CHARACTER JOBZ, UPLO
216 INTEGER INFO, KD, LDAB, LDZ, LIWORK, LRWORK, LWORK, N
217* ..
218* .. Array Arguments ..
219 INTEGER IWORK( * )
220 DOUBLE PRECISION RWORK( * ), W( * )
221 COMPLEX*16 AB( LDAB, * ), WORK( * ), Z( LDZ, * )
222* ..
223*
224* =====================================================================
225*
226* .. Parameters ..
227 DOUBLE PRECISION ZERO, ONE
228 parameter( zero = 0.0d0, one = 1.0d0 )
229 COMPLEX*16 CZERO, CONE
230 parameter( czero = ( 0.0d0, 0.0d0 ),
231 $ cone = ( 1.0d0, 0.0d0 ) )
232* ..
233* .. Local Scalars ..
234 LOGICAL LOWER, LQUERY, WANTZ
235 INTEGER IINFO, IMAX, INDE, INDWK2, INDWRK, ISCALE,
236 $ LIWMIN, LLRWK, LLWK2, LRWMIN, LWMIN
237 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
238 $ SMLNUM
239* ..
240* .. External Functions ..
241 LOGICAL LSAME
242 DOUBLE PRECISION DLAMCH, ZLANHB
243 EXTERNAL lsame, dlamch, zlanhb
244* ..
245* .. External Subroutines ..
246 EXTERNAL dscal, dsterf, xerbla, zgemm, zhbtrd, zlacpy,
247 $ zlascl, zstedc
248* ..
249* .. Intrinsic Functions ..
250 INTRINSIC sqrt
251* ..
252* .. Executable Statements ..
253*
254* Test the input parameters.
255*
256 wantz = lsame( jobz, 'V' )
257 lower = lsame( uplo, 'L' )
258 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 .OR. lrwork.EQ.-1 )
259*
260 info = 0
261 IF( n.LE.1 ) THEN
262 lwmin = 1
263 lrwmin = 1
264 liwmin = 1
265 ELSE
266 IF( wantz ) THEN
267 lwmin = 2*n**2
268 lrwmin = 1 + 5*n + 2*n**2
269 liwmin = 3 + 5*n
270 ELSE
271 lwmin = n
272 lrwmin = n
273 liwmin = 1
274 END IF
275 END IF
276 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
277 info = -1
278 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
279 info = -2
280 ELSE IF( n.LT.0 ) THEN
281 info = -3
282 ELSE IF( kd.LT.0 ) THEN
283 info = -4
284 ELSE IF( ldab.LT.kd+1 ) THEN
285 info = -6
286 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
287 info = -9
288 END IF
289*
290 IF( info.EQ.0 ) THEN
291 work( 1 ) = lwmin
292 rwork( 1 ) = lrwmin
293 iwork( 1 ) = liwmin
294*
295 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
296 info = -11
297 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
298 info = -13
299 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
300 info = -15
301 END IF
302 END IF
303*
304 IF( info.NE.0 ) THEN
305 CALL xerbla( 'ZHBEVD', -info )
306 RETURN
307 ELSE IF( lquery ) THEN
308 RETURN
309 END IF
310*
311* Quick return if possible
312*
313 IF( n.EQ.0 )
314 $ RETURN
315*
316 IF( n.EQ.1 ) THEN
317 w( 1 ) = dble( ab( 1, 1 ) )
318 IF( wantz )
319 $ z( 1, 1 ) = cone
320 RETURN
321 END IF
322*
323* Get machine constants.
324*
325 safmin = dlamch( 'Safe minimum' )
326 eps = dlamch( 'Precision' )
327 smlnum = safmin / eps
328 bignum = one / smlnum
329 rmin = sqrt( smlnum )
330 rmax = sqrt( bignum )
331*
332* Scale matrix to allowable range, if necessary.
333*
334 anrm = zlanhb( 'M', uplo, n, kd, ab, ldab, rwork )
335 iscale = 0
336 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
337 iscale = 1
338 sigma = rmin / anrm
339 ELSE IF( anrm.GT.rmax ) THEN
340 iscale = 1
341 sigma = rmax / anrm
342 END IF
343 IF( iscale.EQ.1 ) THEN
344 IF( lower ) THEN
345 CALL zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
346 ELSE
347 CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
348 END IF
349 END IF
350*
351* Call ZHBTRD to reduce Hermitian band matrix to tridiagonal form.
352*
353 inde = 1
354 indwrk = inde + n
355 indwk2 = 1 + n*n
356 llwk2 = lwork - indwk2 + 1
357 llrwk = lrwork - indwrk + 1
358 CALL zhbtrd( jobz, uplo, n, kd, ab, ldab, w, rwork( inde ), z,
359 $ ldz, work, iinfo )
360*
361* For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEDC.
362*
363 IF( .NOT.wantz ) THEN
364 CALL dsterf( n, w, rwork( inde ), info )
365 ELSE
366 CALL zstedc( 'I', n, w, rwork( inde ), work, n, work( indwk2 ),
367 $ llwk2, rwork( indwrk ), llrwk, iwork, liwork,
368 $ info )
369 CALL zgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
370 $ work( indwk2 ), n )
371 CALL zlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
372 END IF
373*
374* If matrix was scaled, then rescale eigenvalues appropriately.
375*
376 IF( iscale.EQ.1 ) THEN
377 IF( info.EQ.0 ) THEN
378 imax = n
379 ELSE
380 imax = info - 1
381 END IF
382 CALL dscal( imax, one / sigma, w, 1 )
383 END IF
384*
385 work( 1 ) = lwmin
386 rwork( 1 ) = lrwmin
387 iwork( 1 ) = liwmin
388 RETURN
389*
390* End of ZHBEVD
391*
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 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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
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,...
Definition zlanhb.f:132
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:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
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: