LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine zhegv_2stage ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( * )  W,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  INFO 
)

ZHEGV_2STAGE

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

Purpose:
 ZHEGV_2STAGE computes all the eigenvalues, and optionally, the eigenvectors
 of a complex generalized Hermitian-definite eigenproblem, of the form
 A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
 Here A and B are assumed to be Hermitian and B is also
 positive definite.
 This routine use the 2stage technique for the reduction to tridiagonal
 which showed higher performance on recent architecture and for large  
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the problem type to be solved:
          = 1:  A*x = (lambda)*B*x
          = 2:  A*B*x = (lambda)*x
          = 3:  B*A*x = (lambda)*x
[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 triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA, N)
          On entry, the Hermitian matrix A.  If UPLO = 'U', the
          leading N-by-N upper triangular part of A contains the
          upper triangular part of the matrix A.  If UPLO = 'L',
          the leading N-by-N lower triangular part of A contains
          the lower triangular part of the matrix A.

          On exit, if JOBZ = 'V', then if INFO = 0, A contains the
          matrix Z of eigenvectors.  The eigenvectors are normalized
          as follows:
          if ITYPE = 1 or 2, Z**H*B*Z = I;
          if ITYPE = 3, Z**H*inv(B)*Z = I.
          If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
          or the lower triangle (if UPLO='L') of A, including the
          diagonal, is destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX*16 array, dimension (LDB, N)
          On entry, the Hermitian positive definite matrix B.
          If UPLO = 'U', the leading N-by-N upper triangular part of B
          contains the upper triangular part of the matrix B.
          If UPLO = 'L', the leading N-by-N lower triangular part of B
          contains the lower triangular part of the matrix B.

          On exit, if INFO <= N, the part of B containing the matrix is
          overwritten by the triangular factor U or L from the Cholesky
          factorization B = U**H*U or B = L*L**H.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is COMPLEX*16 array, dimension (MAX(1,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 = max(stage1,stage2) + (KD+1)*N + N
                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                                               + max(2*KD*KD, KD*NTHREADS) 
                                               + (KD+1)*N + N
                                   where KD is the blocking size of the reduction,
                                   FACTOPTNB is the blocking used by the QR or LQ
                                   algorithm, usually FACTOPTNB=128 is a good choice
                                   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]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:  ZPOTRF or ZHEEV returned an error code:
             <= N:  if INFO = i, ZHEEV failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not converge to zero;
             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
                    minor of order i of B is not positive definite.
                    The factorization of B could not be completed and
                    no eigenvalues or eigenvectors were computed.
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 234 of file zhegv_2stage.f.

234 *
235  IMPLICIT NONE
236 *
237 * -- LAPACK driver routine (version 3.7.0) --
238 * -- LAPACK is a software package provided by Univ. of Tennessee, --
239 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
240 * December 2016
241 *
242 * .. Scalar Arguments ..
243  CHARACTER jobz, uplo
244  INTEGER info, itype, lda, ldb, lwork, n
245 * ..
246 * .. Array Arguments ..
247  DOUBLE PRECISION rwork( * ), w( * )
248  COMPLEX*16 a( lda, * ), b( ldb, * ), work( * )
249 * ..
250 *
251 * =====================================================================
252 *
253 * .. Parameters ..
254  COMPLEX*16 one
255  parameter ( one = ( 1.0d+0, 0.0d+0 ) )
256 * ..
257 * .. Local Scalars ..
258  LOGICAL lquery, upper, wantz
259  CHARACTER trans
260  INTEGER neig, lwmin, lhtrd, lwtrd, kd, ib
261 * ..
262 * .. External Functions ..
263  LOGICAL lsame
264  INTEGER ilaenv
265  EXTERNAL lsame, ilaenv
266 * ..
267 * .. External Subroutines ..
268  EXTERNAL xerbla, zhegst, zpotrf, ztrmm, ztrsm,
269  $ zheev_2stage
270 * ..
271 * .. Intrinsic Functions ..
272  INTRINSIC max
273 * ..
274 * .. Executable Statements ..
275 *
276 * Test the input parameters.
277 *
278  wantz = lsame( jobz, 'V' )
279  upper = lsame( uplo, 'U' )
280  lquery = ( lwork.EQ.-1 )
281 *
282  info = 0
283  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
284  info = -1
285  ELSE IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
286  info = -2
287  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
288  info = -3
289  ELSE IF( n.LT.0 ) THEN
290  info = -4
291  ELSE IF( lda.LT.max( 1, n ) ) THEN
292  info = -6
293  ELSE IF( ldb.LT.max( 1, n ) ) THEN
294  info = -8
295  END IF
296 *
297  IF( info.EQ.0 ) THEN
298  kd = ilaenv( 17, 'ZHETRD_2STAGE', jobz, n, -1, -1, -1 )
299  ib = ilaenv( 18, 'ZHETRD_2STAGE', jobz, n, kd, -1, -1 )
300  lhtrd = ilaenv( 19, 'ZHETRD_2STAGE', jobz, n, kd, ib, -1 )
301  lwtrd = ilaenv( 20, 'ZHETRD_2STAGE', jobz, n, kd, ib, -1 )
302  lwmin = n + lhtrd + lwtrd
303  work( 1 ) = lwmin
304 *
305  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
306  info = -11
307  END IF
308  END IF
309 *
310  IF( info.NE.0 ) THEN
311  CALL xerbla( 'ZHEGV_2STAGE ', -info )
312  RETURN
313  ELSE IF( lquery ) THEN
314  RETURN
315  END IF
316 *
317 * Quick return if possible
318 *
319  IF( n.EQ.0 )
320  $ RETURN
321 *
322 * Form a Cholesky factorization of B.
323 *
324  CALL zpotrf( uplo, n, b, ldb, info )
325  IF( info.NE.0 ) THEN
326  info = n + info
327  RETURN
328  END IF
329 *
330 * Transform problem to standard eigenvalue problem and solve.
331 *
332  CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
333  CALL zheev_2stage( jobz, uplo, n, a, lda, w,
334  $ work, lwork, rwork, info )
335 *
336  IF( wantz ) THEN
337 *
338 * Backtransform eigenvectors to the original problem.
339 *
340  neig = n
341  IF( info.GT.0 )
342  $ neig = info - 1
343  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
344 *
345 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
346 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
347 *
348  IF( upper ) THEN
349  trans = 'N'
350  ELSE
351  trans = 'C'
352  END IF
353 *
354  CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
355  $ b, ldb, a, lda )
356 *
357  ELSE IF( itype.EQ.3 ) THEN
358 *
359 * For B*A*x=(lambda)*x;
360 * backtransform eigenvectors: x = L*y or U**H *y
361 *
362  IF( upper ) THEN
363  trans = 'C'
364  ELSE
365  trans = 'N'
366  END IF
367 *
368  CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
369  $ b, ldb, a, lda )
370  END IF
371  END IF
372 *
373  work( 1 ) = lwmin
374 *
375  RETURN
376 *
377 * End of ZHEGV_2STAGE
378 *
subroutine zhegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
ZHEGST
Definition: zhegst.f:129
subroutine zpotrf(UPLO, N, A, LDA, INFO)
ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
Definition: zpotrf.f:102
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine ztrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRMM
Definition: ztrmm.f:179
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zheev_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, INFO)
ZHEEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE mat...
Definition: zheev_2stage.f:191
subroutine ztrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRSM
Definition: ztrsm.f:182

Here is the call graph for this function:

Here is the caller graph for this function: