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

◆ ssyev_2stage()

subroutine ssyev_2stage ( character jobz,
character uplo,
integer n,
real, dimension( lda, * ) a,
integer lda,
real, dimension( * ) w,
real, dimension( * ) work,
integer lwork,
integer info )

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

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

Purpose:
!>
!> SSYEV_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 REAL 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 REAL array, dimension (N)
!>          If INFO = 0, the eigenvalues in ascending order.
!> 
[out]WORK
!>          WORK is REAL 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.
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 ssyev_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 REAL A( LDA, * ), W( * ), WORK( * )
194* ..
195*
196* =====================================================================
197*
198* .. Parameters ..
199 REAL ZERO, ONE
200 parameter( zero = 0.0e0, one = 1.0e0 )
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 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
207 $ SMLNUM
208* ..
209* .. External Functions ..
210 LOGICAL LSAME
211 INTEGER ILAENV2STAGE
212 REAL SLAMCH, SLANSY, SROUNDUP_LWORK
213 EXTERNAL lsame, slamch, slansy, ilaenv2stage,
215* ..
216* .. External Subroutines ..
217 EXTERNAL slascl, sorgtr, sscal, ssteqr,
218 $ ssterf,
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC max, sqrt
223* ..
224* .. Executable Statements ..
225*
226* Test the input parameters.
227*
228 wantz = lsame( jobz, 'V' )
229 lower = lsame( uplo, 'L' )
230 lquery = ( lwork.EQ.-1 )
231*
232 info = 0
233 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
234 info = -1
235 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
236 info = -2
237 ELSE IF( n.LT.0 ) THEN
238 info = -3
239 ELSE IF( lda.LT.max( 1, n ) ) THEN
240 info = -5
241 END IF
242*
243 IF( info.EQ.0 ) THEN
244 kd = ilaenv2stage( 1, 'SSYTRD_2STAGE', jobz, n, -1, -1,
245 $ -1 )
246 ib = ilaenv2stage( 2, 'SSYTRD_2STAGE', jobz, n, kd, -1,
247 $ -1 )
248 lhtrd = ilaenv2stage( 3, 'SSYTRD_2STAGE', jobz, n, kd, ib,
249 $ -1 )
250 lwtrd = ilaenv2stage( 4, 'SSYTRD_2STAGE', jobz, n, kd, ib,
251 $ -1 )
252 lwmin = 2*n + lhtrd + lwtrd
253 work( 1 ) = real( lwmin )
254*
255 IF( lwork.LT.lwmin .AND. .NOT.lquery )
256 $ info = -8
257 END IF
258*
259 IF( info.NE.0 ) THEN
260 CALL xerbla( 'SSYEV_2STAGE ', -info )
261 RETURN
262 ELSE IF( lquery ) THEN
263 RETURN
264 END IF
265*
266* Quick return if possible
267*
268 IF( n.EQ.0 ) THEN
269 RETURN
270 END IF
271*
272 IF( n.EQ.1 ) THEN
273 w( 1 ) = a( 1, 1 )
274 work( 1 ) = 2
275 IF( wantz )
276 $ a( 1, 1 ) = one
277 RETURN
278 END IF
279*
280* Get machine constants.
281*
282 safmin = slamch( 'Safe minimum' )
283 eps = slamch( 'Precision' )
284 smlnum = safmin / eps
285 bignum = one / smlnum
286 rmin = sqrt( smlnum )
287 rmax = sqrt( bignum )
288*
289* Scale matrix to allowable range, if necessary.
290*
291 anrm = slansy( 'M', uplo, n, a, lda, work )
292 iscale = 0
293 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
294 iscale = 1
295 sigma = rmin / anrm
296 ELSE IF( anrm.GT.rmax ) THEN
297 iscale = 1
298 sigma = rmax / anrm
299 END IF
300 IF( iscale.EQ.1 )
301 $ CALL slascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
302*
303* Call SSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
304*
305 inde = 1
306 indtau = inde + n
307 indhous = indtau + n
308 indwrk = indhous + lhtrd
309 llwork = lwork - indwrk + 1
310*
311 CALL ssytrd_2stage( jobz, uplo, n, a, lda, w, work( inde ),
312 $ work( indtau ), work( indhous ), lhtrd,
313 $ work( indwrk ), llwork, iinfo )
314*
315* For eigenvalues only, call SSTERF. For eigenvectors, first call
316* SORGTR to generate the orthogonal matrix, then call SSTEQR.
317*
318 IF( .NOT.wantz ) THEN
319 CALL ssterf( n, w, work( inde ), info )
320 ELSE
321* Not available in this release, and argument checking should not
322* let it getting here
323 RETURN
324 CALL sorgtr( uplo, n, a, lda, work( indtau ),
325 $ work( indwrk ),
326 $ llwork, iinfo )
327 CALL ssteqr( jobz, n, w, work( inde ), a, lda,
328 $ work( indtau ),
329 $ info )
330 END IF
331*
332* If matrix was scaled, then rescale eigenvalues appropriately.
333*
334 IF( iscale.EQ.1 ) THEN
335 IF( info.EQ.0 ) THEN
336 imax = n
337 ELSE
338 imax = info - 1
339 END IF
340 CALL sscal( imax, one / sigma, w, 1 )
341 END IF
342*
343* Set WORK(1) to optimal workspace size.
344*
345 work( 1 ) = sroundup_lwork(lwmin)
346*
347 RETURN
348*
349* End of SSYEV_2STAGE
350*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
SSYTRD_2STAGE
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:120
subroutine slascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
SLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition slascl.f:142
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine ssteqr(compz, n, d, e, z, ldz, work, info)
SSTEQR
Definition ssteqr.f:129
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84
subroutine sorgtr(uplo, n, a, lda, tau, work, lwork, info)
SORGTR
Definition sorgtr.f:121
Here is the call graph for this function:
Here is the caller graph for this function: