LAPACK  3.10.0
LAPACK: Linear Algebra PACKage
dsygv_2stage.f
Go to the documentation of this file.
1 *> \brief \b DSYGV_2STAGE
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 DSYGV_2STAGE + dependencies
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsygv_2stage.f">
13 *> [TGZ]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsygv_2stage.f">
15 *> [ZIP]</a>
16 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsygv_2stage.f">
17 *> [TXT]</a>
18 *> \endhtmlonly
19 *
20 * Definition:
21 * ===========
22 *
23 * SUBROUTINE DSYGV_2STAGE( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W,
24 * WORK, LWORK, INFO )
25 *
26 * IMPLICIT NONE
27 *
28 * .. Scalar Arguments ..
29 * CHARACTER JOBZ, UPLO
30 * INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
31 * ..
32 * .. Array Arguments ..
33 * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
34 * ..
35 *
36 *
37 *> \par Purpose:
38 * =============
39 *>
40 *> \verbatim
41 *>
42 *> DSYGV_2STAGE computes all the eigenvalues, and optionally, the eigenvectors
43 *> of a real generalized symmetric-definite eigenproblem, of the form
44 *> A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
45 *> Here A and B are assumed to be symmetric and B is also
46 *> positive definite.
47 *> This routine use the 2stage technique for the reduction to tridiagonal
48 *> which showed higher performance on recent architecture and for large
49 *> sizes N>2000.
50 *> \endverbatim
51 *
52 * Arguments:
53 * ==========
54 *
55 *> \param[in] ITYPE
56 *> \verbatim
57 *> ITYPE is INTEGER
58 *> Specifies the problem type to be solved:
59 *> = 1: A*x = (lambda)*B*x
60 *> = 2: A*B*x = (lambda)*x
61 *> = 3: B*A*x = (lambda)*x
62 *> \endverbatim
63 *>
64 *> \param[in] JOBZ
65 *> \verbatim
66 *> JOBZ is CHARACTER*1
67 *> = 'N': Compute eigenvalues only;
68 *> = 'V': Compute eigenvalues and eigenvectors.
69 *> Not available in this release.
70 *> \endverbatim
71 *>
72 *> \param[in] UPLO
73 *> \verbatim
74 *> UPLO is CHARACTER*1
75 *> = 'U': Upper triangles of A and B are stored;
76 *> = 'L': Lower triangles of A and B are stored.
77 *> \endverbatim
78 *>
79 *> \param[in] N
80 *> \verbatim
81 *> N is INTEGER
82 *> The order of the matrices A and B. N >= 0.
83 *> \endverbatim
84 *>
85 *> \param[in,out] A
86 *> \verbatim
87 *> A is DOUBLE PRECISION array, dimension (LDA, N)
88 *> On entry, the symmetric matrix A. If UPLO = 'U', the
89 *> leading N-by-N upper triangular part of A contains the
90 *> upper triangular part of the matrix A. If UPLO = 'L',
91 *> the leading N-by-N lower triangular part of A contains
92 *> the lower triangular part of the matrix A.
93 *>
94 *> On exit, if JOBZ = 'V', then if INFO = 0, A contains the
95 *> matrix Z of eigenvectors. The eigenvectors are normalized
96 *> as follows:
97 *> if ITYPE = 1 or 2, Z**T*B*Z = I;
98 *> if ITYPE = 3, Z**T*inv(B)*Z = I.
99 *> If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
100 *> or the lower triangle (if UPLO='L') of A, including the
101 *> diagonal, is destroyed.
102 *> \endverbatim
103 *>
104 *> \param[in] LDA
105 *> \verbatim
106 *> LDA is INTEGER
107 *> The leading dimension of the array A. LDA >= max(1,N).
108 *> \endverbatim
109 *>
110 *> \param[in,out] B
111 *> \verbatim
112 *> B is DOUBLE PRECISION array, dimension (LDB, N)
113 *> On entry, the symmetric positive definite matrix B.
114 *> If UPLO = 'U', the leading N-by-N upper triangular part of B
115 *> contains the upper triangular part of the matrix B.
116 *> If UPLO = 'L', the leading N-by-N lower triangular part of B
117 *> contains the lower triangular part of the matrix B.
118 *>
119 *> On exit, if INFO <= N, the part of B containing the matrix is
120 *> overwritten by the triangular factor U or L from the Cholesky
121 *> factorization B = U**T*U or B = L*L**T.
122 *> \endverbatim
123 *>
124 *> \param[in] LDB
125 *> \verbatim
126 *> LDB is INTEGER
127 *> The leading dimension of the array B. LDB >= max(1,N).
128 *> \endverbatim
129 *>
130 *> \param[out] W
131 *> \verbatim
132 *> W is DOUBLE PRECISION array, dimension (N)
133 *> If INFO = 0, the eigenvalues in ascending order.
134 *> \endverbatim
135 *>
136 *> \param[out] WORK
137 *> \verbatim
138 *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
139 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
140 *> \endverbatim
141 *>
142 *> \param[in] LWORK
143 *> \verbatim
144 *> LWORK is INTEGER
145 *> The length of the array WORK. LWORK >= 1, when N <= 1;
146 *> otherwise
147 *> If JOBZ = 'N' and N > 1, LWORK must be queried.
148 *> LWORK = MAX(1, dimension) where
149 *> dimension = max(stage1,stage2) + (KD+1)*N + 2*N
150 *> = N*KD + N*max(KD+1,FACTOPTNB)
151 *> + max(2*KD*KD, KD*NTHREADS)
152 *> + (KD+1)*N + 2*N
153 *> where KD is the blocking size of the reduction,
154 *> FACTOPTNB is the blocking used by the QR or LQ
155 *> algorithm, usually FACTOPTNB=128 is a good choice
156 *> NTHREADS is the number of threads used when
157 *> openMP compilation is enabled, otherwise =1.
158 *> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
159 *>
160 *> If LWORK = -1, then a workspace query is assumed; the routine
161 *> only calculates the optimal size of the WORK array, returns
162 *> this value as the first entry of the WORK array, and no error
163 *> message related to LWORK is issued by XERBLA.
164 *> \endverbatim
165 *>
166 *> \param[out] INFO
167 *> \verbatim
168 *> INFO is INTEGER
169 *> = 0: successful exit
170 *> < 0: if INFO = -i, the i-th argument had an illegal value
171 *> > 0: DPOTRF or DSYEV returned an error code:
172 *> <= N: if INFO = i, DSYEV failed to converge;
173 *> i off-diagonal elements of an intermediate
174 *> tridiagonal form did not converge to zero;
175 *> > N: if INFO = N + i, for 1 <= i <= N, then the leading
176 *> minor of order i of B is not positive definite.
177 *> The factorization of B could not be completed and
178 *> no eigenvalues or eigenvectors were computed.
179 *> \endverbatim
180 *
181 * Authors:
182 * ========
183 *
184 *> \author Univ. of Tennessee
185 *> \author Univ. of California Berkeley
186 *> \author Univ. of Colorado Denver
187 *> \author NAG Ltd.
188 *
189 *> \ingroup doubleSYeigen
190 *
191 *> \par Further Details:
192 * =====================
193 *>
194 *> \verbatim
195 *>
196 *> All details about the 2stage techniques are available in:
197 *>
198 *> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
199 *> Parallel reduction to condensed forms for symmetric eigenvalue problems
200 *> using aggregated fine-grained and memory-aware kernels. In Proceedings
201 *> of 2011 International Conference for High Performance Computing,
202 *> Networking, Storage and Analysis (SC '11), New York, NY, USA,
203 *> Article 8 , 11 pages.
204 *> http://doi.acm.org/10.1145/2063384.2063394
205 *>
206 *> A. Haidar, J. Kurzak, P. Luszczek, 2013.
207 *> An improved parallel singular value algorithm and its implementation
208 *> for multicore hardware, In Proceedings of 2013 International Conference
209 *> for High Performance Computing, Networking, Storage and Analysis (SC '13).
210 *> Denver, Colorado, USA, 2013.
211 *> Article 90, 12 pages.
212 *> http://doi.acm.org/10.1145/2503210.2503292
213 *>
214 *> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
215 *> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
216 *> calculations based on fine-grained memory aware tasks.
217 *> International Journal of High Performance Computing Applications.
218 *> Volume 28 Issue 2, Pages 196-209, May 2014.
219 *> http://hpc.sagepub.com/content/28/2/196
220 *>
221 *> \endverbatim
222 *
223 * =====================================================================
224  SUBROUTINE dsygv_2stage( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W,
225  $ WORK, LWORK, INFO )
226 *
227  IMPLICIT NONE
228 *
229 * -- LAPACK driver routine --
230 * -- LAPACK is a software package provided by Univ. of Tennessee, --
231 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
232 *
233 * .. Scalar Arguments ..
234  CHARACTER JOBZ, UPLO
235  INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
236 * ..
237 * .. Array Arguments ..
238  DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
239 * ..
240 *
241 * =====================================================================
242 *
243 * .. Parameters ..
244  DOUBLE PRECISION ONE
245  parameter( one = 1.0d+0 )
246 * ..
247 * .. Local Scalars ..
248  LOGICAL LQUERY, UPPER, WANTZ
249  CHARACTER TRANS
250  INTEGER NEIG, LWMIN, LHTRD, LWTRD, KD, IB
251 * ..
252 * .. External Functions ..
253  LOGICAL LSAME
254  INTEGER ILAENV2STAGE
255  EXTERNAL lsame, ilaenv2stage
256 * ..
257 * .. External Subroutines ..
258  EXTERNAL dpotrf, dsygst, dtrmm, dtrsm, xerbla,
259  $ dsyev_2stage
260 * ..
261 * .. Intrinsic Functions ..
262  INTRINSIC max
263 * ..
264 * .. Executable Statements ..
265 *
266 * Test the input parameters.
267 *
268  wantz = lsame( jobz, 'V' )
269  upper = lsame( uplo, 'U' )
270  lquery = ( lwork.EQ.-1 )
271 *
272  info = 0
273  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
274  info = -1
275  ELSE IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
276  info = -2
277  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
278  info = -3
279  ELSE IF( n.LT.0 ) THEN
280  info = -4
281  ELSE IF( lda.LT.max( 1, n ) ) THEN
282  info = -6
283  ELSE IF( ldb.LT.max( 1, n ) ) THEN
284  info = -8
285  END IF
286 *
287  IF( info.EQ.0 ) THEN
288  kd = ilaenv2stage( 1, 'DSYTRD_2STAGE', jobz, n, -1, -1, -1 )
289  ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1, -1 )
290  lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
291  lwtrd = ilaenv2stage( 4, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
292  lwmin = 2*n + lhtrd + lwtrd
293  work( 1 ) = lwmin
294 *
295  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
296  info = -11
297  END IF
298  END IF
299 *
300  IF( info.NE.0 ) THEN
301  CALL xerbla( 'DSYGV_2STAGE ', -info )
302  RETURN
303  ELSE IF( lquery ) THEN
304  RETURN
305  END IF
306 *
307 * Quick return if possible
308 *
309  IF( n.EQ.0 )
310  $ RETURN
311 *
312 * Form a Cholesky factorization of B.
313 *
314  CALL dpotrf( uplo, n, b, ldb, info )
315  IF( info.NE.0 ) THEN
316  info = n + info
317  RETURN
318  END IF
319 *
320 * Transform problem to standard eigenvalue problem and solve.
321 *
322  CALL dsygst( itype, uplo, n, a, lda, b, ldb, info )
323  CALL dsyev_2stage( jobz, uplo, n, a, lda, w, work, lwork, info )
324 *
325  IF( wantz ) THEN
326 *
327 * Backtransform eigenvectors to the original problem.
328 *
329  neig = n
330  IF( info.GT.0 )
331  $ neig = info - 1
332  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
333 *
334 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
335 * backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
336 *
337  IF( upper ) THEN
338  trans = 'N'
339  ELSE
340  trans = 'T'
341  END IF
342 *
343  CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
344  $ b, ldb, a, lda )
345 *
346  ELSE IF( itype.EQ.3 ) THEN
347 *
348 * For B*A*x=(lambda)*x;
349 * backtransform eigenvectors: x = L*y or U**T*y
350 *
351  IF( upper ) THEN
352  trans = 'T'
353  ELSE
354  trans = 'N'
355  END IF
356 *
357  CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
358  $ b, ldb, a, lda )
359  END IF
360  END IF
361 *
362  work( 1 ) = lwmin
363  RETURN
364 *
365 * End of DSYGV_2STAGE
366 *
367  END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRSM
Definition: dtrsm.f:181
subroutine dtrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRMM
Definition: dtrmm.f:177
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:107
subroutine dsygst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
DSYGST
Definition: dsygst.f:127
subroutine dsygv_2stage(ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK, LWORK, INFO)
DSYGV_2STAGE
Definition: dsygv_2stage.f:226
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