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