LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine chbev_2stage ( character  JOBZ,
character  UPLO,
integer  N,
integer  KD,
complex, dimension( ldab, * )  AB,
integer  LDAB,
real, dimension( * )  W,
complex, dimension( ldz, * )  Z,
integer  LDZ,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  INFO 
)

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

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

Purpose:
 CHBEV_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
 a complex Hermitian 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 COMPLEX array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the Hermitian 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 COMPLEX 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 COMPLEX 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
                                   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 sizes of the WORK, RWORK and
          IWORK arrays, returns these values as the first entries of
          the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]RWORK
          RWORK 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
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 213 of file chbev_2stage.f.

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