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

◆ dsygv_2stage()

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
!> 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 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
!>                    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 222 of file dsygv_2stage.f.

225*
226 IMPLICIT NONE
227*
228* -- LAPACK driver routine --
229* -- LAPACK is a software package provided by Univ. of Tennessee, --
230* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
231*
232* .. Scalar Arguments ..
233 CHARACTER JOBZ, UPLO
234 INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
235* ..
236* .. Array Arguments ..
237 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
238* ..
239*
240* =====================================================================
241*
242* .. Parameters ..
243 DOUBLE PRECISION ONE
244 parameter( one = 1.0d+0 )
245* ..
246* .. Local Scalars ..
247 LOGICAL LQUERY, UPPER, WANTZ
248 CHARACTER TRANS
249 INTEGER NEIG, LWMIN, LHTRD, LWTRD, KD, IB
250* ..
251* .. External Functions ..
252 LOGICAL LSAME
253 INTEGER ILAENV2STAGE
254 EXTERNAL lsame, ilaenv2stage
255* ..
256* .. External Subroutines ..
257 EXTERNAL dpotrf, dsygst, dtrmm, dtrsm,
258 $ xerbla,
260* ..
261* .. Intrinsic Functions ..
262 INTRINSIC max
263* ..
264* .. Executable Statements ..
265*
266* Test the input parameters.
267*
268 wantz = lsame( jobz, 'V' )
269 upper = lsame( uplo, 'U' )
270 lquery = ( lwork.EQ.-1 )
271*
272 info = 0
273 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
274 info = -1
275 ELSE IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
276 info = -2
277 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
278 info = -3
279 ELSE IF( n.LT.0 ) THEN
280 info = -4
281 ELSE IF( lda.LT.max( 1, n ) ) THEN
282 info = -6
283 ELSE IF( ldb.LT.max( 1, n ) ) THEN
284 info = -8
285 END IF
286*
287 IF( info.EQ.0 ) THEN
288 kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz, n, -1, -1,
289 $ -1 )
290 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1,
291 $ -1 )
292 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib,
293 $ -1 )
294 lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz, n, kd, ib,
295 $ -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( 'DSYGV_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 dpotrf( 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 dsygst( itype, uplo, n, a, lda, b, ldb, info )
327 CALL dsyev_2stage( jobz, uplo, n, a, lda, w, work, lwork,
328 $ info )
329*
330 IF( wantz ) THEN
331*
332* Backtransform eigenvectors to the original problem.
333*
334 neig = n
335 IF( info.GT.0 )
336 $ neig = info - 1
337 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
338*
339* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
340* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
341*
342 IF( upper ) THEN
343 trans = 'N'
344 ELSE
345 trans = 'T'
346 END IF
347*
348 CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, neig,
349 $ one,
350 $ b, ldb, a, lda )
351*
352 ELSE IF( itype.EQ.3 ) THEN
353*
354* For B*A*x=(lambda)*x;
355* backtransform eigenvectors: x = L*y or U**T*y
356*
357 IF( upper ) THEN
358 trans = 'T'
359 ELSE
360 trans = 'N'
361 END IF
362*
363 CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, neig,
364 $ one,
365 $ b, ldb, a, lda )
366 END IF
367 END IF
368*
369 work( 1 ) = lwmin
370 RETURN
371*
372* End of DSYGV_2STAGE
373*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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 matr...
subroutine dsygst(itype, uplo, n, a, lda, b, ldb, info)
DSYGST
Definition dsygst.f:125
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:105
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: