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

SSYGV_2STAGE

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

Purpose:
 SSYGV_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 REAL 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 REAL 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 REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is REAL 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:  SPOTRF or SSYEV returned an error code:
             <= N:  if INFO = i, SSYEV 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 ssygv_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  REAL a( lda, * ), b( ldb, * ), w( * ), work( * )
242 * ..
243 *
244 * =====================================================================
245 *
246 * .. Parameters ..
247  REAL one
248  parameter ( one = 1.0e+0 )
249 * ..
250 * .. Local Scalars ..
251  LOGICAL lquery, upper, wantz
252  CHARACTER trans
253  INTEGER neig,
254  $ llwork, lwmin, lhtrd, lwtrd, kd, ib, indhous
255 * ..
256 * .. External Functions ..
257  LOGICAL lsame
258  INTEGER ilaenv
259  EXTERNAL lsame, ilaenv
260 * ..
261 * .. External Subroutines ..
262  EXTERNAL spotrf, ssygst, strmm, strsm, xerbla,
263  $ ssyev_2stage
264 * ..
265 * .. Intrinsic Functions ..
266  INTRINSIC max
267 * ..
268 * .. Executable Statements ..
269 *
270 * Test the input parameters.
271 *
272  wantz = lsame( jobz, 'V' )
273  upper = lsame( uplo, 'U' )
274  lquery = ( lwork.EQ.-1 )
275 *
276  info = 0
277  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
278  info = -1
279  ELSE IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
280  info = -2
281  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
282  info = -3
283  ELSE IF( n.LT.0 ) THEN
284  info = -4
285  ELSE IF( lda.LT.max( 1, n ) ) THEN
286  info = -6
287  ELSE IF( ldb.LT.max( 1, n ) ) THEN
288  info = -8
289  END IF
290 *
291  IF( info.EQ.0 ) THEN
292  kd = ilaenv( 17, 'SSYTRD_2STAGE', jobz, n, -1, -1, -1 )
293  ib = ilaenv( 18, 'SSYTRD_2STAGE', jobz, n, kd, -1, -1 )
294  lhtrd = ilaenv( 19, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
295  lwtrd = ilaenv( 20, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
296  lwmin = 2*n + lhtrd + lwtrd
297  work( 1 ) = lwmin
298 *
299  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
300  info = -11
301  END IF
302  END IF
303 *
304  IF( info.NE.0 ) THEN
305  CALL xerbla( 'SSYGV_2STAGE ', -info )
306  RETURN
307  ELSE IF( lquery ) THEN
308  RETURN
309  END IF
310 *
311 * Quick return if possible
312 *
313  IF( n.EQ.0 )
314  $ RETURN
315 *
316 * Form a Cholesky factorization of B.
317 *
318  CALL spotrf( uplo, n, b, ldb, info )
319  IF( info.NE.0 ) THEN
320  info = n + info
321  RETURN
322  END IF
323 *
324 * Transform problem to standard eigenvalue problem and solve.
325 *
326  CALL ssygst( itype, uplo, n, a, lda, b, ldb, info )
327  CALL ssyev_2stage( jobz, uplo, n, a, lda, w, work, lwork, info )
328 *
329  IF( wantz ) THEN
330 *
331 * Backtransform eigenvectors to the original problem.
332 *
333  neig = n
334  IF( info.GT.0 )
335  $ neig = info - 1
336  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
337 *
338 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
339 * backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
340 *
341  IF( upper ) THEN
342  trans = 'N'
343  ELSE
344  trans = 'T'
345  END IF
346 *
347  CALL strsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
348  $ b, ldb, a, lda )
349 *
350  ELSE IF( itype.EQ.3 ) THEN
351 *
352 * For B*A*x=(lambda)*x;
353 * backtransform eigenvectors: x = L*y or U**T*y
354 *
355  IF( upper ) THEN
356  trans = 'T'
357  ELSE
358  trans = 'N'
359  END IF
360 *
361  CALL strmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
362  $ b, ldb, a, lda )
363  END IF
364  END IF
365 *
366  work( 1 ) = lwmin
367  RETURN
368 *
369 * End of SSYGV_2STAGE
370 *
subroutine ssyev_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO)
SSYEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...
Definition: ssyev_2stage.f:185
subroutine strsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
STRSM
Definition: strsm.f:183
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:109
subroutine strmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
STRMM
Definition: strmm.f:179
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine ssygst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
SSYGST
Definition: ssygst.f:129

Here is the call graph for this function:

Here is the caller graph for this function: