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

◆ zhbev()

subroutine zhbev ( 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,
double precision, dimension( * )  rwork,
integer  info 
)

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

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

Purpose:
 ZHBEV computes all the eigenvalues and, optionally, eigenvectors of
 a complex Hermitian band matrix A.
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 (N)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(1,3*N-2))
[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 150 of file zhbev.f.

152*
153* -- LAPACK driver routine --
154* -- LAPACK is a software package provided by Univ. of Tennessee, --
155* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
156*
157* .. Scalar Arguments ..
158 CHARACTER JOBZ, UPLO
159 INTEGER INFO, KD, LDAB, LDZ, N
160* ..
161* .. Array Arguments ..
162 DOUBLE PRECISION RWORK( * ), W( * )
163 COMPLEX*16 AB( LDAB, * ), WORK( * ), Z( LDZ, * )
164* ..
165*
166* =====================================================================
167*
168* .. Parameters ..
169 DOUBLE PRECISION ZERO, ONE
170 parameter( zero = 0.0d0, one = 1.0d0 )
171* ..
172* .. Local Scalars ..
173 LOGICAL LOWER, WANTZ
174 INTEGER IINFO, IMAX, INDE, INDRWK, ISCALE
175 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
176 $ SMLNUM
177* ..
178* .. External Functions ..
179 LOGICAL LSAME
180 DOUBLE PRECISION DLAMCH, ZLANHB
181 EXTERNAL lsame, dlamch, zlanhb
182* ..
183* .. External Subroutines ..
184 EXTERNAL dscal, dsterf, xerbla, zhbtrd, zlascl, zsteqr
185* ..
186* .. Intrinsic Functions ..
187 INTRINSIC sqrt
188* ..
189* .. Executable Statements ..
190*
191* Test the input parameters.
192*
193 wantz = lsame( jobz, 'V' )
194 lower = lsame( uplo, 'L' )
195*
196 info = 0
197 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
198 info = -1
199 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
200 info = -2
201 ELSE IF( n.LT.0 ) THEN
202 info = -3
203 ELSE IF( kd.LT.0 ) THEN
204 info = -4
205 ELSE IF( ldab.LT.kd+1 ) THEN
206 info = -6
207 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
208 info = -9
209 END IF
210*
211 IF( info.NE.0 ) THEN
212 CALL xerbla( 'ZHBEV ', -info )
213 RETURN
214 END IF
215*
216* Quick return if possible
217*
218 IF( n.EQ.0 )
219 $ RETURN
220*
221 IF( n.EQ.1 ) THEN
222 IF( lower ) THEN
223 w( 1 ) = dble( ab( 1, 1 ) )
224 ELSE
225 w( 1 ) = dble( ab( kd+1, 1 ) )
226 END IF
227 IF( wantz )
228 $ z( 1, 1 ) = one
229 RETURN
230 END IF
231*
232* Get machine constants.
233*
234 safmin = dlamch( 'Safe minimum' )
235 eps = dlamch( 'Precision' )
236 smlnum = safmin / eps
237 bignum = one / smlnum
238 rmin = sqrt( smlnum )
239 rmax = sqrt( bignum )
240*
241* Scale matrix to allowable range, if necessary.
242*
243 anrm = zlanhb( 'M', uplo, n, kd, ab, ldab, rwork )
244 iscale = 0
245 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
246 iscale = 1
247 sigma = rmin / anrm
248 ELSE IF( anrm.GT.rmax ) THEN
249 iscale = 1
250 sigma = rmax / anrm
251 END IF
252 IF( iscale.EQ.1 ) THEN
253 IF( lower ) THEN
254 CALL zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
255 ELSE
256 CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
257 END IF
258 END IF
259*
260* Call ZHBTRD to reduce Hermitian band matrix to tridiagonal form.
261*
262 inde = 1
263 CALL zhbtrd( jobz, uplo, n, kd, ab, ldab, w, rwork( inde ), z,
264 $ ldz, work, iinfo )
265*
266* For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEQR.
267*
268 IF( .NOT.wantz ) THEN
269 CALL dsterf( n, w, rwork( inde ), info )
270 ELSE
271 indrwk = inde + n
272 CALL zsteqr( jobz, n, w, rwork( inde ), z, ldz,
273 $ rwork( indrwk ), info )
274 END IF
275*
276* If matrix was scaled, then rescale eigenvalues appropriately.
277*
278 IF( iscale.EQ.1 ) THEN
279 IF( info.EQ.0 ) THEN
280 imax = n
281 ELSE
282 imax = info - 1
283 END IF
284 CALL dscal( imax, one / sigma, w, 1 )
285 END IF
286*
287 RETURN
288*
289* End of ZHBEV
290*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
ZHBTRD
Definition zhbtrd.f:163
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 zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
Definition zsteqr.f:132
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: