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