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

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

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

Purpose:
 SSBEV computes all the eigenvalues and, optionally, eigenvectors of
 a real symmetric 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 REAL array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the symmetric 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 REAL 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 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 148 of file ssbev.f.

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