LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsbev.f
Go to the documentation of this file.
1*> \brief <b> DSBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DSBEV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbev.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbev.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbev.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DSBEV( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK,
22* INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER JOBZ, UPLO
26* INTEGER INFO, KD, LDAB, LDZ, N
27* ..
28* .. Array Arguments ..
29* DOUBLE PRECISION AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> DSBEV computes all the eigenvalues and, optionally, eigenvectors of
39*> a real symmetric band matrix A.
40*> \endverbatim
41*
42* Arguments:
43* ==========
44*
45*> \param[in] JOBZ
46*> \verbatim
47*> JOBZ is CHARACTER*1
48*> = 'N': Compute eigenvalues only;
49*> = 'V': Compute eigenvalues and eigenvectors.
50*> \endverbatim
51*>
52*> \param[in] UPLO
53*> \verbatim
54*> UPLO is CHARACTER*1
55*> = 'U': Upper triangle of A is stored;
56*> = 'L': Lower triangle of A is stored.
57*> \endverbatim
58*>
59*> \param[in] N
60*> \verbatim
61*> N is INTEGER
62*> The order of the matrix A. N >= 0.
63*> \endverbatim
64*>
65*> \param[in] KD
66*> \verbatim
67*> KD is INTEGER
68*> The number of superdiagonals of the matrix A if UPLO = 'U',
69*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
70*> \endverbatim
71*>
72*> \param[in,out] AB
73*> \verbatim
74*> AB is DOUBLE PRECISION array, dimension (LDAB, N)
75*> On entry, the upper or lower triangle of the symmetric band
76*> matrix A, stored in the first KD+1 rows of the array. The
77*> j-th column of A is stored in the j-th column of the array AB
78*> as follows:
79*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
80*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
81*>
82*> On exit, AB is overwritten by values generated during the
83*> reduction to tridiagonal form. If UPLO = 'U', the first
84*> superdiagonal and the diagonal of the tridiagonal matrix T
85*> are returned in rows KD and KD+1 of AB, and if UPLO = 'L',
86*> the diagonal and first subdiagonal of T are returned in the
87*> first two rows of AB.
88*> \endverbatim
89*>
90*> \param[in] LDAB
91*> \verbatim
92*> LDAB is INTEGER
93*> The leading dimension of the array AB. LDAB >= KD + 1.
94*> \endverbatim
95*>
96*> \param[out] W
97*> \verbatim
98*> W is DOUBLE PRECISION array, dimension (N)
99*> If INFO = 0, the eigenvalues in ascending order.
100*> \endverbatim
101*>
102*> \param[out] Z
103*> \verbatim
104*> Z is DOUBLE PRECISION array, dimension (LDZ, N)
105*> If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
106*> eigenvectors of the matrix A, with the i-th column of Z
107*> holding the eigenvector associated with W(i).
108*> If JOBZ = 'N', then Z is not referenced.
109*> \endverbatim
110*>
111*> \param[in] LDZ
112*> \verbatim
113*> LDZ is INTEGER
114*> The leading dimension of the array Z. LDZ >= 1, and if
115*> JOBZ = 'V', LDZ >= max(1,N).
116*> \endverbatim
117*>
118*> \param[out] WORK
119*> \verbatim
120*> WORK is DOUBLE PRECISION array, dimension (max(1,3*N-2))
121*> \endverbatim
122*>
123*> \param[out] INFO
124*> \verbatim
125*> INFO is INTEGER
126*> = 0: successful exit
127*> < 0: if INFO = -i, the i-th argument had an illegal value
128*> > 0: if INFO = i, the algorithm failed to converge; i
129*> off-diagonal elements of an intermediate tridiagonal
130*> form did not converge to zero.
131*> \endverbatim
132*
133* Authors:
134* ========
135*
136*> \author Univ. of Tennessee
137*> \author Univ. of California Berkeley
138*> \author Univ. of Colorado Denver
139*> \author NAG Ltd.
140*
141*> \ingroup hbev
142*
143* =====================================================================
144 SUBROUTINE dsbev( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK,
145 $ INFO )
146*
147* -- LAPACK driver routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 CHARACTER JOBZ, UPLO
153 INTEGER INFO, KD, LDAB, LDZ, N
154* ..
155* .. Array Arguments ..
156 DOUBLE PRECISION AB( LDAB, * ), W( * ), WORK( * ), Z( LDZ, * )
157* ..
158*
159* =====================================================================
160*
161* .. Parameters ..
162 DOUBLE PRECISION ZERO, ONE
163 parameter( zero = 0.0d0, one = 1.0d0 )
164* ..
165* .. Local Scalars ..
166 LOGICAL LOWER, WANTZ
167 INTEGER IINFO, IMAX, INDE, INDWRK, ISCALE
168 DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
169 $ smlnum
170* ..
171* .. External Functions ..
172 LOGICAL LSAME
173 DOUBLE PRECISION DLAMCH, DLANSB
174 EXTERNAL lsame, dlamch, dlansb
175* ..
176* .. External Subroutines ..
177 EXTERNAL dlascl, dsbtrd, dscal, dsteqr, dsterf, xerbla
178* ..
179* .. Intrinsic Functions ..
180 INTRINSIC sqrt
181* ..
182* .. Executable Statements ..
183*
184* Test the input parameters.
185*
186 wantz = lsame( jobz, 'V' )
187 lower = lsame( uplo, 'L' )
188*
189 info = 0
190 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
191 info = -1
192 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
193 info = -2
194 ELSE IF( n.LT.0 ) THEN
195 info = -3
196 ELSE IF( kd.LT.0 ) THEN
197 info = -4
198 ELSE IF( ldab.LT.kd+1 ) THEN
199 info = -6
200 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
201 info = -9
202 END IF
203*
204 IF( info.NE.0 ) THEN
205 CALL xerbla( 'DSBEV ', -info )
206 RETURN
207 END IF
208*
209* Quick return if possible
210*
211 IF( n.EQ.0 )
212 $ RETURN
213*
214 IF( n.EQ.1 ) THEN
215 IF( lower ) THEN
216 w( 1 ) = ab( 1, 1 )
217 ELSE
218 w( 1 ) = ab( kd+1, 1 )
219 END IF
220 IF( wantz )
221 $ z( 1, 1 ) = one
222 RETURN
223 END IF
224*
225* Get machine constants.
226*
227 safmin = dlamch( 'Safe minimum' )
228 eps = dlamch( 'Precision' )
229 smlnum = safmin / eps
230 bignum = one / smlnum
231 rmin = sqrt( smlnum )
232 rmax = sqrt( bignum )
233*
234* Scale matrix to allowable range, if necessary.
235*
236 anrm = dlansb( 'M', uplo, n, kd, ab, ldab, work )
237 iscale = 0
238 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
239 iscale = 1
240 sigma = rmin / anrm
241 ELSE IF( anrm.GT.rmax ) THEN
242 iscale = 1
243 sigma = rmax / anrm
244 END IF
245 IF( iscale.EQ.1 ) THEN
246 IF( lower ) THEN
247 CALL dlascl( 'B', kd, kd, one, sigma, n, n, ab, ldab, info )
248 ELSE
249 CALL dlascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab, info )
250 END IF
251 END IF
252*
253* Call DSBTRD to reduce symmetric band matrix to tridiagonal form.
254*
255 inde = 1
256 indwrk = inde + n
257 CALL dsbtrd( jobz, uplo, n, kd, ab, ldab, w, work( inde ), z, ldz,
258 $ work( indwrk ), iinfo )
259*
260* For eigenvalues only, call DSTERF. For eigenvectors, call SSTEQR.
261*
262 IF( .NOT.wantz ) THEN
263 CALL dsterf( n, w, work( inde ), info )
264 ELSE
265 CALL dsteqr( jobz, n, w, work( inde ), z, ldz, work( indwrk ),
266 $ info )
267 END IF
268*
269* If matrix was scaled, then rescale eigenvalues appropriately.
270*
271 IF( iscale.EQ.1 ) THEN
272 IF( info.EQ.0 ) THEN
273 imax = n
274 ELSE
275 imax = info - 1
276 END IF
277 CALL dscal( imax, one / sigma, w, 1 )
278 END IF
279*
280 RETURN
281*
282* End of DSBEV
283*
284 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsbev(jobz, uplo, n, kd, ab, ldab, w, z, ldz, work, info)
DSBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
Definition dsbev.f:146
subroutine dsbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
DSBTRD
Definition dsbtrd.f:163
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:143
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:131
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86