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