LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine dsyev_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  INFO 
)

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

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

Purpose:
 DSYEV_2STAGE computes all eigenvalues and, optionally, eigenvectors of a
 real symmetric matrix A using the 2stage technique for
 the reduction to tridiagonal.
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 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:  if INFO = i, the algorithm failed to converge; i
                off-diagonal elements of an intermediate tridiagonal
                form did not converge to zero.
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 185 of file dsyev_2stage.f.

185 *
186  IMPLICIT NONE
187 *
188 * -- LAPACK driver routine (version 3.7.0) --
189 * -- LAPACK is a software package provided by Univ. of Tennessee, --
190 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191 * December 2016
192 *
193 * .. Scalar Arguments ..
194  CHARACTER jobz, uplo
195  INTEGER info, lda, lwork, n
196 * ..
197 * .. Array Arguments ..
198  DOUBLE PRECISION a( lda, * ), w( * ), work( * )
199 * ..
200 *
201 * =====================================================================
202 *
203 * .. Parameters ..
204  DOUBLE PRECISION zero, one
205  parameter ( zero = 0.0d0, one = 1.0d0 )
206 * ..
207 * .. Local Scalars ..
208  LOGICAL lower, lquery, wantz
209  INTEGER iinfo, imax, inde, indtau, indwrk, iscale,
210  $ llwork, lwmin, lhtrd, lwtrd, kd, ib, indhous
211  DOUBLE PRECISION anrm, bignum, eps, rmax, rmin, safmin, sigma,
212  $ smlnum
213 * ..
214 * .. External Functions ..
215  LOGICAL lsame
216  INTEGER ilaenv
217  DOUBLE PRECISION dlamch, dlansy
218  EXTERNAL lsame, ilaenv, dlamch, dlansy
219 * ..
220 * .. External Subroutines ..
221  EXTERNAL dlascl, dorgtr, dscal, dsteqr, dsterf,
223 * ..
224 * .. Intrinsic Functions ..
225  INTRINSIC max, sqrt
226 * ..
227 * .. Executable Statements ..
228 *
229 * Test the input parameters.
230 *
231  wantz = lsame( jobz, 'V' )
232  lower = lsame( uplo, 'L' )
233  lquery = ( lwork.EQ.-1 )
234 *
235  info = 0
236  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
237  info = -1
238  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
239  info = -2
240  ELSE IF( n.LT.0 ) THEN
241  info = -3
242  ELSE IF( lda.LT.max( 1, n ) ) THEN
243  info = -5
244  END IF
245 *
246  IF( info.EQ.0 ) THEN
247  kd = ilaenv( 17, 'DSYTRD_2STAGE', jobz, n, -1, -1, -1 )
248  ib = ilaenv( 18, 'DSYTRD_2STAGE', jobz, n, kd, -1, -1 )
249  lhtrd = ilaenv( 19, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
250  lwtrd = ilaenv( 20, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
251  lwmin = 2*n + lhtrd + lwtrd
252  work( 1 ) = lwmin
253 *
254  IF( lwork.LT.lwmin .AND. .NOT.lquery )
255  $ info = -8
256  END IF
257 *
258  IF( info.NE.0 ) THEN
259  CALL xerbla( 'DSYEV_2STAGE ', -info )
260  RETURN
261  ELSE IF( lquery ) THEN
262  RETURN
263  END IF
264 *
265 * Quick return if possible
266 *
267  IF( n.EQ.0 ) THEN
268  RETURN
269  END IF
270 *
271  IF( n.EQ.1 ) THEN
272  w( 1 ) = a( 1, 1 )
273  work( 1 ) = 2
274  IF( wantz )
275  $ a( 1, 1 ) = one
276  RETURN
277  END IF
278 *
279 * Get machine constants.
280 *
281  safmin = dlamch( 'Safe minimum' )
282  eps = dlamch( 'Precision' )
283  smlnum = safmin / eps
284  bignum = one / smlnum
285  rmin = sqrt( smlnum )
286  rmax = sqrt( bignum )
287 *
288 * Scale matrix to allowable range, if necessary.
289 *
290  anrm = dlansy( 'M', uplo, n, a, lda, work )
291  iscale = 0
292  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
293  iscale = 1
294  sigma = rmin / anrm
295  ELSE IF( anrm.GT.rmax ) THEN
296  iscale = 1
297  sigma = rmax / anrm
298  END IF
299  IF( iscale.EQ.1 )
300  $ CALL dlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
301 *
302 * Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
303 *
304  inde = 1
305  indtau = inde + n
306  indhous = indtau + n
307  indwrk = indhous + lhtrd
308  llwork = lwork - indwrk + 1
309 *
310  CALL dsytrd_2stage( jobz, uplo, n, a, lda, w, work( inde ),
311  $ work( indtau ), work( indhous ), lhtrd,
312  $ work( indwrk ), llwork, iinfo )
313 *
314 * For eigenvalues only, call DSTERF. For eigenvectors, first call
315 * DORGTR to generate the orthogonal matrix, then call DSTEQR.
316 *
317  IF( .NOT.wantz ) THEN
318  CALL dsterf( n, w, work( inde ), info )
319  ELSE
320 * Not available in this release, and agrument checking should not
321 * let it getting here
322  RETURN
323  CALL dorgtr( uplo, n, a, lda, work( indtau ), work( indwrk ),
324  $ llwork, iinfo )
325  CALL dsteqr( jobz, n, w, work( inde ), a, lda, work( indtau ),
326  $ info )
327  END IF
328 *
329 * If matrix was scaled, then rescale eigenvalues appropriately.
330 *
331  IF( iscale.EQ.1 ) THEN
332  IF( info.EQ.0 ) THEN
333  imax = n
334  ELSE
335  imax = info - 1
336  END IF
337  CALL dscal( imax, one / sigma, w, 1 )
338  END IF
339 *
340 * Set WORK(1) to optimal workspace size.
341 *
342  work( 1 ) = lwmin
343 *
344  RETURN
345 *
346 * End of DSYEV_2STAGE
347 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:133
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, or the element of largest absolute value of a real symmetric matrix.
Definition: dlansy.f:124
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
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine dorgtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
DORGTR
Definition: dorgtr.f:125
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: