LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dstegr.f
Go to the documentation of this file.
1*> \brief \b DSTEGR
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DSTEGR + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dstegr.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dstegr.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dstegr.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DSTEGR( JOBZ, RANGE, N, D, E, VL, VU, IL, IU,
22* ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK,
23* LIWORK, INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER JOBZ, RANGE
27* INTEGER IL, INFO, IU, LDZ, LIWORK, LWORK, M, N
28* DOUBLE PRECISION ABSTOL, VL, VU
29* ..
30* .. Array Arguments ..
31* INTEGER ISUPPZ( * ), IWORK( * )
32* DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * )
33* DOUBLE PRECISION Z( LDZ, * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> DSTEGR computes selected eigenvalues and, optionally, eigenvectors
43*> of a real symmetric tridiagonal matrix T. Any such unreduced matrix has
44*> a well defined set of pairwise different real eigenvalues, the corresponding
45*> real eigenvectors are pairwise orthogonal.
46*>
47*> The spectrum may be computed either completely or partially by specifying
48*> either an interval (VL,VU] or a range of indices IL:IU for the desired
49*> eigenvalues.
50*>
51*> DSTEGR is a compatibility wrapper around the improved DSTEMR routine.
52*> See DSTEMR for further details.
53*>
54*> One important change is that the ABSTOL parameter no longer provides any
55*> benefit and hence is no longer used.
56*>
57*> Note : DSTEGR and DSTEMR work only on machines which follow
58*> IEEE-754 floating-point standard in their handling of infinities and
59*> NaNs. Normal execution may create these exceptional values and hence
60*> may abort due to a floating point exception in environments which
61*> do not conform to the IEEE-754 standard.
62*> \endverbatim
63*
64* Arguments:
65* ==========
66*
67*> \param[in] JOBZ
68*> \verbatim
69*> JOBZ is CHARACTER*1
70*> = 'N': Compute eigenvalues only;
71*> = 'V': Compute eigenvalues and eigenvectors.
72*> \endverbatim
73*>
74*> \param[in] RANGE
75*> \verbatim
76*> RANGE is CHARACTER*1
77*> = 'A': all eigenvalues will be found.
78*> = 'V': all eigenvalues in the half-open interval (VL,VU]
79*> will be found.
80*> = 'I': the IL-th through IU-th eigenvalues will be found.
81*> \endverbatim
82*>
83*> \param[in] N
84*> \verbatim
85*> N is INTEGER
86*> The order of the matrix. N >= 0.
87*> \endverbatim
88*>
89*> \param[in,out] D
90*> \verbatim
91*> D is DOUBLE PRECISION array, dimension (N)
92*> On entry, the N diagonal elements of the tridiagonal matrix
93*> T. On exit, D is overwritten.
94*> \endverbatim
95*>
96*> \param[in,out] E
97*> \verbatim
98*> E is DOUBLE PRECISION array, dimension (N)
99*> On entry, the (N-1) subdiagonal elements of the tridiagonal
100*> matrix T in elements 1 to N-1 of E. E(N) need not be set on
101*> input, but is used internally as workspace.
102*> On exit, E is overwritten.
103*> \endverbatim
104*>
105*> \param[in] VL
106*> \verbatim
107*> VL is DOUBLE PRECISION
108*>
109*> If RANGE='V', the lower bound of the interval to
110*> be searched for eigenvalues. VL < VU.
111*> Not referenced if RANGE = 'A' or 'I'.
112*> \endverbatim
113*>
114*> \param[in] VU
115*> \verbatim
116*> VU is DOUBLE PRECISION
117*>
118*> If RANGE='V', the upper bound of the interval to
119*> be searched for eigenvalues. VL < VU.
120*> Not referenced if RANGE = 'A' or 'I'.
121*> \endverbatim
122*>
123*> \param[in] IL
124*> \verbatim
125*> IL is INTEGER
126*>
127*> If RANGE='I', the index of the
128*> smallest eigenvalue to be returned.
129*> 1 <= IL <= IU <= N, if N > 0.
130*> Not referenced if RANGE = 'A' or 'V'.
131*> \endverbatim
132*>
133*> \param[in] IU
134*> \verbatim
135*> IU is INTEGER
136*>
137*> If RANGE='I', the index of the
138*> largest eigenvalue to be returned.
139*> 1 <= IL <= IU <= N, if N > 0.
140*> Not referenced if RANGE = 'A' or 'V'.
141*> \endverbatim
142*>
143*> \param[in] ABSTOL
144*> \verbatim
145*> ABSTOL is DOUBLE PRECISION
146*> Unused. Was the absolute error tolerance for the
147*> eigenvalues/eigenvectors in previous versions.
148*> \endverbatim
149*>
150*> \param[out] M
151*> \verbatim
152*> M is INTEGER
153*> The total number of eigenvalues found. 0 <= M <= N.
154*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
155*> \endverbatim
156*>
157*> \param[out] W
158*> \verbatim
159*> W is DOUBLE PRECISION array, dimension (N)
160*> The first M elements contain the selected eigenvalues in
161*> ascending order.
162*> \endverbatim
163*>
164*> \param[out] Z
165*> \verbatim
166*> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M) )
167*> If JOBZ = 'V', and if INFO = 0, then the first M columns of Z
168*> contain the orthonormal eigenvectors of the matrix T
169*> corresponding to the selected eigenvalues, with the i-th
170*> column of Z holding the eigenvector associated with W(i).
171*> If JOBZ = 'N', then Z is not referenced.
172*> Note: the user must ensure that at least max(1,M) columns are
173*> supplied in the array Z; if RANGE = 'V', the exact value of M
174*> is not known in advance and an upper bound must be used.
175*> Supplying N columns is always safe.
176*> \endverbatim
177*>
178*> \param[in] LDZ
179*> \verbatim
180*> LDZ is INTEGER
181*> The leading dimension of the array Z. LDZ >= 1, and if
182*> JOBZ = 'V', then LDZ >= max(1,N).
183*> \endverbatim
184*>
185*> \param[out] ISUPPZ
186*> \verbatim
187*> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
188*> The support of the eigenvectors in Z, i.e., the indices
189*> indicating the nonzero elements in Z. The i-th computed eigenvector
190*> is nonzero only in elements ISUPPZ( 2*i-1 ) through
191*> ISUPPZ( 2*i ). This is relevant in the case when the matrix
192*> is split. ISUPPZ is only accessed when JOBZ is 'V' and N > 0.
193*> \endverbatim
194*>
195*> \param[out] WORK
196*> \verbatim
197*> WORK is DOUBLE PRECISION array, dimension (LWORK)
198*> On exit, if INFO = 0, WORK(1) returns the optimal
199*> (and minimal) LWORK.
200*> \endverbatim
201*>
202*> \param[in] LWORK
203*> \verbatim
204*> LWORK is INTEGER
205*> The dimension of the array WORK. LWORK >= max(1,18*N)
206*> if JOBZ = 'V', and LWORK >= max(1,12*N) if JOBZ = 'N'.
207*> If LWORK = -1, then a workspace query is assumed; the routine
208*> only calculates the optimal size of the WORK array, returns
209*> this value as the first entry of the WORK array, and no error
210*> message related to LWORK is issued by XERBLA.
211*> \endverbatim
212*>
213*> \param[out] IWORK
214*> \verbatim
215*> IWORK is INTEGER array, dimension (LIWORK)
216*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
217*> \endverbatim
218*>
219*> \param[in] LIWORK
220*> \verbatim
221*> LIWORK is INTEGER
222*> The dimension of the array IWORK. LIWORK >= max(1,10*N)
223*> if the eigenvectors are desired, and LIWORK >= max(1,8*N)
224*> if only the eigenvalues are to be computed.
225*> If LIWORK = -1, then a workspace query is assumed; the
226*> routine only calculates the optimal size of the IWORK array,
227*> returns this value as the first entry of the IWORK array, and
228*> no error message related to LIWORK is issued by XERBLA.
229*> \endverbatim
230*>
231*> \param[out] INFO
232*> \verbatim
233*> INFO is INTEGER
234*> On exit, INFO
235*> = 0: successful exit
236*> < 0: if INFO = -i, the i-th argument had an illegal value
237*> > 0: if INFO = 1X, internal error in DLARRE,
238*> if INFO = 2X, internal error in DLARRV.
239*> Here, the digit X = ABS( IINFO ) < 10, where IINFO is
240*> the nonzero error code returned by DLARRE or
241*> DLARRV, respectively.
242*> \endverbatim
243*
244* Authors:
245* ========
246*
247*> \author Univ. of Tennessee
248*> \author Univ. of California Berkeley
249*> \author Univ. of Colorado Denver
250*> \author NAG Ltd.
251*
252*> \ingroup stegr
253*
254*> \par Contributors:
255* ==================
256*>
257*> Inderjit Dhillon, IBM Almaden, USA \n
258*> Osni Marques, LBNL/NERSC, USA \n
259*> Christof Voemel, LBNL/NERSC, USA \n
260*
261* =====================================================================
262 SUBROUTINE dstegr( JOBZ, RANGE, N, D, E, VL, VU, IL, IU,
263 $ ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK,
264 $ LIWORK, INFO )
265*
266* -- LAPACK computational routine --
267* -- LAPACK is a software package provided by Univ. of Tennessee, --
268* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
269*
270* .. Scalar Arguments ..
271 CHARACTER JOBZ, RANGE
272 INTEGER IL, INFO, IU, LDZ, LIWORK, LWORK, M, N
273 DOUBLE PRECISION ABSTOL, VL, VU
274* ..
275* .. Array Arguments ..
276 INTEGER ISUPPZ( * ), IWORK( * )
277 DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * )
278 DOUBLE PRECISION Z( LDZ, * )
279* ..
280*
281* =====================================================================
282*
283* .. Local Scalars ..
284 LOGICAL TRYRAC
285* ..
286* .. External Subroutines ..
287 EXTERNAL dstemr
288* ..
289* .. Executable Statements ..
290 info = 0
291 tryrac = .false.
292
293 CALL dstemr( jobz, range, n, d, e, vl, vu, il, iu,
294 $ m, w, z, ldz, n, isuppz, tryrac, work, lwork,
295 $ iwork, liwork, info )
296*
297* End of DSTEGR
298*
299 END
subroutine dstegr(jobz, range, n, d, e, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, iwork, liwork, info)
DSTEGR
Definition dstegr.f:265
subroutine dstemr(jobz, range, n, d, e, vl, vu, il, iu, m, w, z, ldz, nzc, isuppz, tryrac, work, lwork, iwork, liwork, info)
DSTEMR
Definition dstemr.f:322