LAPACK 3.12.0
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 181 of file ssyev_2stage.f.

183*
184 IMPLICIT NONE
185*
186* -- LAPACK driver routine --
187* -- LAPACK is a software package provided by Univ. of Tennessee, --
188* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189*
190* .. Scalar Arguments ..
191 CHARACTER JOBZ, UPLO
192 INTEGER INFO, LDA, LWORK, N
193* ..
194* .. Array Arguments ..
195 REAL A( LDA, * ), W( * ), WORK( * )
196* ..
197*
198* =====================================================================
199*
200* .. Parameters ..
201 REAL ZERO, ONE
202 parameter( zero = 0.0e0, one = 1.0e0 )
203* ..
204* .. Local Scalars ..
205 LOGICAL LOWER, LQUERY, WANTZ
206 INTEGER IINFO, IMAX, INDE, INDTAU, INDWRK, ISCALE,
207 $ LLWORK, LWMIN, LHTRD, LWTRD, KD, IB, INDHOUS
208 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
209 $ SMLNUM
210* ..
211* .. External Functions ..
212 LOGICAL LSAME
213 INTEGER ILAENV2STAGE
214 REAL SLAMCH, SLANSY, SROUNDUP_LWORK
215 EXTERNAL lsame, slamch, slansy, ilaenv2stage,
217* ..
218* .. External Subroutines ..
219 EXTERNAL slascl, sorgtr, sscal, ssteqr, ssterf,
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC max, sqrt
224* ..
225* .. Executable Statements ..
226*
227* Test the input parameters.
228*
229 wantz = lsame( jobz, 'V' )
230 lower = lsame( uplo, 'L' )
231 lquery = ( lwork.EQ.-1 )
232*
233 info = 0
234 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
235 info = -1
236 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
237 info = -2
238 ELSE IF( n.LT.0 ) THEN
239 info = -3
240 ELSE IF( lda.LT.max( 1, n ) ) THEN
241 info = -5
242 END IF
243*
244 IF( info.EQ.0 ) THEN
245 kd = ilaenv2stage( 1, 'SSYTRD_2STAGE', jobz, n, -1, -1, -1 )
246 ib = ilaenv2stage( 2, 'SSYTRD_2STAGE', jobz, n, kd, -1, -1 )
247 lhtrd = ilaenv2stage( 3, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
248 lwtrd = ilaenv2stage( 4, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
249 lwmin = 2*n + lhtrd + lwtrd
250 work( 1 ) = lwmin
251*
252 IF( lwork.LT.lwmin .AND. .NOT.lquery )
253 $ info = -8
254 END IF
255*
256 IF( info.NE.0 ) THEN
257 CALL xerbla( 'SSYEV_2STAGE ', -info )
258 RETURN
259 ELSE IF( lquery ) THEN
260 RETURN
261 END IF
262*
263* Quick return if possible
264*
265 IF( n.EQ.0 ) THEN
266 RETURN
267 END IF
268*
269 IF( n.EQ.1 ) THEN
270 w( 1 ) = a( 1, 1 )
271 work( 1 ) = 2
272 IF( wantz )
273 $ a( 1, 1 ) = one
274 RETURN
275 END IF
276*
277* Get machine constants.
278*
279 safmin = slamch( 'Safe minimum' )
280 eps = slamch( 'Precision' )
281 smlnum = safmin / eps
282 bignum = one / smlnum
283 rmin = sqrt( smlnum )
284 rmax = sqrt( bignum )
285*
286* Scale matrix to allowable range, if necessary.
287*
288 anrm = slansy( 'M', uplo, n, a, lda, work )
289 iscale = 0
290 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
291 iscale = 1
292 sigma = rmin / anrm
293 ELSE IF( anrm.GT.rmax ) THEN
294 iscale = 1
295 sigma = rmax / anrm
296 END IF
297 IF( iscale.EQ.1 )
298 $ CALL slascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
299*
300* Call SSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
301*
302 inde = 1
303 indtau = inde + n
304 indhous = indtau + n
305 indwrk = indhous + lhtrd
306 llwork = lwork - indwrk + 1
307*
308 CALL ssytrd_2stage( jobz, uplo, n, a, lda, w, work( inde ),
309 $ work( indtau ), work( indhous ), lhtrd,
310 $ work( indwrk ), llwork, iinfo )
311*
312* For eigenvalues only, call SSTERF. For eigenvectors, first call
313* SORGTR to generate the orthogonal matrix, then call SSTEQR.
314*
315 IF( .NOT.wantz ) THEN
316 CALL ssterf( n, w, work( inde ), info )
317 ELSE
318* Not available in this release, and argument checking should not
319* let it getting here
320 RETURN
321 CALL sorgtr( uplo, n, a, lda, work( indtau ), work( indwrk ),
322 $ llwork, iinfo )
323 CALL ssteqr( jobz, n, w, work( inde ), a, lda, work( indtau ),
324 $ info )
325 END IF
326*
327* If matrix was scaled, then rescale eigenvalues appropriately.
328*
329 IF( iscale.EQ.1 ) THEN
330 IF( info.EQ.0 ) THEN
331 imax = n
332 ELSE
333 imax = info - 1
334 END IF
335 CALL sscal( imax, one / sigma, w, 1 )
336 END IF
337*
338* Set WORK(1) to optimal workspace size.
339*
340 work( 1 ) = sroundup_lwork(lwmin)
341*
342 RETURN
343*
344* End of SSYEV_2STAGE
345*
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:122
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:143
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:131
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
subroutine sorgtr(uplo, n, a, lda, tau, work, lwork, info)
SORGTR
Definition sorgtr.f:123
Here is the call graph for this function:
Here is the caller graph for this function: