LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dsbgv.f
Go to the documentation of this file.
1 *> \brief \b DSBGV
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DSBGV + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbgv.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbgv.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbgv.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DSBGV( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z,
22 * LDZ, WORK, INFO )
23 *
24 * .. Scalar Arguments ..
25 * CHARACTER JOBZ, UPLO
26 * INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, N
27 * ..
28 * .. Array Arguments ..
29 * DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), W( * ),
30 * $ WORK( * ), Z( LDZ, * )
31 * ..
32 *
33 *
34 *> \par Purpose:
35 * =============
36 *>
37 *> \verbatim
38 *>
39 *> DSBGV computes all the eigenvalues, and optionally, the eigenvectors
40 *> of a real generalized symmetric-definite banded eigenproblem, of
41 *> the form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric
42 *> and banded, and B is also positive definite.
43 *> \endverbatim
44 *
45 * Arguments:
46 * ==========
47 *
48 *> \param[in] JOBZ
49 *> \verbatim
50 *> JOBZ is CHARACTER*1
51 *> = 'N': Compute eigenvalues only;
52 *> = 'V': Compute eigenvalues and eigenvectors.
53 *> \endverbatim
54 *>
55 *> \param[in] UPLO
56 *> \verbatim
57 *> UPLO is CHARACTER*1
58 *> = 'U': Upper triangles of A and B are stored;
59 *> = 'L': Lower triangles of A and B are stored.
60 *> \endverbatim
61 *>
62 *> \param[in] N
63 *> \verbatim
64 *> N is INTEGER
65 *> The order of the matrices A and B. N >= 0.
66 *> \endverbatim
67 *>
68 *> \param[in] KA
69 *> \verbatim
70 *> KA is INTEGER
71 *> The number of superdiagonals of the matrix A if UPLO = 'U',
72 *> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
73 *> \endverbatim
74 *>
75 *> \param[in] KB
76 *> \verbatim
77 *> KB is INTEGER
78 *> The number of superdiagonals of the matrix B if UPLO = 'U',
79 *> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
80 *> \endverbatim
81 *>
82 *> \param[in,out] AB
83 *> \verbatim
84 *> AB is DOUBLE PRECISION array, dimension (LDAB, N)
85 *> On entry, the upper or lower triangle of the symmetric band
86 *> matrix A, stored in the first ka+1 rows of the array. The
87 *> j-th column of A is stored in the j-th column of the array AB
88 *> as follows:
89 *> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
90 *> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
91 *>
92 *> On exit, the contents of AB are destroyed.
93 *> \endverbatim
94 *>
95 *> \param[in] LDAB
96 *> \verbatim
97 *> LDAB is INTEGER
98 *> The leading dimension of the array AB. LDAB >= KA+1.
99 *> \endverbatim
100 *>
101 *> \param[in,out] BB
102 *> \verbatim
103 *> BB is DOUBLE PRECISION array, dimension (LDBB, N)
104 *> On entry, the upper or lower triangle of the symmetric band
105 *> matrix B, stored in the first kb+1 rows of the array. The
106 *> j-th column of B is stored in the j-th column of the array BB
107 *> as follows:
108 *> if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
109 *> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
110 *>
111 *> On exit, the factor S from the split Cholesky factorization
112 *> B = S**T*S, as returned by DPBSTF.
113 *> \endverbatim
114 *>
115 *> \param[in] LDBB
116 *> \verbatim
117 *> LDBB is INTEGER
118 *> The leading dimension of the array BB. LDBB >= KB+1.
119 *> \endverbatim
120 *>
121 *> \param[out] W
122 *> \verbatim
123 *> W is DOUBLE PRECISION array, dimension (N)
124 *> If INFO = 0, the eigenvalues in ascending order.
125 *> \endverbatim
126 *>
127 *> \param[out] Z
128 *> \verbatim
129 *> Z is DOUBLE PRECISION array, dimension (LDZ, N)
130 *> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
131 *> eigenvectors, with the i-th column of Z holding the
132 *> eigenvector associated with W(i). The eigenvectors are
133 *> normalized so that Z**T*B*Z = I.
134 *> If JOBZ = 'N', then Z is not referenced.
135 *> \endverbatim
136 *>
137 *> \param[in] LDZ
138 *> \verbatim
139 *> LDZ is INTEGER
140 *> The leading dimension of the array Z. LDZ >= 1, and if
141 *> JOBZ = 'V', LDZ >= N.
142 *> \endverbatim
143 *>
144 *> \param[out] WORK
145 *> \verbatim
146 *> WORK is DOUBLE PRECISION array, dimension (3*N)
147 *> \endverbatim
148 *>
149 *> \param[out] INFO
150 *> \verbatim
151 *> INFO is INTEGER
152 *> = 0: successful exit
153 *> < 0: if INFO = -i, the i-th argument had an illegal value
154 *> > 0: if INFO = i, and i is:
155 *> <= N: the algorithm failed to converge:
156 *> i off-diagonal elements of an intermediate
157 *> tridiagonal form did not converge to zero;
158 *> > N: if INFO = N + i, for 1 <= i <= N, then DPBSTF
159 *> returned INFO = i: B is not positive definite.
160 *> The factorization of B could not be completed and
161 *> no eigenvalues or eigenvectors were computed.
162 *> \endverbatim
163 *
164 * Authors:
165 * ========
166 *
167 *> \author Univ. of Tennessee
168 *> \author Univ. of California Berkeley
169 *> \author Univ. of Colorado Denver
170 *> \author NAG Ltd.
171 *
172 *> \date November 2015
173 *
174 *> \ingroup doubleOTHEReigen
175 *
176 * =====================================================================
177  SUBROUTINE dsbgv( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z,
178  $ ldz, work, info )
179 *
180 * -- LAPACK driver routine (version 3.6.0) --
181 * -- LAPACK is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * November 2015
184 *
185 * .. Scalar Arguments ..
186  CHARACTER JOBZ, UPLO
187  INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, N
188 * ..
189 * .. Array Arguments ..
190  DOUBLE PRECISION AB( ldab, * ), BB( ldbb, * ), W( * ),
191  $ work( * ), z( ldz, * )
192 * ..
193 *
194 * =====================================================================
195 *
196 * .. Local Scalars ..
197  LOGICAL UPPER, WANTZ
198  CHARACTER VECT
199  INTEGER IINFO, INDE, INDWRK
200 * ..
201 * .. External Functions ..
202  LOGICAL LSAME
203  EXTERNAL lsame
204 * ..
205 * .. External Subroutines ..
206  EXTERNAL dpbstf, dsbgst, dsbtrd, dsteqr, dsterf, xerbla
207 * ..
208 * .. Executable Statements ..
209 *
210 * Test the input parameters.
211 *
212  wantz = lsame( jobz, 'V' )
213  upper = lsame( uplo, 'U' )
214 *
215  info = 0
216  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
217  info = -1
218  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
219  info = -2
220  ELSE IF( n.LT.0 ) THEN
221  info = -3
222  ELSE IF( ka.LT.0 ) THEN
223  info = -4
224  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
225  info = -5
226  ELSE IF( ldab.LT.ka+1 ) THEN
227  info = -7
228  ELSE IF( ldbb.LT.kb+1 ) THEN
229  info = -9
230  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
231  info = -12
232  END IF
233  IF( info.NE.0 ) THEN
234  CALL xerbla( 'DSBGV ', -info )
235  RETURN
236  END IF
237 *
238 * Quick return if possible
239 *
240  IF( n.EQ.0 )
241  $ RETURN
242 *
243 * Form a split Cholesky factorization of B.
244 *
245  CALL dpbstf( uplo, n, kb, bb, ldbb, info )
246  IF( info.NE.0 ) THEN
247  info = n + info
248  RETURN
249  END IF
250 *
251 * Transform problem to standard eigenvalue problem.
252 *
253  inde = 1
254  indwrk = inde + n
255  CALL dsbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
256  $ work( indwrk ), iinfo )
257 *
258 * Reduce to tridiagonal form.
259 *
260  IF( wantz ) THEN
261  vect = 'U'
262  ELSE
263  vect = 'N'
264  END IF
265  CALL dsbtrd( vect, uplo, n, ka, ab, ldab, w, work( inde ), z, ldz,
266  $ work( indwrk ), iinfo )
267 *
268 * For eigenvalues only, call DSTERF. For eigenvectors, call SSTEQR.
269 *
270  IF( .NOT.wantz ) THEN
271  CALL dsterf( n, w, work( inde ), info )
272  ELSE
273  CALL dsteqr( jobz, n, w, work( inde ), z, ldz, work( indwrk ),
274  $ info )
275  END IF
276  RETURN
277 *
278 * End of DSBGV
279 *
280  END
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine dsbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
DSBTRD
Definition: dsbtrd.f:165
subroutine dpbstf(UPLO, N, KD, AB, LDAB, INFO)
DPBSTF
Definition: dpbstf.f:154
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:133
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dsbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, INFO)
DSBGST
Definition: dsbgst.f:161
subroutine dsbgv(JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z, LDZ, WORK, INFO)
DSBGV
Definition: dsbgv.f:179