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

◆ 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.
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 202 of file ssbev_2stage.f.

204*
205 IMPLICIT NONE
206*
207* -- LAPACK driver routine --
208* -- LAPACK is a software package provided by Univ. of Tennessee, --
209* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
210*
211* .. Scalar Arguments ..
212 CHARACTER JOBZ, UPLO
213 INTEGER INFO, KD, LDAB, LDZ, N, LWORK
214* ..
215* .. Array Arguments ..
216 REAL AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
217* ..
218*
219* =====================================================================
220*
221* .. Parameters ..
222 REAL ZERO, ONE
223 parameter( zero = 0.0e0, one = 1.0e0 )
224* ..
225* .. Local Scalars ..
226 LOGICAL LOWER, WANTZ, LQUERY
227 INTEGER IINFO, IMAX, INDE, INDWRK, ISCALE,
228 $ LLWORK, LWMIN, LHTRD, LWTRD, IB, INDHOUS
229 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
230 $ SMLNUM
231* ..
232* .. External Functions ..
233 LOGICAL LSAME
234 INTEGER ILAENV2STAGE
235 REAL SLAMCH, SLANSB, SROUNDUP_LWORK
236 EXTERNAL lsame, slamch, slansb, ilaenv2stage,
238* ..
239* .. External Subroutines ..
240 EXTERNAL slascl, sscal, ssteqr, ssterf, xerbla,
242* ..
243* .. Intrinsic Functions ..
244 INTRINSIC sqrt
245* ..
246* .. Executable Statements ..
247*
248* Test the input parameters.
249*
250 wantz = lsame( jobz, 'V' )
251 lower = lsame( uplo, 'L' )
252 lquery = ( lwork.EQ.-1 )
253*
254 info = 0
255 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
256 info = -1
257 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
258 info = -2
259 ELSE IF( n.LT.0 ) THEN
260 info = -3
261 ELSE IF( kd.LT.0 ) THEN
262 info = -4
263 ELSE IF( ldab.LT.kd+1 ) THEN
264 info = -6
265 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
266 info = -9
267 END IF
268*
269 IF( info.EQ.0 ) THEN
270 IF( n.LE.1 ) THEN
271 lwmin = 1
272 work( 1 ) = sroundup_lwork(lwmin)
273 ELSE
274 ib = ilaenv2stage( 2, 'SSYTRD_SB2ST', jobz,
275 $ n, kd, -1, -1 )
276 lhtrd = ilaenv2stage( 3, 'SSYTRD_SB2ST', jobz,
277 $ n, kd, ib, -1 )
278 lwtrd = ilaenv2stage( 4, 'SSYTRD_SB2ST', jobz,
279 $ n, kd, ib, -1 )
280 lwmin = n + lhtrd + lwtrd
281 work( 1 ) = sroundup_lwork(lwmin)
282 ENDIF
283*
284 IF( lwork.LT.lwmin .AND. .NOT.lquery )
285 $ info = -11
286 END IF
287*
288 IF( info.NE.0 ) THEN
289 CALL xerbla( 'SSBEV_2STAGE ', -info )
290 RETURN
291 ELSE IF( lquery ) THEN
292 RETURN
293 END IF
294*
295* Quick return if possible
296*
297 IF( n.EQ.0 )
298 $ RETURN
299*
300 IF( n.EQ.1 ) THEN
301 IF( lower ) THEN
302 w( 1 ) = ab( 1, 1 )
303 ELSE
304 w( 1 ) = ab( kd+1, 1 )
305 END IF
306 IF( wantz )
307 $ z( 1, 1 ) = one
308 RETURN
309 END IF
310*
311* Get machine constants.
312*
313 safmin = slamch( 'Safe minimum' )
314 eps = slamch( 'Precision' )
315 smlnum = safmin / eps
316 bignum = one / smlnum
317 rmin = sqrt( smlnum )
318 rmax = sqrt( bignum )
319*
320* Scale matrix to allowable range, if necessary.
321*
322 anrm = slansb( 'M', uplo, n, kd, ab, ldab, work )
323 iscale = 0
324 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
325 iscale = 1
326 sigma = rmin / anrm
327 ELSE IF( anrm.GT.rmax ) THEN
328 iscale = 1
329 sigma = rmax / anrm
330 END IF
331 IF( iscale.EQ.1 ) THEN
332 IF( lower ) THEN
333 CALL slascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
334 ELSE
335 CALL slascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
336 END IF
337 END IF
338*
339* Call SSYTRD_SB2ST to reduce symmetric band matrix to tridiagonal form.
340*
341 inde = 1
342 indhous = inde + n
343 indwrk = indhous + lhtrd
344 llwork = lwork - indwrk + 1
345*
346 CALL ssytrd_sb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
347 $ work( inde ), work( indhous ), lhtrd,
348 $ work( indwrk ), llwork, iinfo )
349*
350* For eigenvalues only, call SSTERF. For eigenvectors, call SSTEQR.
351*
352 IF( .NOT.wantz ) THEN
353 CALL ssterf( n, w, work( inde ), info )
354 ELSE
355 CALL ssteqr( jobz, n, w, work( inde ), z, ldz, work( indwrk ),
356 $ info )
357 END IF
358*
359* If matrix was scaled, then rescale eigenvalues appropriately.
360*
361 IF( iscale.EQ.1 ) THEN
362 IF( info.EQ.0 ) THEN
363 imax = n
364 ELSE
365 imax = info - 1
366 END IF
367 CALL sscal( imax, one / sigma, w, 1 )
368 END IF
369*
370* Set WORK(1) to optimal workspace size.
371*
372 work( 1 ) = sroundup_lwork(lwmin)
373*
374 RETURN
375*
376* End of SSBEV_2STAGE
377*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
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 ssteqr(compz, n, d, e, z, ldz, work, info)
SSTEQR
Definition ssteqr.f:131
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: