LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssyev_2stage.f
Go to the documentation of this file.
1*> \brief <b> SSYEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices</b>
2*
3* @generated from dsyev_2stage.f, fortran d -> s, Sat Nov 5 23:55:51 2016
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> \htmlonly
11*> Download SSYEV_2STAGE + dependencies
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssyev_2stage.f">
13*> [TGZ]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssyev_2stage.f">
15*> [ZIP]</a>
16*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssyev_2stage.f">
17*> [TXT]</a>
18*> \endhtmlonly
19*
20* Definition:
21* ===========
22*
23* SUBROUTINE SSYEV_2STAGE( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK,
24* INFO )
25*
26* IMPLICIT NONE
27*
28* .. Scalar Arguments ..
29* CHARACTER JOBZ, UPLO
30* INTEGER INFO, LDA, LWORK, N
31* ..
32* .. Array Arguments ..
33* REAL A( LDA, * ), W( * ), WORK( * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> SSYEV_2STAGE computes all eigenvalues and, optionally, eigenvectors of a
43*> real symmetric matrix A using the 2stage technique for
44*> the reduction to tridiagonal.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] JOBZ
51*> \verbatim
52*> JOBZ is CHARACTER*1
53*> = 'N': Compute eigenvalues only;
54*> = 'V': Compute eigenvalues and eigenvectors.
55*> Not available in this release.
56*> \endverbatim
57*>
58*> \param[in] UPLO
59*> \verbatim
60*> UPLO is CHARACTER*1
61*> = 'U': Upper triangle of A is stored;
62*> = 'L': Lower triangle of A is stored.
63*> \endverbatim
64*>
65*> \param[in] N
66*> \verbatim
67*> N is INTEGER
68*> The order of the matrix A. N >= 0.
69*> \endverbatim
70*>
71*> \param[in,out] A
72*> \verbatim
73*> A is REAL array, dimension (LDA, N)
74*> On entry, the symmetric matrix A. If UPLO = 'U', the
75*> leading N-by-N upper triangular part of A contains the
76*> upper triangular part of the matrix A. If UPLO = 'L',
77*> the leading N-by-N lower triangular part of A contains
78*> the lower triangular part of the matrix A.
79*> On exit, if JOBZ = 'V', then if INFO = 0, A contains the
80*> orthonormal eigenvectors of the matrix A.
81*> If JOBZ = 'N', then on exit the lower triangle (if UPLO='L')
82*> or the upper triangle (if UPLO='U') of A, including the
83*> diagonal, is destroyed.
84*> \endverbatim
85*>
86*> \param[in] LDA
87*> \verbatim
88*> LDA is INTEGER
89*> The leading dimension of the array A. LDA >= max(1,N).
90*> \endverbatim
91*>
92*> \param[out] W
93*> \verbatim
94*> W is REAL array, dimension (N)
95*> If INFO = 0, the eigenvalues in ascending order.
96*> \endverbatim
97*>
98*> \param[out] WORK
99*> \verbatim
100*> WORK is REAL array, dimension LWORK
101*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
102*> \endverbatim
103*>
104*> \param[in] LWORK
105*> \verbatim
106*> LWORK is INTEGER
107*> The length of the array WORK. LWORK >= 1, when N <= 1;
108*> otherwise
109*> If JOBZ = 'N' and N > 1, LWORK must be queried.
110*> LWORK = MAX(1, dimension) where
111*> dimension = max(stage1,stage2) + (KD+1)*N + 2*N
112*> = N*KD + N*max(KD+1,FACTOPTNB)
113*> + max(2*KD*KD, KD*NTHREADS)
114*> + (KD+1)*N + 2*N
115*> where KD is the blocking size of the reduction,
116*> FACTOPTNB is the blocking used by the QR or LQ
117*> algorithm, usually FACTOPTNB=128 is a good choice
118*> NTHREADS is the number of threads used when
119*> openMP compilation is enabled, otherwise =1.
120*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
121*>
122*> If LWORK = -1, then a workspace query is assumed; the routine
123*> only calculates the optimal size of the WORK array, returns
124*> this value as the first entry of the WORK array, and no error
125*> message related to LWORK is issued by XERBLA.
126*> \endverbatim
127*>
128*> \param[out] INFO
129*> \verbatim
130*> INFO is INTEGER
131*> = 0: successful exit
132*> < 0: if INFO = -i, the i-th argument had an illegal value
133*> > 0: if INFO = i, the algorithm failed to converge; i
134*> off-diagonal elements of an intermediate tridiagonal
135*> form did not converge to zero.
136*> \endverbatim
137*
138* Authors:
139* ========
140*
141*> \author Univ. of Tennessee
142*> \author Univ. of California Berkeley
143*> \author Univ. of Colorado Denver
144*> \author NAG Ltd.
145*
146*> \ingroup heev_2stage
147*
148*> \par Further Details:
149* =====================
150*>
151*> \verbatim
152*>
153*> All details about the 2stage techniques are available in:
154*>
155*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
156*> Parallel reduction to condensed forms for symmetric eigenvalue problems
157*> using aggregated fine-grained and memory-aware kernels. In Proceedings
158*> of 2011 International Conference for High Performance Computing,
159*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
160*> Article 8 , 11 pages.
161*> http://doi.acm.org/10.1145/2063384.2063394
162*>
163*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
164*> An improved parallel singular value algorithm and its implementation
165*> for multicore hardware, In Proceedings of 2013 International Conference
166*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
167*> Denver, Colorado, USA, 2013.
168*> Article 90, 12 pages.
169*> http://doi.acm.org/10.1145/2503210.2503292
170*>
171*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
172*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
173*> calculations based on fine-grained memory aware tasks.
174*> International Journal of High Performance Computing Applications.
175*> Volume 28 Issue 2, Pages 196-209, May 2014.
176*> http://hpc.sagepub.com/content/28/2/196
177*>
178*> \endverbatim
179*
180* =====================================================================
181 SUBROUTINE ssyev_2stage( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK,
182 $ INFO )
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,
216 $ sroundup_lwork
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*
346 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssyev_2stage(jobz, uplo, n, a, lda, w, work, lwork, info)
SSYEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matr...
subroutine ssytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
SSYTRD_2STAGE
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
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