LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ chbev_2stage()

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
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 213 of file chbev_2stage.f.

213 *
214  IMPLICIT NONE
215 *
216 * -- LAPACK driver routine (version 3.8.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 * November 2017
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 ILAENV2STAGE
246  REAL SLAMCH, CLANHB
247  EXTERNAL lsame, slamch, clanhb, ilaenv2stage
248 * ..
249 * .. External Subroutines ..
250  EXTERNAL sscal, ssterf, xerbla, clascl, csteqr,
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 = ilaenv2stage( 2, 'CHETRD_HB2ST', jobz,
285  $ n, kd, -1, -1 )
286  lhtrd = ilaenv2stage( 3, 'CHETRD_HB2ST', jobz,
287  $ n, kd, ib, -1 )
288  lwtrd = ilaenv2stage( 4, 'CHETRD_HB2ST', jobz,
289  $ n, kd, ib, -1 )
290  lwmin = lhtrd + lwtrd
291  work( 1 ) = lwmin
292  ENDIF
293 *
294  IF( lwork.LT.lwmin .AND. .NOT.lquery )
295  $ info = -11
296  END IF
297 *
298  IF( info.NE.0 ) THEN
299  CALL xerbla( 'CHBEV_2STAGE ', -info )
300  RETURN
301  ELSE IF( lquery ) THEN
302  RETURN
303  END IF
304 *
305 * Quick return if possible
306 *
307  IF( n.EQ.0 )
308  $ RETURN
309 *
310  IF( n.EQ.1 ) THEN
311  IF( lower ) THEN
312  w( 1 ) = real( ab( 1, 1 ) )
313  ELSE
314  w( 1 ) = real( ab( kd+1, 1 ) )
315  END IF
316  IF( wantz )
317  $ z( 1, 1 ) = one
318  RETURN
319  END IF
320 *
321 * Get machine constants.
322 *
323  safmin = slamch( 'Safe minimum' )
324  eps = slamch( 'Precision' )
325  smlnum = safmin / eps
326  bignum = one / smlnum
327  rmin = sqrt( smlnum )
328  rmax = sqrt( bignum )
329 *
330 * Scale matrix to allowable range, if necessary.
331 *
332  anrm = clanhb( 'M', uplo, n, kd, ab, ldab, rwork )
333  iscale = 0
334  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
335  iscale = 1
336  sigma = rmin / anrm
337  ELSE IF( anrm.GT.rmax ) THEN
338  iscale = 1
339  sigma = rmax / anrm
340  END IF
341  IF( iscale.EQ.1 ) THEN
342  IF( lower ) THEN
343  CALL clascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
344  ELSE
345  CALL clascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
346  END IF
347  END IF
348 *
349 * Call CHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form.
350 *
351  inde = 1
352  indhous = 1
353  indwrk = indhous + lhtrd
354  llwork = lwork - indwrk + 1
355 *
356  CALL chetrd_hb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
357  $ rwork( inde ), work( indhous ), lhtrd,
358  $ work( indwrk ), llwork, iinfo )
359 *
360 * For eigenvalues only, call SSTERF. For eigenvectors, call CSTEQR.
361 *
362  IF( .NOT.wantz ) THEN
363  CALL ssterf( n, w, rwork( inde ), info )
364  ELSE
365  indrwk = inde + n
366  CALL csteqr( jobz, n, w, rwork( inde ), z, ldz,
367  $ rwork( indrwk ), info )
368  END IF
369 *
370 * If matrix was scaled, then rescale eigenvalues appropriately.
371 *
372  IF( iscale.EQ.1 ) THEN
373  IF( info.EQ.0 ) THEN
374  imax = n
375  ELSE
376  imax = info - 1
377  END IF
378  CALL sscal( imax, one / sigma, w, 1 )
379  END IF
380 *
381 * Set WORK(1) to optimal workspace size.
382 *
383  work( 1 ) = lwmin
384 *
385  RETURN
386 *
387 * End of CHBEV_2STAGE
388 *
Here is the call graph for this function:
Here is the caller graph for this function:
csteqr
subroutine csteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
CSTEQR
Definition: csteqr.f:134
clanhb
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,...
Definition: clanhb.f:134
ssterf
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
sscal
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:81
chetrd_2stage
subroutine chetrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
CHETRD_2STAGE
Definition: chetrd_2stage.f:226
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
chetrd_hb2st
subroutine chetrd_hb2st(STAGE1, VECT, UPLO, N, KD, AB, LDAB, D, E, HOUS, LHOUS, WORK, LWORK, INFO)
CHBTRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
Definition: chetrd_hb2st.F:232
ilaenv2stage
integer function ilaenv2stage(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV2STAGE
Definition: ilaenv2stage.f:151
clascl
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