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