LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ ssbev_2stage()

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
November 2017
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.8.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 * November 2017
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 ILAENV2STAGE
238  REAL SLAMCH, SLANSB
239  EXTERNAL lsame, slamch, slansb, ilaenv2stage
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 = ilaenv2stage( 2, 'SSYTRD_SB2ST', jobz,
277  $ n, kd, -1, -1 )
278  lhtrd = ilaenv2stage( 3, 'SSYTRD_SB2ST', jobz,
279  $ n, kd, ib, -1 )
280  lwtrd = ilaenv2stage( 4, 'SSYTRD_SB2ST', jobz,
281  $ n, kd, ib, -1 )
282  lwmin = n + lhtrd + lwtrd
283  work( 1 ) = lwmin
284  ENDIF
285 *
286  IF( lwork.LT.lwmin .AND. .NOT.lquery )
287  $ info = -11
288  END IF
289 *
290  IF( info.NE.0 ) THEN
291  CALL xerbla( 'SSBEV_2STAGE ', -info )
292  RETURN
293  ELSE IF( lquery ) THEN
294  RETURN
295  END IF
296 *
297 * Quick return if possible
298 *
299  IF( n.EQ.0 )
300  $ RETURN
301 *
302  IF( n.EQ.1 ) THEN
303  IF( lower ) THEN
304  w( 1 ) = ab( 1, 1 )
305  ELSE
306  w( 1 ) = ab( kd+1, 1 )
307  END IF
308  IF( wantz )
309  $ z( 1, 1 ) = one
310  RETURN
311  END IF
312 *
313 * Get machine constants.
314 *
315  safmin = slamch( 'Safe minimum' )
316  eps = slamch( 'Precision' )
317  smlnum = safmin / eps
318  bignum = one / smlnum
319  rmin = sqrt( smlnum )
320  rmax = sqrt( bignum )
321 *
322 * Scale matrix to allowable range, if necessary.
323 *
324  anrm = slansb( 'M', uplo, n, kd, ab, ldab, work )
325  iscale = 0
326  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
327  iscale = 1
328  sigma = rmin / anrm
329  ELSE IF( anrm.GT.rmax ) THEN
330  iscale = 1
331  sigma = rmax / anrm
332  END IF
333  IF( iscale.EQ.1 ) THEN
334  IF( lower ) THEN
335  CALL slascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
336  ELSE
337  CALL slascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
338  END IF
339  END IF
340 *
341 * Call SSYTRD_SB2ST to reduce symmetric band matrix to tridiagonal form.
342 *
343  inde = 1
344  indhous = inde + n
345  indwrk = indhous + lhtrd
346  llwork = lwork - indwrk + 1
347 *
348  CALL ssytrd_sb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
349  $ work( inde ), work( indhous ), lhtrd,
350  $ work( indwrk ), llwork, iinfo )
351 *
352 * For eigenvalues only, call SSTERF. For eigenvectors, call SSTEQR.
353 *
354  IF( .NOT.wantz ) THEN
355  CALL ssterf( n, w, work( inde ), info )
356  ELSE
357  CALL ssteqr( jobz, n, w, work( inde ), z, ldz, work( indwrk ),
358  $ info )
359  END IF
360 *
361 * If matrix was scaled, then rescale eigenvalues appropriately.
362 *
363  IF( iscale.EQ.1 ) THEN
364  IF( info.EQ.0 ) THEN
365  imax = n
366  ELSE
367  imax = info - 1
368  END IF
369  CALL sscal( imax, one / sigma, w, 1 )
370  END IF
371 *
372 * Set WORK(1) to optimal workspace size.
373 *
374  work( 1 ) = lwmin
375 *
376  RETURN
377 *
378 * End of SSBEV_2STAGE
379 *
Here is the call graph for this function:
Here is the caller graph for this function:
slansb
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:131
slascl
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
ssterf
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
sscal
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:81
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
lsame
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
slamch
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:70
ssytrd_sb2st
subroutine ssytrd_sb2st(STAGE1, VECT, UPLO, N, KD, AB, LDAB, D, E, HOUS, LHOUS, WORK, LWORK, INFO)
SSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric tridiagonal form T
Definition: ssytrd_sb2st.F:232
ssteqr
subroutine ssteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
SSTEQR
Definition: ssteqr.f:133
ilaenv2stage
integer function ilaenv2stage(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV2STAGE
Definition: ilaenv2stage.f:151