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

◆ zhegv_2stage()

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
!> 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*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
!>                    principal minor of order i of B is not positive.
!>                    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 228 of file zhegv_2stage.f.

231*
232 IMPLICIT NONE
233*
234* -- LAPACK driver routine --
235* -- LAPACK is a software package provided by Univ. of Tennessee, --
236* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
237*
238* .. Scalar Arguments ..
239 CHARACTER JOBZ, UPLO
240 INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
241* ..
242* .. Array Arguments ..
243 DOUBLE PRECISION RWORK( * ), W( * )
244 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
245* ..
246*
247* =====================================================================
248*
249* .. Parameters ..
250 COMPLEX*16 ONE
251 parameter( one = ( 1.0d+0, 0.0d+0 ) )
252* ..
253* .. Local Scalars ..
254 LOGICAL LQUERY, UPPER, WANTZ
255 CHARACTER TRANS
256 INTEGER NEIG, LWMIN, LHTRD, LWTRD, KD, IB
257* ..
258* .. External Functions ..
259 LOGICAL LSAME
260 INTEGER ILAENV2STAGE
261 EXTERNAL lsame, ilaenv2stage
262* ..
263* .. External Subroutines ..
264 EXTERNAL xerbla, zhegst, zpotrf, ztrmm,
265 $ ztrsm,
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, 'ZHETRD_2STAGE', jobz, n, -1, -1,
296 $ -1 )
297 ib = ilaenv2stage( 2, 'ZHETRD_2STAGE', jobz, n, kd, -1,
298 $ -1 )
299 lhtrd = ilaenv2stage( 3, 'ZHETRD_2STAGE', jobz, n, kd, ib,
300 $ -1 )
301 lwtrd = ilaenv2stage( 4, 'ZHETRD_2STAGE', jobz, n, kd, ib,
302 $ -1 )
303 lwmin = n + lhtrd + lwtrd
304 work( 1 ) = lwmin
305*
306 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
307 info = -11
308 END IF
309 END IF
310*
311 IF( info.NE.0 ) THEN
312 CALL xerbla( 'ZHEGV_2STAGE ', -info )
313 RETURN
314 ELSE IF( lquery ) THEN
315 RETURN
316 END IF
317*
318* Quick return if possible
319*
320 IF( n.EQ.0 )
321 $ RETURN
322*
323* Form a Cholesky factorization of B.
324*
325 CALL zpotrf( uplo, n, b, ldb, info )
326 IF( info.NE.0 ) THEN
327 info = n + info
328 RETURN
329 END IF
330*
331* Transform problem to standard eigenvalue problem and solve.
332*
333 CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
334 CALL zheev_2stage( jobz, uplo, n, a, lda, w,
335 $ work, lwork, rwork, info )
336*
337 IF( wantz ) THEN
338*
339* Backtransform eigenvectors to the original problem.
340*
341 neig = n
342 IF( info.GT.0 )
343 $ neig = info - 1
344 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
345*
346* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
347* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
348*
349 IF( upper ) THEN
350 trans = 'N'
351 ELSE
352 trans = 'C'
353 END IF
354*
355 CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, neig,
356 $ one,
357 $ b, ldb, a, lda )
358*
359 ELSE IF( itype.EQ.3 ) THEN
360*
361* For B*A*x=(lambda)*x;
362* backtransform eigenvectors: x = L*y or U**H *y
363*
364 IF( upper ) THEN
365 trans = 'C'
366 ELSE
367 trans = 'N'
368 END IF
369*
370 CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, neig,
371 $ one,
372 $ b, ldb, a, lda )
373 END IF
374 END IF
375*
376 work( 1 ) = lwmin
377*
378 RETURN
379*
380* End of ZHEGV_2STAGE
381*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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 matr...
subroutine zhegst(itype, uplo, n, a, lda, b, ldb, info)
ZHEGST
Definition zhegst.f:126
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:105
subroutine ztrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRMM
Definition ztrmm.f:177
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180
Here is the call graph for this function:
Here is the caller graph for this function: