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