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

◆ ssbevd()

subroutine ssbevd ( 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  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

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

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

Purpose:
 SSBEVD computes all the eigenvalues and, optionally, eigenvectors of
 a real symmetric 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 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 (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 > 2, LWORK must be at least 2*N.
          If JOBZ  = 'V' and N > 2, LWORK must be at least
                         ( 1 + 5*N + 2*N**2 ).

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK 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 the array IWORK.
          If JOBZ  = 'N' or N <= 1, LIWORK must be at least 1.
          If JOBZ  = 'V' and N > 2, 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 and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK 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 185 of file ssbevd.f.

187*
188* -- LAPACK driver routine --
189* -- LAPACK is a software package provided by Univ. of Tennessee, --
190* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191*
192* .. Scalar Arguments ..
193 CHARACTER JOBZ, UPLO
194 INTEGER INFO, KD, LDAB, LDZ, LIWORK, LWORK, N
195* ..
196* .. Array Arguments ..
197 INTEGER IWORK( * )
198 REAL AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
199* ..
200*
201* =====================================================================
202*
203* .. Parameters ..
204 REAL ZERO, ONE
205 parameter( zero = 0.0e+0, one = 1.0e+0 )
206* ..
207* .. Local Scalars ..
208 LOGICAL LOWER, LQUERY, WANTZ
209 INTEGER IINFO, INDE, INDWK2, INDWRK, ISCALE, LIWMIN,
210 $ LLWRK2, LWMIN
211 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
212 $ SMLNUM
213* ..
214* .. External Functions ..
215 LOGICAL LSAME
216 REAL SLAMCH, SLANSB, SROUNDUP_LWORK
218* ..
219* .. External Subroutines ..
220 EXTERNAL sgemm, slacpy, slascl, ssbtrd, sscal, sstedc,
221 $ ssterf, xerbla
222* ..
223* .. Intrinsic Functions ..
224 INTRINSIC sqrt
225* ..
226* .. Executable Statements ..
227*
228* Test the input parameters.
229*
230 wantz = lsame( jobz, 'V' )
231 lower = lsame( uplo, 'L' )
232 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
233*
234 info = 0
235 IF( n.LE.1 ) THEN
236 liwmin = 1
237 lwmin = 1
238 ELSE
239 IF( wantz ) THEN
240 liwmin = 3 + 5*n
241 lwmin = 1 + 5*n + 2*n**2
242 ELSE
243 liwmin = 1
244 lwmin = 2*n
245 END IF
246 END IF
247 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
248 info = -1
249 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
250 info = -2
251 ELSE IF( n.LT.0 ) THEN
252 info = -3
253 ELSE IF( kd.LT.0 ) THEN
254 info = -4
255 ELSE IF( ldab.LT.kd+1 ) THEN
256 info = -6
257 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
258 info = -9
259 END IF
260*
261 IF( info.EQ.0 ) THEN
262 work( 1 ) = sroundup_lwork(lwmin)
263 iwork( 1 ) = liwmin
264*
265 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
266 info = -11
267 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
268 info = -13
269 END IF
270 END IF
271*
272 IF( info.NE.0 ) THEN
273 CALL xerbla( 'SSBEVD', -info )
274 RETURN
275 ELSE IF( lquery ) THEN
276 RETURN
277 END IF
278*
279* Quick return if possible
280*
281 IF( n.EQ.0 )
282 $ RETURN
283*
284 IF( n.EQ.1 ) THEN
285 w( 1 ) = ab( 1, 1 )
286 IF( wantz )
287 $ z( 1, 1 ) = one
288 RETURN
289 END IF
290*
291* Get machine constants.
292*
293 safmin = slamch( 'Safe minimum' )
294 eps = slamch( 'Precision' )
295 smlnum = safmin / eps
296 bignum = one / smlnum
297 rmin = sqrt( smlnum )
298 rmax = sqrt( bignum )
299*
300* Scale matrix to allowable range, if necessary.
301*
302 anrm = slansb( 'M', uplo, n, kd, ab, ldab, work )
303 iscale = 0
304 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
305 iscale = 1
306 sigma = rmin / anrm
307 ELSE IF( anrm.GT.rmax ) THEN
308 iscale = 1
309 sigma = rmax / anrm
310 END IF
311 IF( iscale.EQ.1 ) THEN
312 IF( lower ) THEN
313 CALL slascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
314 ELSE
315 CALL slascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
316 END IF
317 END IF
318*
319* Call SSBTRD to reduce symmetric band matrix to tridiagonal form.
320*
321 inde = 1
322 indwrk = inde + n
323 indwk2 = indwrk + n*n
324 llwrk2 = lwork - indwk2 + 1
325 CALL ssbtrd( jobz, uplo, n, kd, ab, ldab, w, work( inde ), z, ldz,
326 $ work( indwrk ), iinfo )
327*
328* For eigenvalues only, call SSTERF. For eigenvectors, call SSTEDC.
329*
330 IF( .NOT.wantz ) THEN
331 CALL ssterf( n, w, work( inde ), info )
332 ELSE
333 CALL sstedc( 'I', n, w, work( inde ), work( indwrk ), n,
334 $ work( indwk2 ), llwrk2, iwork, liwork, info )
335 CALL sgemm( 'N', 'N', n, n, n, one, z, ldz, work( indwrk ), n,
336 $ zero, work( indwk2 ), n )
337 CALL slacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
338 END IF
339*
340* If matrix was scaled, then rescale eigenvalues appropriately.
341*
342 IF( iscale.EQ.1 )
343 $ CALL sscal( n, one / sigma, w, 1 )
344*
345 work( 1 ) = sroundup_lwork(lwmin)
346 iwork( 1 ) = liwmin
347 RETURN
348*
349* End of SSBEVD
350*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine ssbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
SSBTRD
Definition ssbtrd.f:163
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
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,...
Definition slansb.f:129
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:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sstedc(compz, n, d, e, z, ldz, work, lwork, iwork, liwork, info)
SSTEDC
Definition sstedc.f:182
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
Here is the call graph for this function:
Here is the caller graph for this function: