LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ chegv_2stage()

subroutine chegv_2stage ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  W,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CHEGV_2STAGE

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

Purpose:
 CHEGV_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
 sizes N>2000.
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 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 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 REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is COMPLEX 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 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:  CPOTRF or CHEEV returned an error code:
             <= N:  if INFO = i, CHEEV 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.
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 230 of file chegv_2stage.f.

232 *
233  IMPLICIT NONE
234 *
235 * -- LAPACK driver routine --
236 * -- LAPACK is a software package provided by Univ. of Tennessee, --
237 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
238 *
239 * .. Scalar Arguments ..
240  CHARACTER JOBZ, UPLO
241  INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
242 * ..
243 * .. Array Arguments ..
244  REAL RWORK( * ), W( * )
245  COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
246 * ..
247 *
248 * =====================================================================
249 *
250 * .. Parameters ..
251  COMPLEX ONE
252  parameter( one = ( 1.0e+0, 0.0e+0 ) )
253 * ..
254 * .. Local Scalars ..
255  LOGICAL LQUERY, UPPER, WANTZ
256  CHARACTER TRANS
257  INTEGER NEIG, LWMIN, LHTRD, LWTRD, KD, IB
258 * ..
259 * .. External Functions ..
260  LOGICAL LSAME
261  INTEGER ILAENV2STAGE
262  EXTERNAL lsame, ilaenv2stage
263 * ..
264 * .. External Subroutines ..
265  EXTERNAL xerbla, chegst, cpotrf, ctrmm, ctrsm,
266  $ cheev_2stage
267 * ..
268 * .. Intrinsic Functions ..
269  INTRINSIC max
270 * ..
271 * .. Executable Statements ..
272 *
273 * Test the input parameters.
274 *
275  wantz = lsame( jobz, 'V' )
276  upper = lsame( uplo, 'U' )
277  lquery = ( lwork.EQ.-1 )
278 *
279  info = 0
280  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
281  info = -1
282  ELSE IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
283  info = -2
284  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
285  info = -3
286  ELSE IF( n.LT.0 ) THEN
287  info = -4
288  ELSE IF( lda.LT.max( 1, n ) ) THEN
289  info = -6
290  ELSE IF( ldb.LT.max( 1, n ) ) THEN
291  info = -8
292  END IF
293 *
294  IF( info.EQ.0 ) THEN
295  kd = ilaenv2stage( 1, 'CHETRD_2STAGE', jobz, n, -1, -1, -1 )
296  ib = ilaenv2stage( 2, 'CHETRD_2STAGE', jobz, n, kd, -1, -1 )
297  lhtrd = ilaenv2stage( 3, 'CHETRD_2STAGE', jobz, n, kd, ib, -1 )
298  lwtrd = ilaenv2stage( 4, 'CHETRD_2STAGE', jobz, n, kd, ib, -1 )
299  lwmin = n + lhtrd + lwtrd
300  work( 1 ) = lwmin
301 *
302  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
303  info = -11
304  END IF
305  END IF
306 *
307  IF( info.NE.0 ) THEN
308  CALL xerbla( 'CHEGV_2STAGE ', -info )
309  RETURN
310  ELSE IF( lquery ) THEN
311  RETURN
312  END IF
313 *
314 * Quick return if possible
315 *
316  IF( n.EQ.0 )
317  $ RETURN
318 *
319 * Form a Cholesky factorization of B.
320 *
321  CALL cpotrf( uplo, n, b, ldb, info )
322  IF( info.NE.0 ) THEN
323  info = n + info
324  RETURN
325  END IF
326 *
327 * Transform problem to standard eigenvalue problem and solve.
328 *
329  CALL chegst( itype, uplo, n, a, lda, b, ldb, info )
330  CALL cheev_2stage( jobz, uplo, n, a, lda, w,
331  $ work, lwork, rwork, info )
332 *
333  IF( wantz ) THEN
334 *
335 * Backtransform eigenvectors to the original problem.
336 *
337  neig = n
338  IF( info.GT.0 )
339  $ neig = info - 1
340  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
341 *
342 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
343 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
344 *
345  IF( upper ) THEN
346  trans = 'N'
347  ELSE
348  trans = 'C'
349  END IF
350 *
351  CALL ctrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
352  $ b, ldb, a, lda )
353 *
354  ELSE IF( itype.EQ.3 ) THEN
355 *
356 * For B*A*x=(lambda)*x;
357 * backtransform eigenvectors: x = L*y or U**H *y
358 *
359  IF( upper ) THEN
360  trans = 'C'
361  ELSE
362  trans = 'N'
363  END IF
364 *
365  CALL ctrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
366  $ b, ldb, a, lda )
367  END IF
368  END IF
369 *
370  work( 1 ) = lwmin
371 *
372  RETURN
373 *
374 * End of CHEGV_2STAGE
375 *
integer function ilaenv2stage(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV2STAGE
Definition: ilaenv2stage.f:149
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine ctrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRMM
Definition: ctrmm.f:177
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:180
subroutine chegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
CHEGST
Definition: chegst.f:128
subroutine cheev_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, INFO)
CHEEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matr...
Definition: cheev_2stage.f:189
subroutine cpotrf(UPLO, N, A, LDA, INFO)
CPOTRF
Definition: cpotrf.f:107
Here is the call graph for this function:
Here is the caller graph for this function: