LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine chbev ( character  JOBZ,
character  UPLO,
integer  N,
integer  KD,
complex, dimension( ldab, * )  AB,
integer  LDAB,
real, dimension( * )  W,
complex, dimension( ldz, * )  Z,
integer  LDZ,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

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

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

Purpose:
 CHBEV 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 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 REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX 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 array, dimension (N)
[out]RWORK
          RWORK is REAL 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.
Date
November 2011

Definition at line 154 of file chbev.f.

154 *
155 * -- LAPACK driver routine (version 3.4.0) --
156 * -- LAPACK is a software package provided by Univ. of Tennessee, --
157 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
158 * November 2011
159 *
160 * .. Scalar Arguments ..
161  CHARACTER jobz, uplo
162  INTEGER info, kd, ldab, ldz, n
163 * ..
164 * .. Array Arguments ..
165  REAL rwork( * ), w( * )
166  COMPLEX ab( ldab, * ), work( * ), z( ldz, * )
167 * ..
168 *
169 * =====================================================================
170 *
171 * .. Parameters ..
172  REAL zero, one
173  parameter ( zero = 0.0e0, one = 1.0e0 )
174 * ..
175 * .. Local Scalars ..
176  LOGICAL lower, wantz
177  INTEGER iinfo, imax, inde, indrwk, iscale
178  REAL anrm, bignum, eps, rmax, rmin, safmin, sigma,
179  $ smlnum
180 * ..
181 * .. External Functions ..
182  LOGICAL lsame
183  REAL clanhb, slamch
184  EXTERNAL lsame, clanhb, slamch
185 * ..
186 * .. External Subroutines ..
187  EXTERNAL chbtrd, clascl, csteqr, sscal, ssterf, xerbla
188 * ..
189 * .. Intrinsic Functions ..
190  INTRINSIC sqrt
191 * ..
192 * .. Executable Statements ..
193 *
194 * Test the input parameters.
195 *
196  wantz = lsame( jobz, 'V' )
197  lower = lsame( uplo, 'L' )
198 *
199  info = 0
200  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
201  info = -1
202  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
203  info = -2
204  ELSE IF( n.LT.0 ) THEN
205  info = -3
206  ELSE IF( kd.LT.0 ) THEN
207  info = -4
208  ELSE IF( ldab.LT.kd+1 ) THEN
209  info = -6
210  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
211  info = -9
212  END IF
213 *
214  IF( info.NE.0 ) THEN
215  CALL xerbla( 'CHBEV ', -info )
216  RETURN
217  END IF
218 *
219 * Quick return if possible
220 *
221  IF( n.EQ.0 )
222  $ RETURN
223 *
224  IF( n.EQ.1 ) THEN
225  IF( lower ) THEN
226  w( 1 ) = ab( 1, 1 )
227  ELSE
228  w( 1 ) = ab( kd+1, 1 )
229  END IF
230  IF( wantz )
231  $ z( 1, 1 ) = one
232  RETURN
233  END IF
234 *
235 * Get machine constants.
236 *
237  safmin = slamch( 'Safe minimum' )
238  eps = slamch( 'Precision' )
239  smlnum = safmin / eps
240  bignum = one / smlnum
241  rmin = sqrt( smlnum )
242  rmax = sqrt( bignum )
243 *
244 * Scale matrix to allowable range, if necessary.
245 *
246  anrm = clanhb( 'M', uplo, n, kd, ab, ldab, rwork )
247  iscale = 0
248  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
249  iscale = 1
250  sigma = rmin / anrm
251  ELSE IF( anrm.GT.rmax ) THEN
252  iscale = 1
253  sigma = rmax / anrm
254  END IF
255  IF( iscale.EQ.1 ) THEN
256  IF( lower ) THEN
257  CALL clascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
258  ELSE
259  CALL clascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
260  END IF
261  END IF
262 *
263 * Call CHBTRD to reduce Hermitian band matrix to tridiagonal form.
264 *
265  inde = 1
266  CALL chbtrd( jobz, uplo, n, kd, ab, ldab, w, rwork( inde ), z,
267  $ ldz, work, iinfo )
268 *
269 * For eigenvalues only, call SSTERF. For eigenvectors, call CSTEQR.
270 *
271  IF( .NOT.wantz ) THEN
272  CALL ssterf( n, w, rwork( inde ), info )
273  ELSE
274  indrwk = inde + n
275  CALL csteqr( jobz, n, w, rwork( inde ), z, ldz,
276  $ rwork( indrwk ), info )
277  END IF
278 *
279 * If matrix was scaled, then rescale eigenvalues appropriately.
280 *
281  IF( iscale.EQ.1 ) THEN
282  IF( info.EQ.0 ) THEN
283  imax = n
284  ELSE
285  imax = info - 1
286  END IF
287  CALL sscal( imax, one / sigma, w, 1 )
288  END IF
289 *
290  RETURN
291 *
292 * End of CHBEV
293 *
subroutine clascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: clascl.f:145
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine chbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
CHBTRD
Definition: chbtrd.f:165
subroutine csteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
CSTEQR
Definition: csteqr.f:134
real function clanhb(NORM, UPLO, N, K, AB, LDAB, WORK)
CLANHB 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: clanhb.f:134
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
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: