LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsbgvd.f
Go to the documentation of this file.
1*> \brief \b DSBGVD
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DSBGVD + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbgvd.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbgvd.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbgvd.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DSBGVD( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
20* Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER JOBZ, UPLO
24* INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LWORK, N
25* ..
26* .. Array Arguments ..
27* INTEGER IWORK( * )
28* DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), W( * ),
29* $ WORK( * ), Z( LDZ, * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> DSBGVD computes all the eigenvalues, and optionally, the eigenvectors
39*> of a real generalized symmetric-definite banded eigenproblem, of the
40*> form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric and
41*> banded, and B is also positive definite. If eigenvectors are
42*> desired, it uses a divide and conquer algorithm.
43*>
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] JOBZ
50*> \verbatim
51*> JOBZ is CHARACTER*1
52*> = 'N': Compute eigenvalues only;
53*> = 'V': Compute eigenvalues and eigenvectors.
54*> \endverbatim
55*>
56*> \param[in] UPLO
57*> \verbatim
58*> UPLO is CHARACTER*1
59*> = 'U': Upper triangles of A and B are stored;
60*> = 'L': Lower triangles of A and B are stored.
61*> \endverbatim
62*>
63*> \param[in] N
64*> \verbatim
65*> N is INTEGER
66*> The order of the matrices A and B. N >= 0.
67*> \endverbatim
68*>
69*> \param[in] KA
70*> \verbatim
71*> KA is INTEGER
72*> The number of superdiagonals of the matrix A if UPLO = 'U',
73*> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
74*> \endverbatim
75*>
76*> \param[in] KB
77*> \verbatim
78*> KB is INTEGER
79*> The number of superdiagonals of the matrix B if UPLO = 'U',
80*> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
81*> \endverbatim
82*>
83*> \param[in,out] AB
84*> \verbatim
85*> AB is DOUBLE PRECISION array, dimension (LDAB, N)
86*> On entry, the upper or lower triangle of the symmetric band
87*> matrix A, stored in the first ka+1 rows of the array. The
88*> j-th column of A is stored in the j-th column of the array AB
89*> as follows:
90*> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
91*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
92*>
93*> On exit, the contents of AB are destroyed.
94*> \endverbatim
95*>
96*> \param[in] LDAB
97*> \verbatim
98*> LDAB is INTEGER
99*> The leading dimension of the array AB. LDAB >= KA+1.
100*> \endverbatim
101*>
102*> \param[in,out] BB
103*> \verbatim
104*> BB is DOUBLE PRECISION array, dimension (LDBB, N)
105*> On entry, the upper or lower triangle of the symmetric band
106*> matrix B, stored in the first kb+1 rows of the array. The
107*> j-th column of B is stored in the j-th column of the array BB
108*> as follows:
109*> if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
110*> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
111*>
112*> On exit, the factor S from the split Cholesky factorization
113*> B = S**T*S, as returned by DPBSTF.
114*> \endverbatim
115*>
116*> \param[in] LDBB
117*> \verbatim
118*> LDBB is INTEGER
119*> The leading dimension of the array BB. LDBB >= KB+1.
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] Z
129*> \verbatim
130*> Z is DOUBLE PRECISION array, dimension (LDZ, N)
131*> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
132*> eigenvectors, with the i-th column of Z holding the
133*> eigenvector associated with W(i). The eigenvectors are
134*> normalized so Z**T*B*Z = I.
135*> If JOBZ = 'N', then Z is not referenced.
136*> \endverbatim
137*>
138*> \param[in] LDZ
139*> \verbatim
140*> LDZ is INTEGER
141*> The leading dimension of the array Z. LDZ >= 1, and if
142*> JOBZ = 'V', LDZ >= max(1,N).
143*> \endverbatim
144*>
145*> \param[out] WORK
146*> \verbatim
147*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
148*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
149*> \endverbatim
150*>
151*> \param[in] LWORK
152*> \verbatim
153*> LWORK is INTEGER
154*> The dimension of the array WORK.
155*> If N <= 1, LWORK >= 1.
156*> If JOBZ = 'N' and N > 1, LWORK >= 2*N.
157*> If JOBZ = 'V' and N > 1, LWORK >= 1 + 5*N + 2*N**2.
158*>
159*> If LWORK = -1, then a workspace query is assumed; the routine
160*> only calculates the optimal sizes of the WORK and IWORK
161*> arrays, returns these values as the first entries of the WORK
162*> and IWORK arrays, and no error message related to LWORK or
163*> LIWORK is issued by XERBLA.
164*> \endverbatim
165*>
166*> \param[out] IWORK
167*> \verbatim
168*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
169*> On exit, if LIWORK > 0, IWORK(1) returns the optimal LIWORK.
170*> \endverbatim
171*>
172*> \param[in] LIWORK
173*> \verbatim
174*> LIWORK is INTEGER
175*> The dimension of the array IWORK.
176*> If JOBZ = 'N' or N <= 1, LIWORK >= 1.
177*> If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.
178*>
179*> If LIWORK = -1, then a workspace query is assumed; the
180*> routine only calculates the optimal sizes of the WORK and
181*> IWORK arrays, returns these values as the first entries of
182*> the WORK and IWORK arrays, and no error message related to
183*> LWORK or LIWORK is issued by XERBLA.
184*> \endverbatim
185*>
186*> \param[out] INFO
187*> \verbatim
188*> INFO is INTEGER
189*> = 0: successful exit
190*> < 0: if INFO = -i, the i-th argument had an illegal value
191*> > 0: if INFO = i, and i is:
192*> <= N: the algorithm failed to converge:
193*> i off-diagonal elements of an intermediate
194*> tridiagonal form did not converge to zero;
195*> > N: if INFO = N + i, for 1 <= i <= N, then DPBSTF
196*> returned INFO = i: B is not positive definite.
197*> The factorization of B could not be completed and
198*> no eigenvalues or eigenvectors were computed.
199*> \endverbatim
200*
201* Authors:
202* ========
203*
204*> \author Univ. of Tennessee
205*> \author Univ. of California Berkeley
206*> \author Univ. of Colorado Denver
207*> \author NAG Ltd.
208*
209*> \ingroup hbgvd
210*
211*> \par Contributors:
212* ==================
213*>
214*> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
215*
216* =====================================================================
217 SUBROUTINE dsbgvd( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB,
218 $ W,
219 $ Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO )
220*
221* -- LAPACK driver routine --
222* -- LAPACK is a software package provided by Univ. of Tennessee, --
223* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
224*
225* .. Scalar Arguments ..
226 CHARACTER JOBZ, UPLO
227 INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LWORK, N
228* ..
229* .. Array Arguments ..
230 INTEGER IWORK( * )
231 DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), W( * ),
232 $ work( * ), z( ldz, * )
233* ..
234*
235* =====================================================================
236*
237* .. Parameters ..
238 DOUBLE PRECISION ONE, ZERO
239 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
240* ..
241* .. Local Scalars ..
242 LOGICAL LQUERY, UPPER, WANTZ
243 CHARACTER VECT
244 INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLWRK2,
245 $ lwmin
246* ..
247* .. External Functions ..
248 LOGICAL LSAME
249 EXTERNAL LSAME
250* ..
251* .. External Subroutines ..
252 EXTERNAL dgemm, dlacpy, dpbstf, dsbgst, dsbtrd,
253 $ dstedc,
254 $ dsterf, xerbla
255* ..
256* .. Executable Statements ..
257*
258* Test the input parameters.
259*
260 wantz = lsame( jobz, 'V' )
261 upper = lsame( uplo, 'U' )
262 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
263*
264 info = 0
265 IF( n.LE.1 ) THEN
266 liwmin = 1
267 lwmin = 1
268 ELSE IF( wantz ) THEN
269 liwmin = 3 + 5*n
270 lwmin = 1 + 5*n + 2*n**2
271 ELSE
272 liwmin = 1
273 lwmin = 2*n
274 END IF
275*
276 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
277 info = -1
278 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
279 info = -2
280 ELSE IF( n.LT.0 ) THEN
281 info = -3
282 ELSE IF( ka.LT.0 ) THEN
283 info = -4
284 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
285 info = -5
286 ELSE IF( ldab.LT.ka+1 ) THEN
287 info = -7
288 ELSE IF( ldbb.LT.kb+1 ) THEN
289 info = -9
290 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
291 info = -12
292 END IF
293*
294 IF( info.EQ.0 ) THEN
295 work( 1 ) = lwmin
296 iwork( 1 ) = liwmin
297*
298 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
299 info = -14
300 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
301 info = -16
302 END IF
303 END IF
304*
305 IF( info.NE.0 ) THEN
306 CALL xerbla( 'DSBGVD', -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 split Cholesky factorization of B.
318*
319 CALL dpbstf( uplo, n, kb, bb, ldbb, info )
320 IF( info.NE.0 ) THEN
321 info = n + info
322 RETURN
323 END IF
324*
325* Transform problem to standard eigenvalue problem.
326*
327 inde = 1
328 indwrk = inde + n
329 indwk2 = indwrk + n*n
330 llwrk2 = lwork - indwk2 + 1
331 CALL dsbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
332 $ work, iinfo )
333*
334* Reduce to tridiagonal form.
335*
336 IF( wantz ) THEN
337 vect = 'U'
338 ELSE
339 vect = 'N'
340 END IF
341 CALL dsbtrd( vect, uplo, n, ka, ab, ldab, w, work( inde ), z,
342 $ ldz,
343 $ work( indwrk ), iinfo )
344*
345* For eigenvalues only, call DSTERF. For eigenvectors, call SSTEDC.
346*
347 IF( .NOT.wantz ) THEN
348 CALL dsterf( n, w, work( inde ), info )
349 ELSE
350 CALL dstedc( 'I', n, w, work( inde ), work( indwrk ), n,
351 $ work( indwk2 ), llwrk2, iwork, liwork, info )
352 CALL dgemm( 'N', 'N', n, n, n, one, z, ldz, work( indwrk ),
353 $ n,
354 $ zero, work( indwk2 ), n )
355 CALL dlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
356 END IF
357*
358 work( 1 ) = lwmin
359 iwork( 1 ) = liwmin
360*
361 RETURN
362*
363* End of DSBGVD
364*
365 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dsbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, info)
DSBGST
Definition dsbgst.f:158
subroutine dsbgvd(jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, w, z, ldz, work, lwork, iwork, liwork, info)
DSBGVD
Definition dsbgvd.f:220
subroutine dsbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
DSBTRD
Definition dsbtrd.f:161
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
subroutine dpbstf(uplo, n, kd, ab, ldab, info)
DPBSTF
Definition dpbstf.f:150
subroutine dstedc(compz, n, d, e, z, ldz, work, lwork, iwork, liwork, info)
DSTEDC
Definition dstedc.f:180
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84