LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine ssbev_2stage ( 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  INFO 
)

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

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

Purpose:
 SSBEV_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
 a real symmetric band matrix A using the 2stage technique for
 the reduction to tridiagonal.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
                  Not available in this release.
[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 length of the array WORK. LWORK >= 1, when N <= 1;
          otherwise  
          If JOBZ = 'N' and N > 1, LWORK must be queried.
                                   LWORK = MAX(1, dimension) where
                                   dimension = (2KD+1)*N + KD*NTHREADS + N
                                   where KD is the size of the band.
                                   NTHREADS is the number of threads used when
                                   openMP compilation is enabled, otherwise =1.
          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK 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.
Date
December 2016
Further Details:
  All details about the 2stage techniques are available in:

  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
  Parallel reduction to condensed forms for symmetric eigenvalue problems
  using aggregated fine-grained and memory-aware kernels. In Proceedings
  of 2011 International Conference for High Performance Computing,
  Networking, Storage and Analysis (SC '11), New York, NY, USA,
  Article 8 , 11 pages.
  http://doi.acm.org/10.1145/2063384.2063394

  A. Haidar, J. Kurzak, P. Luszczek, 2013.
  An improved parallel singular value algorithm and its implementation 
  for multicore hardware, In Proceedings of 2013 International Conference
  for High Performance Computing, Networking, Storage and Analysis (SC '13).
  Denver, Colorado, USA, 2013.
  Article 90, 12 pages.
  http://doi.acm.org/10.1145/2503210.2503292

  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
  A novel hybrid CPU-GPU generalized eigensolver for electronic structure 
  calculations based on fine-grained memory aware tasks.
  International Journal of High Performance Computing Applications.
  Volume 28 Issue 2, Pages 196-209, May 2014.
  http://hpc.sagepub.com/content/28/2/196 

Definition at line 206 of file ssbev_2stage.f.

206 *
207  IMPLICIT NONE
208 *
209 * -- LAPACK driver routine (version 3.7.0) --
210 * -- LAPACK is a software package provided by Univ. of Tennessee, --
211 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
212 * December 2016
213 *
214 * .. Scalar Arguments ..
215  CHARACTER jobz, uplo
216  INTEGER info, kd, ldab, ldz, n, lwork
217 * ..
218 * .. Array Arguments ..
219  REAL ab( ldab, * ), w( * ), work( * ), z( ldz, * )
220 * ..
221 *
222 * =====================================================================
223 *
224 * .. Parameters ..
225  REAL zero, one
226  parameter ( zero = 0.0e0, one = 1.0e0 )
227 * ..
228 * .. Local Scalars ..
229  LOGICAL lower, wantz, lquery
230  INTEGER iinfo, imax, inde, indwrk, iscale,
231  $ llwork, lwmin, lhtrd, lwtrd, ib, indhous
232  REAL anrm, bignum, eps, rmax, rmin, safmin, sigma,
233  $ smlnum
234 * ..
235 * .. External Functions ..
236  LOGICAL lsame
237  INTEGER ilaenv
238  REAL slamch, slansb
239  EXTERNAL lsame, slamch, slansb, ilaenv
240 * ..
241 * .. External Subroutines ..
242  EXTERNAL slascl, sscal, ssteqr, ssterf, xerbla
243  $ ssytrd_sb2st
244 * ..
245 * .. Intrinsic Functions ..
246  INTRINSIC sqrt
247 * ..
248 * .. Executable Statements ..
249 *
250 * Test the input parameters.
251 *
252  wantz = lsame( jobz, 'V' )
253  lower = lsame( uplo, 'L' )
254  lquery = ( lwork.EQ.-1 )
255 *
256  info = 0
257  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
258  info = -1
259  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
260  info = -2
261  ELSE IF( n.LT.0 ) THEN
262  info = -3
263  ELSE IF( kd.LT.0 ) THEN
264  info = -4
265  ELSE IF( ldab.LT.kd+1 ) THEN
266  info = -6
267  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
268  info = -9
269  END IF
270 *
271  IF( info.EQ.0 ) THEN
272  IF( n.LE.1 ) THEN
273  lwmin = 1
274  work( 1 ) = lwmin
275  ELSE
276  ib = ilaenv( 18, 'SSYTRD_SB2ST', jobz, n, kd, -1, -1 )
277  lhtrd = ilaenv( 19, 'SSYTRD_SB2ST', jobz, n, kd, ib, -1 )
278  lwtrd = ilaenv( 20, 'SSYTRD_SB2ST', jobz, n, kd, ib, -1 )
279  lwmin = n + lhtrd + lwtrd
280  work( 1 ) = lwmin
281  ENDIF
282 *
283  IF( lwork.LT.lwmin .AND. .NOT.lquery )
284  $ info = -11
285  END IF
286 *
287  IF( info.NE.0 ) THEN
288  CALL xerbla( 'SSBEV_2STAGE ', -info )
289  RETURN
290  ELSE IF( lquery ) THEN
291  RETURN
292  END IF
293 *
294 * Quick return if possible
295 *
296  IF( n.EQ.0 )
297  $ RETURN
298 *
299  IF( n.EQ.1 ) THEN
300  IF( lower ) THEN
301  w( 1 ) = ab( 1, 1 )
302  ELSE
303  w( 1 ) = ab( kd+1, 1 )
304  END IF
305  IF( wantz )
306  $ z( 1, 1 ) = one
307  RETURN
308  END IF
309 *
310 * Get machine constants.
311 *
312  safmin = slamch( 'Safe minimum' )
313  eps = slamch( 'Precision' )
314  smlnum = safmin / eps
315  bignum = one / smlnum
316  rmin = sqrt( smlnum )
317  rmax = sqrt( bignum )
318 *
319 * Scale matrix to allowable range, if necessary.
320 *
321  anrm = slansb( 'M', uplo, n, kd, ab, ldab, work )
322  iscale = 0
323  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
324  iscale = 1
325  sigma = rmin / anrm
326  ELSE IF( anrm.GT.rmax ) THEN
327  iscale = 1
328  sigma = rmax / anrm
329  END IF
330  IF( iscale.EQ.1 ) THEN
331  IF( lower ) THEN
332  CALL slascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
333  ELSE
334  CALL slascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
335  END IF
336  END IF
337 *
338 * Call SSYTRD_SB2ST to reduce symmetric band matrix to tridiagonal form.
339 *
340  inde = 1
341  indhous = inde + n
342  indwrk = indhous + lhtrd
343  llwork = lwork - indwrk + 1
344 *
345  CALL ssytrd_sb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
346  $ work( inde ), work( indhous ), lhtrd,
347  $ work( indwrk ), llwork, iinfo )
348 *
349 * For eigenvalues only, call SSTERF. For eigenvectors, call SSTEQR.
350 *
351  IF( .NOT.wantz ) THEN
352  CALL ssterf( n, w, work( inde ), info )
353  ELSE
354  CALL ssteqr( jobz, n, w, work( inde ), z, ldz, work( indwrk ),
355  $ info )
356  END IF
357 *
358 * If matrix was scaled, then rescale eigenvalues appropriately.
359 *
360  IF( iscale.EQ.1 ) THEN
361  IF( info.EQ.0 ) THEN
362  imax = n
363  ELSE
364  imax = info - 1
365  END IF
366  CALL sscal( imax, one / sigma, w, 1 )
367  END IF
368 *
369 * Set WORK(1) to optimal workspace size.
370 *
371  work( 1 ) = lwmin
372 *
373  RETURN
374 *
375 * End of SSBEV_2STAGE
376 *
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
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
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
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88

Here is the call graph for this function:

Here is the caller graph for this function: