LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssygv.f
Go to the documentation of this file.
1*> \brief \b SSYGV
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SSYGV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssygv.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssygv.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssygv.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SSYGV( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
22* LWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER JOBZ, UPLO
26* INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
27* ..
28* .. Array Arguments ..
29* REAL A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> SSYGV computes all the eigenvalues, and optionally, the eigenvectors
39*> of a real generalized symmetric-definite eigenproblem, of the form
40*> A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
41*> Here A and B are assumed to be symmetric and B is also
42*> positive definite.
43*> \endverbatim
44*
45* Arguments:
46* ==========
47*
48*> \param[in] ITYPE
49*> \verbatim
50*> ITYPE is INTEGER
51*> Specifies the problem type to be solved:
52*> = 1: A*x = (lambda)*B*x
53*> = 2: A*B*x = (lambda)*x
54*> = 3: B*A*x = (lambda)*x
55*> \endverbatim
56*>
57*> \param[in] JOBZ
58*> \verbatim
59*> JOBZ is CHARACTER*1
60*> = 'N': Compute eigenvalues only;
61*> = 'V': Compute eigenvalues and eigenvectors.
62*> \endverbatim
63*>
64*> \param[in] UPLO
65*> \verbatim
66*> UPLO is CHARACTER*1
67*> = 'U': Upper triangles of A and B are stored;
68*> = 'L': Lower triangles of A and B are stored.
69*> \endverbatim
70*>
71*> \param[in] N
72*> \verbatim
73*> N is INTEGER
74*> The order of the matrices A and B. N >= 0.
75*> \endverbatim
76*>
77*> \param[in,out] A
78*> \verbatim
79*> A is REAL array, dimension (LDA, N)
80*> On entry, the symmetric matrix A. If UPLO = 'U', the
81*> leading N-by-N upper triangular part of A contains the
82*> upper triangular part of the matrix A. If UPLO = 'L',
83*> the leading N-by-N lower triangular part of A contains
84*> the lower triangular part of the matrix A.
85*>
86*> On exit, if JOBZ = 'V', then if INFO = 0, A contains the
87*> matrix Z of eigenvectors. The eigenvectors are normalized
88*> as follows:
89*> if ITYPE = 1 or 2, Z**T*B*Z = I;
90*> if ITYPE = 3, Z**T*inv(B)*Z = I.
91*> If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
92*> or the lower triangle (if UPLO='L') of A, including the
93*> diagonal, is destroyed.
94*> \endverbatim
95*>
96*> \param[in] LDA
97*> \verbatim
98*> LDA is INTEGER
99*> The leading dimension of the array A. LDA >= max(1,N).
100*> \endverbatim
101*>
102*> \param[in,out] B
103*> \verbatim
104*> B is REAL array, dimension (LDB, N)
105*> On entry, the symmetric positive definite matrix B.
106*> If UPLO = 'U', the leading N-by-N upper triangular part of B
107*> contains the upper triangular part of the matrix B.
108*> If UPLO = 'L', the leading N-by-N lower triangular part of B
109*> contains the lower triangular part of the matrix B.
110*>
111*> On exit, if INFO <= N, the part of B containing the matrix is
112*> overwritten by the triangular factor U or L from the Cholesky
113*> factorization B = U**T*U or B = L*L**T.
114*> \endverbatim
115*>
116*> \param[in] LDB
117*> \verbatim
118*> LDB is INTEGER
119*> The leading dimension of the array B. LDB >= max(1,N).
120*> \endverbatim
121*>
122*> \param[out] W
123*> \verbatim
124*> W is REAL array, dimension (N)
125*> If INFO = 0, the eigenvalues in ascending order.
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> WORK is REAL array, dimension (MAX(1,LWORK))
131*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
132*> \endverbatim
133*>
134*> \param[in] LWORK
135*> \verbatim
136*> LWORK is INTEGER
137*> The length of the array WORK. LWORK >= max(1,3*N-1).
138*> For optimal efficiency, LWORK >= (NB+2)*N,
139*> where NB is the blocksize for SSYTRD returned by ILAENV.
140*>
141*> If LWORK = -1, then a workspace query is assumed; the routine
142*> only calculates the optimal size of the WORK array, returns
143*> this value as the first entry of the WORK array, and no error
144*> message related to LWORK is issued by XERBLA.
145*> \endverbatim
146*>
147*> \param[out] INFO
148*> \verbatim
149*> INFO is INTEGER
150*> = 0: successful exit
151*> < 0: if INFO = -i, the i-th argument had an illegal value
152*> > 0: SPOTRF or SSYEV returned an error code:
153*> <= N: if INFO = i, SSYEV failed to converge;
154*> i off-diagonal elements of an intermediate
155*> tridiagonal form did not converge to zero;
156*> > N: if INFO = N + i, for 1 <= i <= N, then the leading
157*> principal minor of order i of B is not positive.
158*> The factorization of B could not be completed and
159*> no eigenvalues or eigenvectors were computed.
160*> \endverbatim
161*
162* Authors:
163* ========
164*
165*> \author Univ. of Tennessee
166*> \author Univ. of California Berkeley
167*> \author Univ. of Colorado Denver
168*> \author NAG Ltd.
169*
170*> \ingroup hegv
171*
172* =====================================================================
173 SUBROUTINE ssygv( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
174 $ LWORK, INFO )
175*
176* -- LAPACK driver routine --
177* -- LAPACK is a software package provided by Univ. of Tennessee, --
178* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
179*
180* .. Scalar Arguments ..
181 CHARACTER JOBZ, UPLO
182 INTEGER INFO, ITYPE, LDA, LDB, LWORK, N
183* ..
184* .. Array Arguments ..
185 REAL A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
186* ..
187*
188* =====================================================================
189*
190* .. Parameters ..
191 REAL ONE
192 parameter( one = 1.0e+0 )
193* ..
194* .. Local Scalars ..
195 LOGICAL LQUERY, UPPER, WANTZ
196 CHARACTER TRANS
197 INTEGER LWKMIN, LWKOPT, NB, NEIG
198* ..
199* .. External Functions ..
200 LOGICAL LSAME
201 INTEGER ILAENV
202 REAL SROUNDUP_LWORK
203 EXTERNAL ilaenv, lsame, sroundup_lwork
204* ..
205* .. External Subroutines ..
206 EXTERNAL spotrf, ssyev, ssygst, strmm, strsm, xerbla
207* ..
208* .. Intrinsic Functions ..
209 INTRINSIC max
210* ..
211* .. Executable Statements ..
212*
213* Test the input parameters.
214*
215 wantz = lsame( jobz, 'V' )
216 upper = lsame( uplo, 'U' )
217 lquery = ( lwork.EQ.-1 )
218*
219 info = 0
220 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
221 info = -1
222 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
223 info = -2
224 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
225 info = -3
226 ELSE IF( n.LT.0 ) THEN
227 info = -4
228 ELSE IF( lda.LT.max( 1, n ) ) THEN
229 info = -6
230 ELSE IF( ldb.LT.max( 1, n ) ) THEN
231 info = -8
232 END IF
233*
234 IF( info.EQ.0 ) THEN
235 lwkmin = max( 1, 3*n - 1 )
236 nb = ilaenv( 1, 'SSYTRD', uplo, n, -1, -1, -1 )
237 lwkopt = max( lwkmin, ( nb + 2 )*n )
238 work( 1 ) = sroundup_lwork(lwkopt)
239*
240 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
241 info = -11
242 END IF
243 END IF
244*
245 IF( info.NE.0 ) THEN
246 CALL xerbla( 'SSYGV ', -info )
247 RETURN
248 ELSE IF( lquery ) THEN
249 RETURN
250 END IF
251*
252* Quick return if possible
253*
254 IF( n.EQ.0 )
255 $ RETURN
256*
257* Form a Cholesky factorization of B.
258*
259 CALL spotrf( uplo, n, b, ldb, info )
260 IF( info.NE.0 ) THEN
261 info = n + info
262 RETURN
263 END IF
264*
265* Transform problem to standard eigenvalue problem and solve.
266*
267 CALL ssygst( itype, uplo, n, a, lda, b, ldb, info )
268 CALL ssyev( jobz, uplo, n, a, lda, w, work, lwork, info )
269*
270 IF( wantz ) THEN
271*
272* Backtransform eigenvectors to the original problem.
273*
274 neig = n
275 IF( info.GT.0 )
276 $ neig = info - 1
277 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
278*
279* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
280* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
281*
282 IF( upper ) THEN
283 trans = 'N'
284 ELSE
285 trans = 'T'
286 END IF
287*
288 CALL strsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
289 $ b, ldb, a, lda )
290*
291 ELSE IF( itype.EQ.3 ) THEN
292*
293* For B*A*x=(lambda)*x;
294* backtransform eigenvectors: x = L*y or U**T*y
295*
296 IF( upper ) THEN
297 trans = 'T'
298 ELSE
299 trans = 'N'
300 END IF
301*
302 CALL strmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
303 $ b, ldb, a, lda )
304 END IF
305 END IF
306*
307 work( 1 ) = sroundup_lwork(lwkopt)
308 RETURN
309*
310* End of SSYGV
311*
312 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssyev(jobz, uplo, n, a, lda, w, work, lwork, info)
SSYEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Definition ssyev.f:132
subroutine ssygst(itype, uplo, n, a, lda, b, ldb, info)
SSYGST
Definition ssygst.f:127
subroutine ssygv(itype, jobz, uplo, n, a, lda, b, ldb, w, work, lwork, info)
SSYGV
Definition ssygv.f:175
subroutine spotrf(uplo, n, a, lda, info)
SPOTRF
Definition spotrf.f:107
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