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

◆ dsyev_2stage()

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 (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:  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.
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 179 of file dsyev_2stage.f.

181*
182 IMPLICIT NONE
183*
184* -- LAPACK driver routine --
185* -- LAPACK is a software package provided by Univ. of Tennessee, --
186* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
187*
188* .. Scalar Arguments ..
189 CHARACTER JOBZ, UPLO
190 INTEGER INFO, LDA, LWORK, N
191* ..
192* .. Array Arguments ..
193 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
194* ..
195*
196* =====================================================================
197*
198* .. Parameters ..
199 DOUBLE PRECISION ZERO, ONE
200 parameter( zero = 0.0d0, one = 1.0d0 )
201* ..
202* .. Local Scalars ..
203 LOGICAL LOWER, LQUERY, WANTZ
204 INTEGER IINFO, IMAX, INDE, INDTAU, INDWRK, ISCALE,
205 $ LLWORK, LWMIN, LHTRD, LWTRD, KD, IB, INDHOUS
206 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
207 $ SMLNUM
208* ..
209* .. External Functions ..
210 LOGICAL LSAME
211 INTEGER ILAENV2STAGE
212 DOUBLE PRECISION DLAMCH, DLANSY
213 EXTERNAL lsame, dlamch, dlansy, ilaenv2stage
214* ..
215* .. External Subroutines ..
216 EXTERNAL dlascl, dorgtr, dscal, dsteqr,
217 $ dsterf,
219* ..
220* .. Intrinsic Functions ..
221 INTRINSIC max, sqrt
222* ..
223* .. Executable Statements ..
224*
225* Test the input parameters.
226*
227 wantz = lsame( jobz, 'V' )
228 lower = lsame( uplo, 'L' )
229 lquery = ( lwork.EQ.-1 )
230*
231 info = 0
232 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
233 info = -1
234 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
235 info = -2
236 ELSE IF( n.LT.0 ) THEN
237 info = -3
238 ELSE IF( lda.LT.max( 1, n ) ) THEN
239 info = -5
240 END IF
241*
242 IF( info.EQ.0 ) THEN
243 kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz, n, -1, -1,
244 $ -1 )
245 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1,
246 $ -1 )
247 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib,
248 $ -1 )
249 lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz, n, kd, ib,
250 $ -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 argument checking should not
321* let it getting here
322 RETURN
323 CALL dorgtr( uplo, n, a, lda, work( indtau ),
324 $ work( indwrk ),
325 $ llwork, iinfo )
326 CALL dsteqr( jobz, n, w, work( inde ), a, lda,
327 $ work( indtau ),
328 $ info )
329 END IF
330*
331* If matrix was scaled, then rescale eigenvalues appropriately.
332*
333 IF( iscale.EQ.1 ) THEN
334 IF( info.EQ.0 ) THEN
335 imax = n
336 ELSE
337 imax = info - 1
338 END IF
339 CALL dscal( imax, one / sigma, w, 1 )
340 END IF
341*
342* Set WORK(1) to optimal workspace size.
343*
344 work( 1 ) = lwmin
345*
346 RETURN
347*
348* End of DSYEV_2STAGE
349*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
DSYTRD_2STAGE
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
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:120
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:142
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:129
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84
subroutine dorgtr(uplo, n, a, lda, tau, work, lwork, info)
DORGTR
Definition dorgtr.f:121
Here is the call graph for this function:
Here is the caller graph for this function: