LAPACK  3.10.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.
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 225 of file dsyevd_2stage.f.

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