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

DSYGV_2STAGE

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

Purpose:
 DSYGV_2STAGE computes all the eigenvalues, and optionally, the eigenvectors
 of a real generalized symmetric-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 symmetric 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 DOUBLE PRECISION array, dimension (LDA, N)
          On entry, the symmetric 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**T*B*Z = I;
          if ITYPE = 3, Z**T*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 DOUBLE PRECISION array, dimension (LDB, N)
          On entry, the symmetric 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**T*U or B = L*L**T.
[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 DOUBLE PRECISION 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 + 2*N
                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                                               + max(2*KD*KD, KD*NTHREADS) 
                                               + (KD+1)*N + 2*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]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  DPOTRF or DSYEV returned an error code:
             <= N:  if INFO = i, DSYEV 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 228 of file dsygv_2stage.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: