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

◆ zhbev_2stage()

subroutine zhbev_2stage ( character  jobz,
character  uplo,
integer  n,
integer  kd,
complex*16, dimension( ldab, * )  ab,
integer  ldab,
double precision, dimension( * )  w,
complex*16, dimension( ldz, * )  z,
integer  ldz,
complex*16, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  rwork,
integer  info 
)

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

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

Purpose:
 ZHBEV_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*16 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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 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*16 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 DOUBLE PRECISION 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.
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 209 of file zhbev_2stage.f.

211*
212 IMPLICIT NONE
213*
214* -- LAPACK driver routine --
215* -- LAPACK is a software package provided by Univ. of Tennessee, --
216* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
217*
218* .. Scalar Arguments ..
219 CHARACTER JOBZ, UPLO
220 INTEGER INFO, KD, LDAB, LDZ, N, LWORK
221* ..
222* .. Array Arguments ..
223 DOUBLE PRECISION RWORK( * ), W( * )
224 COMPLEX*16 AB( LDAB, * ), WORK( * ), Z( LDZ, * )
225* ..
226*
227* =====================================================================
228*
229* .. Parameters ..
230 DOUBLE PRECISION ZERO, ONE
231 parameter( zero = 0.0d0, one = 1.0d0 )
232* ..
233* .. Local Scalars ..
234 LOGICAL LOWER, WANTZ, LQUERY
235 INTEGER IINFO, IMAX, INDE, INDWRK, INDRWK, ISCALE,
236 $ LLWORK, LWMIN, LHTRD, LWTRD, IB, INDHOUS
237 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
238 $ SMLNUM
239* ..
240* .. External Functions ..
241 LOGICAL LSAME
242 INTEGER ILAENV2STAGE
243 DOUBLE PRECISION DLAMCH, ZLANHB
244 EXTERNAL lsame, dlamch, zlanhb, ilaenv2stage
245* ..
246* .. External Subroutines ..
247 EXTERNAL dscal, dsterf, xerbla, zlascl, zsteqr,
249* ..
250* .. Intrinsic Functions ..
251 INTRINSIC dble, sqrt
252* ..
253* .. Executable Statements ..
254*
255* Test the input parameters.
256*
257 wantz = lsame( jobz, 'V' )
258 lower = lsame( uplo, 'L' )
259 lquery = ( lwork.EQ.-1 )
260*
261 info = 0
262 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
263 info = -1
264 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
265 info = -2
266 ELSE IF( n.LT.0 ) THEN
267 info = -3
268 ELSE IF( kd.LT.0 ) THEN
269 info = -4
270 ELSE IF( ldab.LT.kd+1 ) THEN
271 info = -6
272 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
273 info = -9
274 END IF
275*
276 IF( info.EQ.0 ) THEN
277 IF( n.LE.1 ) THEN
278 lwmin = 1
279 work( 1 ) = lwmin
280 ELSE
281 ib = ilaenv2stage( 2, 'ZHETRD_HB2ST', jobz,
282 $ n, kd, -1, -1 )
283 lhtrd = ilaenv2stage( 3, 'ZHETRD_HB2ST', jobz,
284 $ n, kd, ib, -1 )
285 lwtrd = ilaenv2stage( 4, 'ZHETRD_HB2ST', jobz,
286 $ 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( 'ZHBEV_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 ) = dble( ab( 1, 1 ) )
310 ELSE
311 w( 1 ) = dble( 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 = dlamch( 'Safe minimum' )
321 eps = dlamch( '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 = zlanhb( '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 zlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
341 ELSE
342 CALL zlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
343 END IF
344 END IF
345*
346* Call ZHBTRD_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 zhetrd_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 DSTERF. For eigenvectors, call ZSTEQR.
358*
359 IF( .NOT.wantz ) THEN
360 CALL dsterf( n, w, rwork( inde ), info )
361 ELSE
362 indrwk = inde + n
363 CALL zsteqr( 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 dscal( 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 ZHBEV_2STAGE
385*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhetrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
ZHETRD_2STAGE
subroutine zhetrd_hb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
ZHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlanhb(norm, uplo, n, k, ab, ldab, work)
ZLANHB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhb.f:132
subroutine zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition zlascl.f:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
Definition zsteqr.f:132
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
Here is the call graph for this function:
Here is the caller graph for this function: