LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssbev_2stage.f
Go to the documentation of this file.
1*> \brief <b> SSBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* @generated from dsbev_2stage.f, fortran d -> s, Sat Nov 5 23:58:09 2016
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> Download SSBEV_2STAGE + dependencies
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssbev_2stage.f">
12*> [TGZ]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssbev_2stage.f">
14*> [ZIP]</a>
15*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssbev_2stage.f">
16*> [TXT]</a>
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SSBEV_2STAGE( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ,
22* WORK, LWORK, INFO )
23*
24* IMPLICIT NONE
25*
26* .. Scalar Arguments ..
27* CHARACTER JOBZ, UPLO
28* INTEGER INFO, KD, LDAB, LDZ, N, LWORK
29* ..
30* .. Array Arguments ..
31* REAL AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> SSBEV_2STAGE computes all the eigenvalues and, optionally, eigenvectors of
41*> a real symmetric band matrix A using the 2stage technique for
42*> the reduction to tridiagonal.
43*> \endverbatim
44*
45* Arguments:
46* ==========
47*
48*> \param[in] JOBZ
49*> \verbatim
50*> JOBZ is CHARACTER*1
51*> = 'N': Compute eigenvalues only;
52*> = 'V': Compute eigenvalues and eigenvectors.
53*> Not available in this release.
54*> \endverbatim
55*>
56*> \param[in] UPLO
57*> \verbatim
58*> UPLO is CHARACTER*1
59*> = 'U': Upper triangle of A is stored;
60*> = 'L': Lower triangle of A is stored.
61*> \endverbatim
62*>
63*> \param[in] N
64*> \verbatim
65*> N is INTEGER
66*> The order of the matrix A. N >= 0.
67*> \endverbatim
68*>
69*> \param[in] KD
70*> \verbatim
71*> KD is INTEGER
72*> The number of superdiagonals of the matrix A if UPLO = 'U',
73*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
74*> \endverbatim
75*>
76*> \param[in,out] AB
77*> \verbatim
78*> AB is REAL array, dimension (LDAB, N)
79*> On entry, the upper or lower triangle of the symmetric band
80*> matrix A, stored in the first KD+1 rows of the array. The
81*> j-th column of A is stored in the j-th column of the array AB
82*> as follows:
83*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
84*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
85*>
86*> On exit, AB is overwritten by values generated during the
87*> reduction to tridiagonal form. If UPLO = 'U', the first
88*> superdiagonal and the diagonal of the tridiagonal matrix T
89*> are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
90*> the diagonal and first subdiagonal of T are returned in the
91*> first two rows of AB.
92*> \endverbatim
93*>
94*> \param[in] LDAB
95*> \verbatim
96*> LDAB is INTEGER
97*> The leading dimension of the array AB. LDAB >= KD + 1.
98*> \endverbatim
99*>
100*> \param[out] W
101*> \verbatim
102*> W is REAL array, dimension (N)
103*> If INFO = 0, the eigenvalues in ascending order.
104*> \endverbatim
105*>
106*> \param[out] Z
107*> \verbatim
108*> Z is REAL array, dimension (LDZ, N)
109*> If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
110*> eigenvectors of the matrix A, with the i-th column of Z
111*> holding the eigenvector associated with W(i).
112*> If JOBZ = 'N', then Z is not referenced.
113*> \endverbatim
114*>
115*> \param[in] LDZ
116*> \verbatim
117*> LDZ is INTEGER
118*> The leading dimension of the array Z. LDZ >= 1, and if
119*> JOBZ = 'V', LDZ >= max(1,N).
120*> \endverbatim
121*>
122*> \param[out] WORK
123*> \verbatim
124*> WORK is REAL array, dimension LWORK
125*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
126*> \endverbatim
127*>
128*> \param[in] LWORK
129*> \verbatim
130*> LWORK is INTEGER
131*> The length of the array WORK. LWORK >= 1, when N <= 1;
132*> otherwise
133*> If JOBZ = 'N' and N > 1, LWORK must be queried.
134*> LWORK = MAX(1, dimension) where
135*> dimension = (2KD+1)*N + KD*NTHREADS + N
136*> where KD is the size of the band.
137*> NTHREADS is the number of threads used when
138*> openMP compilation is enabled, otherwise =1.
139*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.
140*>
141*> If LWORK = -1, then a workspace query is assumed; the routine
142*> only calculates the optimal size of the WORK array, returns
143*> this value as the first entry of the WORK array, and no error
144*> message related to LWORK is issued by XERBLA.
145*> \endverbatim
146*>
147*> \param[out] INFO
148*> \verbatim
149*> INFO is INTEGER
150*> = 0: successful exit
151*> < 0: if INFO = -i, the i-th argument had an illegal value
152*> > 0: if INFO = i, the algorithm failed to converge; i
153*> off-diagonal elements of an intermediate tridiagonal
154*> form did not converge to zero.
155*> \endverbatim
156*
157* Authors:
158* ========
159*
160*> \author Univ. of Tennessee
161*> \author Univ. of California Berkeley
162*> \author Univ. of Colorado Denver
163*> \author NAG Ltd.
164*
165*> \ingroup hbev_2stage
166*
167*> \par Further Details:
168* =====================
169*>
170*> \verbatim
171*>
172*> All details about the 2stage techniques are available in:
173*>
174*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
175*> Parallel reduction to condensed forms for symmetric eigenvalue problems
176*> using aggregated fine-grained and memory-aware kernels. In Proceedings
177*> of 2011 International Conference for High Performance Computing,
178*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
179*> Article 8 , 11 pages.
180*> http://doi.acm.org/10.1145/2063384.2063394
181*>
182*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
183*> An improved parallel singular value algorithm and its implementation
184*> for multicore hardware, In Proceedings of 2013 International Conference
185*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
186*> Denver, Colorado, USA, 2013.
187*> Article 90, 12 pages.
188*> http://doi.acm.org/10.1145/2503210.2503292
189*>
190*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
191*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
192*> calculations based on fine-grained memory aware tasks.
193*> International Journal of High Performance Computing Applications.
194*> Volume 28 Issue 2, Pages 196-209, May 2014.
195*> http://hpc.sagepub.com/content/28/2/196
196*>
197*> \endverbatim
198*
199* =====================================================================
200 SUBROUTINE ssbev_2stage( JOBZ, UPLO, N, KD, AB, LDAB, W, Z,
201 $ LDZ,
202 $ WORK, LWORK, INFO )
203*
204 IMPLICIT NONE
205*
206* -- LAPACK driver routine --
207* -- LAPACK is a software package provided by Univ. of Tennessee, --
208* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
209*
210* .. Scalar Arguments ..
211 CHARACTER JOBZ, UPLO
212 INTEGER INFO, KD, LDAB, LDZ, N, LWORK
213* ..
214* .. Array Arguments ..
215 REAL AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
216* ..
217*
218* =====================================================================
219*
220* .. Parameters ..
221 REAL ZERO, ONE
222 PARAMETER ( ZERO = 0.0e0, one = 1.0e0 )
223* ..
224* .. Local Scalars ..
225 LOGICAL LOWER, WANTZ, LQUERY
226 INTEGER IINFO, IMAX, INDE, INDWRK, ISCALE,
227 $ llwork, lwmin, lhtrd, lwtrd, ib, indhous
228 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
229 $ SMLNUM
230* ..
231* .. External Functions ..
232 LOGICAL LSAME
233 INTEGER ILAENV2STAGE
234 REAL SLAMCH, SLANSB, SROUNDUP_LWORK
235 EXTERNAL lsame, slamch, slansb, ilaenv2stage,
236 $ sroundup_lwork
237* ..
238* .. External Subroutines ..
239 EXTERNAL slascl, sscal, ssteqr, ssterf,
240 $ xerbla,
242* ..
243* .. Intrinsic Functions ..
244 INTRINSIC sqrt
245* ..
246* .. Executable Statements ..
247*
248* Test the input parameters.
249*
250 wantz = lsame( jobz, 'V' )
251 lower = lsame( uplo, 'L' )
252 lquery = ( lwork.EQ.-1 )
253*
254 info = 0
255 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
256 info = -1
257 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
258 info = -2
259 ELSE IF( n.LT.0 ) THEN
260 info = -3
261 ELSE IF( kd.LT.0 ) THEN
262 info = -4
263 ELSE IF( ldab.LT.kd+1 ) THEN
264 info = -6
265 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
266 info = -9
267 END IF
268*
269 IF( info.EQ.0 ) THEN
270 IF( n.LE.1 ) THEN
271 lwmin = 1
272 work( 1 ) = sroundup_lwork(lwmin)
273 ELSE
274 ib = ilaenv2stage( 2, 'SSYTRD_SB2ST', jobz,
275 $ n, kd, -1, -1 )
276 lhtrd = ilaenv2stage( 3, 'SSYTRD_SB2ST', jobz,
277 $ n, kd, ib, -1 )
278 lwtrd = ilaenv2stage( 4, 'SSYTRD_SB2ST', jobz,
279 $ n, kd, ib, -1 )
280 lwmin = n + lhtrd + lwtrd
281 work( 1 ) = sroundup_lwork(lwmin)
282 ENDIF
283*
284 IF( lwork.LT.lwmin .AND. .NOT.lquery )
285 $ info = -11
286 END IF
287*
288 IF( info.NE.0 ) THEN
289 CALL xerbla( 'SSBEV_2STAGE ', -info )
290 RETURN
291 ELSE IF( lquery ) THEN
292 RETURN
293 END IF
294*
295* Quick return if possible
296*
297 IF( n.EQ.0 )
298 $ RETURN
299*
300 IF( n.EQ.1 ) THEN
301 IF( lower ) THEN
302 w( 1 ) = ab( 1, 1 )
303 ELSE
304 w( 1 ) = ab( kd+1, 1 )
305 END IF
306 IF( wantz )
307 $ z( 1, 1 ) = one
308 RETURN
309 END IF
310*
311* Get machine constants.
312*
313 safmin = slamch( 'Safe minimum' )
314 eps = slamch( 'Precision' )
315 smlnum = safmin / eps
316 bignum = one / smlnum
317 rmin = sqrt( smlnum )
318 rmax = sqrt( bignum )
319*
320* Scale matrix to allowable range, if necessary.
321*
322 anrm = slansb( 'M', uplo, n, kd, ab, ldab, work )
323 iscale = 0
324 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
325 iscale = 1
326 sigma = rmin / anrm
327 ELSE IF( anrm.GT.rmax ) THEN
328 iscale = 1
329 sigma = rmax / anrm
330 END IF
331 IF( iscale.EQ.1 ) THEN
332 IF( lower ) THEN
333 CALL slascl( 'B', kd, kd, one, sigma, n, n, ab, ldab,
334 $ info )
335 ELSE
336 CALL slascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab,
337 $ info )
338 END IF
339 END IF
340*
341* Call SSYTRD_SB2ST to reduce symmetric band matrix to tridiagonal form.
342*
343 inde = 1
344 indhous = inde + n
345 indwrk = indhous + lhtrd
346 llwork = lwork - indwrk + 1
347*
348 CALL ssytrd_sb2st( "N", jobz, uplo, n, kd, ab, ldab, w,
349 $ work( inde ), work( indhous ), lhtrd,
350 $ work( indwrk ), llwork, iinfo )
351*
352* For eigenvalues only, call SSTERF. For eigenvectors, call SSTEQR.
353*
354 IF( .NOT.wantz ) THEN
355 CALL ssterf( n, w, work( inde ), info )
356 ELSE
357 CALL ssteqr( jobz, n, w, work( inde ), z, ldz,
358 $ work( indwrk ),
359 $ info )
360 END IF
361*
362* If matrix was scaled, then rescale eigenvalues appropriately.
363*
364 IF( iscale.EQ.1 ) THEN
365 IF( info.EQ.0 ) THEN
366 imax = n
367 ELSE
368 imax = info - 1
369 END IF
370 CALL sscal( imax, one / sigma, w, 1 )
371 END IF
372*
373* Set WORK(1) to optimal workspace size.
374*
375 work( 1 ) = sroundup_lwork(lwmin)
376*
377 RETURN
378*
379* End of SSBEV_2STAGE
380*
381 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssbev_2stage(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, lwork, info)
SSBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER m...
subroutine ssytrd_sb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
SSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric tridiagonal form T
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
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