LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ dsyevd_2stage()

subroutine dsyevd_2stage ( character  JOBZ,
character  UPLO,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  W,
double precision, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

DSYEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices

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

Purpose:
 DSYEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a
 real symmetric matrix A using the 2stage technique for
 the reduction to tridiagonal. If eigenvectors are desired, it uses a
 divide and conquer algorithm.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
Parameters
[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 triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  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
          orthonormal eigenvectors of the matrix A.
          If JOBZ = 'N', then on exit the lower triangle (if UPLO='L')
          or the upper triangle (if UPLO='U') of A, including the
          diagonal, is destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= 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 (LWORK)
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK must be at least 1.
          If JOBZ = 'N' and N > 1, LWORK must be queried.
                                   LWORK = MAX(1, dimension) where
                                   dimension = max(stage1,stage2) + (KD+1)*N + 2*N+1
                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                                               + max(2*KD*KD, KD*NTHREADS) 
                                               + (KD+1)*N + 2*N+1
                                   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 at least
                                                1 + 6*N + 2*N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK or
          LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If N <= 1,                LIWORK must be at least 1.
          If JOBZ  = 'N' and N > 1, LIWORK must be at least 1.
          If JOBZ  = 'V' and N > 1, LIWORK must be at least 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK or LIWORK 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:  if INFO = i and JOBZ = 'N', then the algorithm failed
                to converge; i off-diagonal elements of an intermediate
                tridiagonal form did not converge to zero;
                if INFO = i and JOBZ = 'V', then the algorithm failed
                to compute an eigenvalue while working on the submatrix
                lying in rows and columns INFO/(N+1) through
                mod(INFO,N+1).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2017
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA
Modified by Francoise Tisseur, University of Tennessee
Modified description of INFO. Sven, 16 Feb 05.
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 229 of file dsyevd_2stage.f.

229 *
230  IMPLICIT NONE
231 *
232 * -- LAPACK driver routine (version 3.8.0) --
233 * -- LAPACK is a software package provided by Univ. of Tennessee, --
234 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
235 * November 2017
236 *
237 * .. Scalar Arguments ..
238  CHARACTER JOBZ, UPLO
239  INTEGER INFO, LDA, LIWORK, LWORK, N
240 * ..
241 * .. Array Arguments ..
242  INTEGER IWORK( * )
243  DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
244 * ..
245 *
246 * =====================================================================
247 *
248 * .. Parameters ..
249  DOUBLE PRECISION ZERO, ONE
250  parameter( zero = 0.0d+0, one = 1.0d+0 )
251 * ..
252 * .. Local Scalars ..
253 *
254  LOGICAL LOWER, LQUERY, WANTZ
255  INTEGER IINFO, INDE, INDTAU, INDWK2, INDWRK, ISCALE,
256  $ LIWMIN, LLWORK, LLWRK2, LWMIN,
257  $ LHTRD, LWTRD, KD, IB, INDHOUS
258  DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
259  $ SMLNUM
260 * ..
261 * .. External Functions ..
262  LOGICAL LSAME
263  INTEGER ILAENV2STAGE
264  DOUBLE PRECISION DLAMCH, DLANSY
265  EXTERNAL lsame, dlamch, dlansy, ilaenv2stage
266 * ..
267 * .. External Subroutines ..
268  EXTERNAL dlacpy, dlascl, dormtr, dscal, dstedc, dsterf,
270 * ..
271 * .. Intrinsic Functions ..
272  INTRINSIC max, sqrt
273 * ..
274 * .. Executable Statements ..
275 *
276 * Test the input parameters.
277 *
278  wantz = lsame( jobz, 'V' )
279  lower = lsame( uplo, 'L' )
280  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
281 *
282  info = 0
283  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
284  info = -1
285  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
286  info = -2
287  ELSE IF( n.LT.0 ) THEN
288  info = -3
289  ELSE IF( lda.LT.max( 1, n ) ) THEN
290  info = -5
291  END IF
292 *
293  IF( info.EQ.0 ) THEN
294  IF( n.LE.1 ) THEN
295  liwmin = 1
296  lwmin = 1
297  ELSE
298  kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz,
299  $ n, -1, -1, -1 )
300  ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz,
301  $ n, kd, -1, -1 )
302  lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz,
303  $ n, kd, ib, -1 )
304  lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz,
305  $ n, kd, ib, -1 )
306  IF( wantz ) THEN
307  liwmin = 3 + 5*n
308  lwmin = 1 + 6*n + 2*n**2
309  ELSE
310  liwmin = 1
311  lwmin = 2*n + 1 + lhtrd + lwtrd
312  END IF
313  END IF
314  work( 1 ) = lwmin
315  iwork( 1 ) = liwmin
316 *
317  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
318  info = -8
319  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
320  info = -10
321  END IF
322  END IF
323 *
324  IF( info.NE.0 ) THEN
325  CALL xerbla( 'DSYEVD_2STAGE', -info )
326  RETURN
327  ELSE IF( lquery ) THEN
328  RETURN
329  END IF
330 *
331 * Quick return if possible
332 *
333  IF( n.EQ.0 )
334  $ RETURN
335 *
336  IF( n.EQ.1 ) THEN
337  w( 1 ) = a( 1, 1 )
338  IF( wantz )
339  $ a( 1, 1 ) = one
340  RETURN
341  END IF
342 *
343 * Get machine constants.
344 *
345  safmin = dlamch( 'Safe minimum' )
346  eps = dlamch( 'Precision' )
347  smlnum = safmin / eps
348  bignum = one / smlnum
349  rmin = sqrt( smlnum )
350  rmax = sqrt( bignum )
351 *
352 * Scale matrix to allowable range, if necessary.
353 *
354  anrm = dlansy( 'M', uplo, n, a, lda, work )
355  iscale = 0
356  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
357  iscale = 1
358  sigma = rmin / anrm
359  ELSE IF( anrm.GT.rmax ) THEN
360  iscale = 1
361  sigma = rmax / anrm
362  END IF
363  IF( iscale.EQ.1 )
364  $ CALL dlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
365 *
366 * Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
367 *
368  inde = 1
369  indtau = inde + n
370  indhous = indtau + n
371  indwrk = indhous + lhtrd
372  llwork = lwork - indwrk + 1
373  indwk2 = indwrk + n*n
374  llwrk2 = lwork - indwk2 + 1
375 *
376  CALL dsytrd_2stage( jobz, uplo, n, a, lda, w, work( inde ),
377  $ work( indtau ), work( indhous ), lhtrd,
378  $ work( indwrk ), llwork, iinfo )
379 *
380 * For eigenvalues only, call DSTERF. For eigenvectors, first call
381 * DSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
382 * tridiagonal matrix, then call DORMTR to multiply it by the
383 * Householder transformations stored in A.
384 *
385  IF( .NOT.wantz ) THEN
386  CALL dsterf( n, w, work( inde ), info )
387  ELSE
388 * Not available in this release, and argument checking should not
389 * let it getting here
390  RETURN
391  CALL dstedc( 'I', n, w, work( inde ), work( indwrk ), n,
392  $ work( indwk2 ), llwrk2, iwork, liwork, info )
393  CALL dormtr( 'L', uplo, 'N', n, n, a, lda, work( indtau ),
394  $ work( indwrk ), n, work( indwk2 ), llwrk2, iinfo )
395  CALL dlacpy( 'A', n, n, work( indwrk ), n, a, lda )
396  END IF
397 *
398 * If matrix was scaled, then rescale eigenvalues appropriately.
399 *
400  IF( iscale.EQ.1 )
401  $ CALL dscal( n, one / sigma, w, 1 )
402 *
403  work( 1 ) = lwmin
404  iwork( 1 ) = liwmin
405 *
406  RETURN
407 *
408 * End of DSYEVD_2STAGE
409 *
Here is the call graph for this function:
Here is the caller graph for this function:
dlascl
subroutine dlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: dlascl.f:145
dlansy
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition: dlansy.f:124
dsterf
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
dsytrd_2stage
subroutine dsytrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
DSYTRD_2STAGE
Definition: dsytrd_2stage.f:226
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
lsame
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
dstedc
subroutine dstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEDC
Definition: dstedc.f:190
dormtr
subroutine dormtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMTR
Definition: dormtr.f:173
ilaenv2stage
integer function ilaenv2stage(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV2STAGE
Definition: ilaenv2stage.f:151
dscal
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:81
dlamch
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:70
dlacpy
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105