LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dsbgvx.f
Go to the documentation of this file.
1 *> \brief \b DSBGVX
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DSBGVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbgvx.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbgvx.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbgvx.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
22 * LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
23 * LDZ, WORK, IWORK, IFAIL, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER JOBZ, RANGE, UPLO
27 * INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
28 * $ N
29 * DOUBLE PRECISION ABSTOL, VL, VU
30 * ..
31 * .. Array Arguments ..
32 * INTEGER IFAIL( * ), IWORK( * )
33 * DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
34 * $ W( * ), WORK( * ), Z( LDZ, * )
35 * ..
36 *
37 *
38 *> \par Purpose:
39 * =============
40 *>
41 *> \verbatim
42 *>
43 *> DSBGVX computes selected eigenvalues, and optionally, eigenvectors
44 *> of a real generalized symmetric-definite banded eigenproblem, of
45 *> the form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric
46 *> and banded, and B is also positive definite. Eigenvalues and
47 *> eigenvectors can be selected by specifying either all eigenvalues,
48 *> a range of values or a range of indices for the desired eigenvalues.
49 *> \endverbatim
50 *
51 * Arguments:
52 * ==========
53 *
54 *> \param[in] JOBZ
55 *> \verbatim
56 *> JOBZ is CHARACTER*1
57 *> = 'N': Compute eigenvalues only;
58 *> = 'V': Compute eigenvalues and eigenvectors.
59 *> \endverbatim
60 *>
61 *> \param[in] RANGE
62 *> \verbatim
63 *> RANGE is CHARACTER*1
64 *> = 'A': all eigenvalues will be found.
65 *> = 'V': all eigenvalues in the half-open interval (VL,VU]
66 *> will be found.
67 *> = 'I': the IL-th through IU-th eigenvalues will be found.
68 *> \endverbatim
69 *>
70 *> \param[in] UPLO
71 *> \verbatim
72 *> UPLO is CHARACTER*1
73 *> = 'U': Upper triangles of A and B are stored;
74 *> = 'L': Lower triangles of A and B are stored.
75 *> \endverbatim
76 *>
77 *> \param[in] N
78 *> \verbatim
79 *> N is INTEGER
80 *> The order of the matrices A and B. N >= 0.
81 *> \endverbatim
82 *>
83 *> \param[in] KA
84 *> \verbatim
85 *> KA is INTEGER
86 *> The number of superdiagonals of the matrix A if UPLO = 'U',
87 *> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
88 *> \endverbatim
89 *>
90 *> \param[in] KB
91 *> \verbatim
92 *> KB is INTEGER
93 *> The number of superdiagonals of the matrix B if UPLO = 'U',
94 *> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
95 *> \endverbatim
96 *>
97 *> \param[in,out] AB
98 *> \verbatim
99 *> AB is DOUBLE PRECISION array, dimension (LDAB, N)
100 *> On entry, the upper or lower triangle of the symmetric band
101 *> matrix A, stored in the first ka+1 rows of the array. The
102 *> j-th column of A is stored in the j-th column of the array AB
103 *> as follows:
104 *> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
105 *> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
106 *>
107 *> On exit, the contents of AB are destroyed.
108 *> \endverbatim
109 *>
110 *> \param[in] LDAB
111 *> \verbatim
112 *> LDAB is INTEGER
113 *> The leading dimension of the array AB. LDAB >= KA+1.
114 *> \endverbatim
115 *>
116 *> \param[in,out] BB
117 *> \verbatim
118 *> BB is DOUBLE PRECISION array, dimension (LDBB, N)
119 *> On entry, the upper or lower triangle of the symmetric band
120 *> matrix B, stored in the first kb+1 rows of the array. The
121 *> j-th column of B is stored in the j-th column of the array BB
122 *> as follows:
123 *> if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
124 *> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
125 *>
126 *> On exit, the factor S from the split Cholesky factorization
127 *> B = S**T*S, as returned by DPBSTF.
128 *> \endverbatim
129 *>
130 *> \param[in] LDBB
131 *> \verbatim
132 *> LDBB is INTEGER
133 *> The leading dimension of the array BB. LDBB >= KB+1.
134 *> \endverbatim
135 *>
136 *> \param[out] Q
137 *> \verbatim
138 *> Q is DOUBLE PRECISION array, dimension (LDQ, N)
139 *> If JOBZ = 'V', the n-by-n matrix used in the reduction of
140 *> A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
141 *> and consequently C to tridiagonal form.
142 *> If JOBZ = 'N', the array Q is not referenced.
143 *> \endverbatim
144 *>
145 *> \param[in] LDQ
146 *> \verbatim
147 *> LDQ is INTEGER
148 *> The leading dimension of the array Q. If JOBZ = 'N',
149 *> LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).
150 *> \endverbatim
151 *>
152 *> \param[in] VL
153 *> \verbatim
154 *> VL is DOUBLE PRECISION
155 *>
156 *> If RANGE='V', the lower bound of the interval to
157 *> be searched for eigenvalues. VL < VU.
158 *> Not referenced if RANGE = 'A' or 'I'.
159 *> \endverbatim
160 *>
161 *> \param[in] VU
162 *> \verbatim
163 *> VU is DOUBLE PRECISION
164 *>
165 *> If RANGE='V', the upper bound of the interval to
166 *> be searched for eigenvalues. VL < VU.
167 *> Not referenced if RANGE = 'A' or 'I'.
168 *> \endverbatim
169 *>
170 *> \param[in] IL
171 *> \verbatim
172 *> IL is INTEGER
173 *>
174 *> If RANGE='I', the index of the
175 *> smallest eigenvalue to be returned.
176 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
177 *> Not referenced if RANGE = 'A' or 'V'.
178 *> \endverbatim
179 *>
180 *> \param[in] IU
181 *> \verbatim
182 *> IU is INTEGER
183 *>
184 *> If RANGE='I', the index of the
185 *> largest eigenvalue to be returned.
186 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
187 *> Not referenced if RANGE = 'A' or 'V'.
188 *> \endverbatim
189 *>
190 *> \param[in] ABSTOL
191 *> \verbatim
192 *> ABSTOL is DOUBLE PRECISION
193 *> The absolute error tolerance for the eigenvalues.
194 *> An approximate eigenvalue is accepted as converged
195 *> when it is determined to lie in an interval [a,b]
196 *> of width less than or equal to
197 *>
198 *> ABSTOL + EPS * max( |a|,|b| ) ,
199 *>
200 *> where EPS is the machine precision. If ABSTOL is less than
201 *> or equal to zero, then EPS*|T| will be used in its place,
202 *> where |T| is the 1-norm of the tridiagonal matrix obtained
203 *> by reducing A to tridiagonal form.
204 *>
205 *> Eigenvalues will be computed most accurately when ABSTOL is
206 *> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
207 *> If this routine returns with INFO>0, indicating that some
208 *> eigenvectors did not converge, try setting ABSTOL to
209 *> 2*DLAMCH('S').
210 *> \endverbatim
211 *>
212 *> \param[out] M
213 *> \verbatim
214 *> M is INTEGER
215 *> The total number of eigenvalues found. 0 <= M <= N.
216 *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
217 *> \endverbatim
218 *>
219 *> \param[out] W
220 *> \verbatim
221 *> W is DOUBLE PRECISION array, dimension (N)
222 *> If INFO = 0, the eigenvalues in ascending order.
223 *> \endverbatim
224 *>
225 *> \param[out] Z
226 *> \verbatim
227 *> Z is DOUBLE PRECISION array, dimension (LDZ, N)
228 *> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
229 *> eigenvectors, with the i-th column of Z holding the
230 *> eigenvector associated with W(i). The eigenvectors are
231 *> normalized so Z**T*B*Z = I.
232 *> If JOBZ = 'N', then Z is not referenced.
233 *> \endverbatim
234 *>
235 *> \param[in] LDZ
236 *> \verbatim
237 *> LDZ is INTEGER
238 *> The leading dimension of the array Z. LDZ >= 1, and if
239 *> JOBZ = 'V', LDZ >= max(1,N).
240 *> \endverbatim
241 *>
242 *> \param[out] WORK
243 *> \verbatim
244 *> WORK is DOUBLE PRECISION array, dimension (7*N)
245 *> \endverbatim
246 *>
247 *> \param[out] IWORK
248 *> \verbatim
249 *> IWORK is INTEGER array, dimension (5*N)
250 *> \endverbatim
251 *>
252 *> \param[out] IFAIL
253 *> \verbatim
254 *> IFAIL is INTEGER array, dimension (M)
255 *> If JOBZ = 'V', then if INFO = 0, the first M elements of
256 *> IFAIL are zero. If INFO > 0, then IFAIL contains the
257 *> indices of the eigenvalues that failed to converge.
258 *> If JOBZ = 'N', then IFAIL is not referenced.
259 *> \endverbatim
260 *>
261 *> \param[out] INFO
262 *> \verbatim
263 *> INFO is INTEGER
264 *> = 0 : successful exit
265 *> < 0 : if INFO = -i, the i-th argument had an illegal value
266 *> <= N: if INFO = i, then i eigenvectors failed to converge.
267 *> Their indices are stored in IFAIL.
268 *> > N : DPBSTF returned an error code; i.e.,
269 *> if INFO = N + i, for 1 <= i <= N, then the leading
270 *> minor of order i of B is not positive definite.
271 *> The factorization of B could not be completed and
272 *> no eigenvalues or eigenvectors were computed.
273 *> \endverbatim
274 *
275 * Authors:
276 * ========
277 *
278 *> \author Univ. of Tennessee
279 *> \author Univ. of California Berkeley
280 *> \author Univ. of Colorado Denver
281 *> \author NAG Ltd.
282 *
283 *> \date June 2016
284 *
285 *> \ingroup doubleOTHEReigen
286 *
287 *> \par Contributors:
288 * ==================
289 *>
290 *> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
291 *
292 * =====================================================================
293  SUBROUTINE dsbgvx( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
294  $ ldbb, q, ldq, vl, vu, il, iu, abstol, m, w, z,
295  $ ldz, work, iwork, ifail, info )
296 *
297 * -- LAPACK driver routine (version 3.6.1) --
298 * -- LAPACK is a software package provided by Univ. of Tennessee, --
299 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
300 * June 2016
301 *
302 * .. Scalar Arguments ..
303  CHARACTER JOBZ, RANGE, UPLO
304  INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
305  $ n
306  DOUBLE PRECISION ABSTOL, VL, VU
307 * ..
308 * .. Array Arguments ..
309  INTEGER IFAIL( * ), IWORK( * )
310  DOUBLE PRECISION AB( ldab, * ), BB( ldbb, * ), Q( ldq, * ),
311  $ w( * ), work( * ), z( ldz, * )
312 * ..
313 *
314 * =====================================================================
315 *
316 * .. Parameters ..
317  DOUBLE PRECISION ZERO, ONE
318  parameter ( zero = 0.0d+0, one = 1.0d+0 )
319 * ..
320 * .. Local Scalars ..
321  LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
322  CHARACTER ORDER, VECT
323  INTEGER I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP,
324  $ indiwo, indwrk, itmp1, j, jj, nsplit
325  DOUBLE PRECISION TMP1
326 * ..
327 * .. External Functions ..
328  LOGICAL LSAME
329  EXTERNAL lsame
330 * ..
331 * .. External Subroutines ..
332  EXTERNAL dcopy, dgemv, dlacpy, dpbstf, dsbgst, dsbtrd,
334 * ..
335 * .. Intrinsic Functions ..
336  INTRINSIC min
337 * ..
338 * .. Executable Statements ..
339 *
340 * Test the input parameters.
341 *
342  wantz = lsame( jobz, 'V' )
343  upper = lsame( uplo, 'U' )
344  alleig = lsame( range, 'A' )
345  valeig = lsame( range, 'V' )
346  indeig = lsame( range, 'I' )
347 *
348  info = 0
349  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
350  info = -1
351  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
352  info = -2
353  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
354  info = -3
355  ELSE IF( n.LT.0 ) THEN
356  info = -4
357  ELSE IF( ka.LT.0 ) THEN
358  info = -5
359  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
360  info = -6
361  ELSE IF( ldab.LT.ka+1 ) THEN
362  info = -8
363  ELSE IF( ldbb.LT.kb+1 ) THEN
364  info = -10
365  ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
366  info = -12
367  ELSE
368  IF( valeig ) THEN
369  IF( n.GT.0 .AND. vu.LE.vl )
370  $ info = -14
371  ELSE IF( indeig ) THEN
372  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
373  info = -15
374  ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
375  info = -16
376  END IF
377  END IF
378  END IF
379  IF( info.EQ.0) THEN
380  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
381  info = -21
382  END IF
383  END IF
384 *
385  IF( info.NE.0 ) THEN
386  CALL xerbla( 'DSBGVX', -info )
387  RETURN
388  END IF
389 *
390 * Quick return if possible
391 *
392  m = 0
393  IF( n.EQ.0 )
394  $ RETURN
395 *
396 * Form a split Cholesky factorization of B.
397 *
398  CALL dpbstf( uplo, n, kb, bb, ldbb, info )
399  IF( info.NE.0 ) THEN
400  info = n + info
401  RETURN
402  END IF
403 *
404 * Transform problem to standard eigenvalue problem.
405 *
406  CALL dsbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
407  $ work, iinfo )
408 *
409 * Reduce symmetric band matrix to tridiagonal form.
410 *
411  indd = 1
412  inde = indd + n
413  indwrk = inde + n
414  IF( wantz ) THEN
415  vect = 'U'
416  ELSE
417  vect = 'N'
418  END IF
419  CALL dsbtrd( vect, uplo, n, ka, ab, ldab, work( indd ),
420  $ work( inde ), q, ldq, work( indwrk ), iinfo )
421 *
422 * If all eigenvalues are desired and ABSTOL is less than or equal
423 * to zero, then call DSTERF or SSTEQR. If this fails for some
424 * eigenvalue, then try DSTEBZ.
425 *
426  test = .false.
427  IF( indeig ) THEN
428  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
429  test = .true.
430  END IF
431  END IF
432  IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
433  CALL dcopy( n, work( indd ), 1, w, 1 )
434  indee = indwrk + 2*n
435  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
436  IF( .NOT.wantz ) THEN
437  CALL dsterf( n, w, work( indee ), info )
438  ELSE
439  CALL dlacpy( 'A', n, n, q, ldq, z, ldz )
440  CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
441  $ work( indwrk ), info )
442  IF( info.EQ.0 ) THEN
443  DO 10 i = 1, n
444  ifail( i ) = 0
445  10 CONTINUE
446  END IF
447  END IF
448  IF( info.EQ.0 ) THEN
449  m = n
450  GO TO 30
451  END IF
452  info = 0
453  END IF
454 *
455 * Otherwise, call DSTEBZ and, if eigenvectors are desired,
456 * call DSTEIN.
457 *
458  IF( wantz ) THEN
459  order = 'B'
460  ELSE
461  order = 'E'
462  END IF
463  indibl = 1
464  indisp = indibl + n
465  indiwo = indisp + n
466  CALL dstebz( range, order, n, vl, vu, il, iu, abstol,
467  $ work( indd ), work( inde ), m, nsplit, w,
468  $ iwork( indibl ), iwork( indisp ), work( indwrk ),
469  $ iwork( indiwo ), info )
470 *
471  IF( wantz ) THEN
472  CALL dstein( n, work( indd ), work( inde ), m, w,
473  $ iwork( indibl ), iwork( indisp ), z, ldz,
474  $ work( indwrk ), iwork( indiwo ), ifail, info )
475 *
476 * Apply transformation matrix used in reduction to tridiagonal
477 * form to eigenvectors returned by DSTEIN.
478 *
479  DO 20 j = 1, m
480  CALL dcopy( n, z( 1, j ), 1, work( 1 ), 1 )
481  CALL dgemv( 'N', n, n, one, q, ldq, work, 1, zero,
482  $ z( 1, j ), 1 )
483  20 CONTINUE
484  END IF
485 *
486  30 CONTINUE
487 *
488 * If eigenvalues are not in order, then sort them, along with
489 * eigenvectors.
490 *
491  IF( wantz ) THEN
492  DO 50 j = 1, m - 1
493  i = 0
494  tmp1 = w( j )
495  DO 40 jj = j + 1, m
496  IF( w( jj ).LT.tmp1 ) THEN
497  i = jj
498  tmp1 = w( jj )
499  END IF
500  40 CONTINUE
501 *
502  IF( i.NE.0 ) THEN
503  itmp1 = iwork( indibl+i-1 )
504  w( i ) = w( j )
505  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
506  w( j ) = tmp1
507  iwork( indibl+j-1 ) = itmp1
508  CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
509  IF( info.NE.0 ) THEN
510  itmp1 = ifail( i )
511  ifail( i ) = ifail( j )
512  ifail( j ) = itmp1
513  END IF
514  END IF
515  50 CONTINUE
516  END IF
517 *
518  RETURN
519 *
520 * End of DSBGVX
521 *
522  END
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine dsbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
DSBTRD
Definition: dsbtrd.f:165
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dpbstf(UPLO, N, KD, AB, LDAB, INFO)
DPBSTF
Definition: dpbstf.f:154
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
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 dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dsbgvx(JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSBGVX
Definition: dsbgvx.f:296
subroutine dsbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, INFO)
DSBGST
Definition: dsbgst.f:161
subroutine dstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEIN
Definition: dstein.f:176