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

◆ ssygv_2stage()

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
 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 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
                    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 224 of file ssygv_2stage.f.

226*
227 IMPLICIT NONE
228*
229* -- LAPACK driver routine --
230* -- LAPACK is a software package provided by Univ. of Tennessee, --
231* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
232*
233* .. Scalar Arguments ..
234 CHARACTER JOBZ, UPLO
235 INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
236* ..
237* .. Array Arguments ..
238 REAL A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
239* ..
240*
241* =====================================================================
242*
243* .. Parameters ..
244 REAL ONE
245 parameter( one = 1.0e+0 )
246* ..
247* .. Local Scalars ..
248 LOGICAL LQUERY, UPPER, WANTZ
249 CHARACTER TRANS
250 INTEGER NEIG, LWMIN, LHTRD, LWTRD, KD, IB
251* ..
252* .. External Functions ..
253 LOGICAL LSAME
254 INTEGER ILAENV2STAGE
255 REAL SROUNDUP_LWORK
257* ..
258* .. External Subroutines ..
259 EXTERNAL spotrf, ssygst, strmm, strsm, xerbla,
261* ..
262* .. Intrinsic Functions ..
263 INTRINSIC max
264* ..
265* .. Executable Statements ..
266*
267* Test the input parameters.
268*
269 wantz = lsame( jobz, 'V' )
270 upper = lsame( uplo, 'U' )
271 lquery = ( lwork.EQ.-1 )
272*
273 info = 0
274 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
275 info = -1
276 ELSE IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
277 info = -2
278 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
279 info = -3
280 ELSE IF( n.LT.0 ) THEN
281 info = -4
282 ELSE IF( lda.LT.max( 1, n ) ) THEN
283 info = -6
284 ELSE IF( ldb.LT.max( 1, n ) ) THEN
285 info = -8
286 END IF
287*
288 IF( info.EQ.0 ) THEN
289 kd = ilaenv2stage( 1, 'SSYTRD_2STAGE', jobz, n, -1, -1, -1 )
290 ib = ilaenv2stage( 2, 'SSYTRD_2STAGE', jobz, n, kd, -1, -1 )
291 lhtrd = ilaenv2stage( 3, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
292 lwtrd = ilaenv2stage( 4, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
293 lwmin = 2*n + lhtrd + lwtrd
294 work( 1 ) = lwmin
295*
296 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
297 info = -11
298 END IF
299 END IF
300*
301 IF( info.NE.0 ) THEN
302 CALL xerbla( 'SSYGV_2STAGE ', -info )
303 RETURN
304 ELSE IF( lquery ) THEN
305 RETURN
306 END IF
307*
308* Quick return if possible
309*
310 IF( n.EQ.0 )
311 $ RETURN
312*
313* Form a Cholesky factorization of B.
314*
315 CALL spotrf( uplo, n, b, ldb, info )
316 IF( info.NE.0 ) THEN
317 info = n + info
318 RETURN
319 END IF
320*
321* Transform problem to standard eigenvalue problem and solve.
322*
323 CALL ssygst( itype, uplo, n, a, lda, b, ldb, info )
324 CALL ssyev_2stage( jobz, uplo, n, a, lda, w, work, lwork, info )
325*
326 IF( wantz ) THEN
327*
328* Backtransform eigenvectors to the original problem.
329*
330 neig = n
331 IF( info.GT.0 )
332 $ neig = info - 1
333 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
334*
335* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
336* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
337*
338 IF( upper ) THEN
339 trans = 'N'
340 ELSE
341 trans = 'T'
342 END IF
343*
344 CALL strsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
345 $ b, ldb, a, lda )
346*
347 ELSE IF( itype.EQ.3 ) THEN
348*
349* For B*A*x=(lambda)*x;
350* backtransform eigenvectors: x = L*y or U**T*y
351*
352 IF( upper ) THEN
353 trans = 'T'
354 ELSE
355 trans = 'N'
356 END IF
357*
358 CALL strmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
359 $ b, ldb, a, lda )
360 END IF
361 END IF
362*
363 work( 1 ) = sroundup_lwork(lwmin)
364 RETURN
365*
366* End of SSYGV_2STAGE
367*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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 matr...
subroutine ssygst(itype, uplo, n, a, lda, b, ldb, info)
SSYGST
Definition ssygst.f:127
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine spotrf(uplo, n, a, lda, info)
SPOTRF
Definition spotrf.f:107
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine strmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRMM
Definition strmm.f:177
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: