LAPACK  3.10.0
LAPACK: Linear Algebra PACKage
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 *> \htmlonly
11 *> Download DSYEV_2STAGE + dependencies
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyev_2stage.f">
13 *> [TGZ]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyev_2stage.f">
15 *> [ZIP]</a>
16 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyev_2stage.f">
17 *> [TXT]</a>
18 *> \endhtmlonly
19 *
20 * Definition:
21 * ===========
22 *
23 * SUBROUTINE DSYEV_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 * DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
34 * ..
35 *
36 *
37 *> \par Purpose:
38 * =============
39 *>
40 *> \verbatim
41 *>
42 *> DSYEV_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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
95 *> If INFO = 0, the eigenvalues in ascending order.
96 *> \endverbatim
97 *>
98 *> \param[out] WORK
99 *> \verbatim
100 *> WORK is DOUBLE PRECISION 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 doubleSYeigen
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 dsyev_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  DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * )
196 * ..
197 *
198 * =====================================================================
199 *
200 * .. Parameters ..
201  DOUBLE PRECISION ZERO, ONE
202  parameter( zero = 0.0d0, one = 1.0d0 )
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  DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
209  $ smlnum
210 * ..
211 * .. External Functions ..
212  LOGICAL LSAME
213  INTEGER ILAENV2STAGE
214  DOUBLE PRECISION DLAMCH, DLANSY
215  EXTERNAL lsame, dlamch, dlansy, ilaenv2stage
216 * ..
217 * .. External Subroutines ..
218  EXTERNAL dlascl, dorgtr, dscal, dsteqr, dsterf,
220 * ..
221 * .. Intrinsic Functions ..
222  INTRINSIC max, sqrt
223 * ..
224 * .. Executable Statements ..
225 *
226 * Test the input parameters.
227 *
228  wantz = lsame( jobz, 'V' )
229  lower = lsame( uplo, 'L' )
230  lquery = ( lwork.EQ.-1 )
231 *
232  info = 0
233  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
234  info = -1
235  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
236  info = -2
237  ELSE IF( n.LT.0 ) THEN
238  info = -3
239  ELSE IF( lda.LT.max( 1, n ) ) THEN
240  info = -5
241  END IF
242 *
243  IF( info.EQ.0 ) THEN
244  kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz, n, -1, -1, -1 )
245  ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1, -1 )
246  lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
247  lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
248  lwmin = 2*n + lhtrd + lwtrd
249  work( 1 ) = lwmin
250 *
251  IF( lwork.LT.lwmin .AND. .NOT.lquery )
252  $ info = -8
253  END IF
254 *
255  IF( info.NE.0 ) THEN
256  CALL xerbla( 'DSYEV_2STAGE ', -info )
257  RETURN
258  ELSE IF( lquery ) THEN
259  RETURN
260  END IF
261 *
262 * Quick return if possible
263 *
264  IF( n.EQ.0 ) THEN
265  RETURN
266  END IF
267 *
268  IF( n.EQ.1 ) THEN
269  w( 1 ) = a( 1, 1 )
270  work( 1 ) = 2
271  IF( wantz )
272  $ a( 1, 1 ) = one
273  RETURN
274  END IF
275 *
276 * Get machine constants.
277 *
278  safmin = dlamch( 'Safe minimum' )
279  eps = dlamch( 'Precision' )
280  smlnum = safmin / eps
281  bignum = one / smlnum
282  rmin = sqrt( smlnum )
283  rmax = sqrt( bignum )
284 *
285 * Scale matrix to allowable range, if necessary.
286 *
287  anrm = dlansy( 'M', uplo, n, a, lda, work )
288  iscale = 0
289  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
290  iscale = 1
291  sigma = rmin / anrm
292  ELSE IF( anrm.GT.rmax ) THEN
293  iscale = 1
294  sigma = rmax / anrm
295  END IF
296  IF( iscale.EQ.1 )
297  $ CALL dlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
298 *
299 * Call DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
300 *
301  inde = 1
302  indtau = inde + n
303  indhous = indtau + n
304  indwrk = indhous + lhtrd
305  llwork = lwork - indwrk + 1
306 *
307  CALL dsytrd_2stage( jobz, uplo, n, a, lda, w, work( inde ),
308  $ work( indtau ), work( indhous ), lhtrd,
309  $ work( indwrk ), llwork, iinfo )
310 *
311 * For eigenvalues only, call DSTERF. For eigenvectors, first call
312 * DORGTR to generate the orthogonal matrix, then call DSTEQR.
313 *
314  IF( .NOT.wantz ) THEN
315  CALL dsterf( n, w, work( inde ), info )
316  ELSE
317 * Not available in this release, and argument checking should not
318 * let it getting here
319  RETURN
320  CALL dorgtr( uplo, n, a, lda, work( indtau ), work( indwrk ),
321  $ llwork, iinfo )
322  CALL dsteqr( jobz, n, w, work( inde ), a, lda, work( indtau ),
323  $ info )
324  END IF
325 *
326 * If matrix was scaled, then rescale eigenvalues appropriately.
327 *
328  IF( iscale.EQ.1 ) THEN
329  IF( info.EQ.0 ) THEN
330  imax = n
331  ELSE
332  imax = info - 1
333  END IF
334  CALL dscal( imax, one / sigma, w, 1 )
335  END IF
336 *
337 * Set WORK(1) to optimal workspace size.
338 *
339  work( 1 ) = lwmin
340 *
341  RETURN
342 *
343 * End of DSYEV_2STAGE
344 *
345  END
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:143
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:131
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:86
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:79
subroutine dorgtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
DORGTR
Definition: dorgtr.f:123
subroutine dsytrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
DSYTRD_2STAGE
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...
Definition: dsyev_2stage.f:183