LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsygv.f
Go to the documentation of this file.
1*> \brief \b DSYGV
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DSYGV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsygv.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsygv.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsygv.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DSYGV( 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* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> DSYGV 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
125*> If INFO = 0, the eigenvalues in ascending order.
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> WORK is DOUBLE PRECISION 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 DSYTRD 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: DPOTRF or DSYEV returned an error code:
153*> <= N: if INFO = i, DSYEV 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 dsygv( 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 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
186* ..
187*
188* =====================================================================
189*
190* .. Parameters ..
191 DOUBLE PRECISION ONE
192 parameter( one = 1.0d+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 EXTERNAL lsame, ilaenv
203* ..
204* .. External Subroutines ..
205 EXTERNAL dpotrf, dsyev, dsygst, dtrmm, dtrsm, xerbla
206* ..
207* .. Intrinsic Functions ..
208 INTRINSIC max
209* ..
210* .. Executable Statements ..
211*
212* Test the input parameters.
213*
214 wantz = lsame( jobz, 'V' )
215 upper = lsame( uplo, 'U' )
216 lquery = ( lwork.EQ.-1 )
217*
218 info = 0
219 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
220 info = -1
221 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
222 info = -2
223 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
224 info = -3
225 ELSE IF( n.LT.0 ) THEN
226 info = -4
227 ELSE IF( lda.LT.max( 1, n ) ) THEN
228 info = -6
229 ELSE IF( ldb.LT.max( 1, n ) ) THEN
230 info = -8
231 END IF
232*
233 IF( info.EQ.0 ) THEN
234 lwkmin = max( 1, 3*n - 1 )
235 nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
236 lwkopt = max( lwkmin, ( nb + 2 )*n )
237 work( 1 ) = lwkopt
238*
239 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
240 info = -11
241 END IF
242 END IF
243*
244 IF( info.NE.0 ) THEN
245 CALL xerbla( 'DSYGV ', -info )
246 RETURN
247 ELSE IF( lquery ) THEN
248 RETURN
249 END IF
250*
251* Quick return if possible
252*
253 IF( n.EQ.0 )
254 $ RETURN
255*
256* Form a Cholesky factorization of B.
257*
258 CALL dpotrf( uplo, n, b, ldb, info )
259 IF( info.NE.0 ) THEN
260 info = n + info
261 RETURN
262 END IF
263*
264* Transform problem to standard eigenvalue problem and solve.
265*
266 CALL dsygst( itype, uplo, n, a, lda, b, ldb, info )
267 CALL dsyev( jobz, uplo, n, a, lda, w, work, lwork, info )
268*
269 IF( wantz ) THEN
270*
271* Backtransform eigenvectors to the original problem.
272*
273 neig = n
274 IF( info.GT.0 )
275 $ neig = info - 1
276 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
277*
278* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
279* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
280*
281 IF( upper ) THEN
282 trans = 'N'
283 ELSE
284 trans = 'T'
285 END IF
286*
287 CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
288 $ b, ldb, a, lda )
289*
290 ELSE IF( itype.EQ.3 ) THEN
291*
292* For B*A*x=(lambda)*x;
293* backtransform eigenvectors: x = L*y or U**T*y
294*
295 IF( upper ) THEN
296 trans = 'T'
297 ELSE
298 trans = 'N'
299 END IF
300*
301 CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
302 $ b, ldb, a, lda )
303 END IF
304 END IF
305*
306 work( 1 ) = lwkopt
307 RETURN
308*
309* End of DSYGV
310*
311 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsyev(jobz, uplo, n, a, lda, w, work, lwork, info)
DSYEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Definition dsyev.f:132
subroutine dsygst(itype, uplo, n, a, lda, b, ldb, info)
DSYGST
Definition dsygst.f:127
subroutine dsygv(itype, jobz, uplo, n, a, lda, b, ldb, w, work, lwork, info)
DSYGV
Definition dsygv.f:175
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:107
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181