LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
sdrvst2stg.f
Go to the documentation of this file.
1 *> \brief \b SDRVST2STG
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE SDRVST2STG( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12 * NOUNIT, A, LDA, D1, D2, D3, D4, EVEIGS, WA1,
13 * WA2, WA3, U, LDU, V, TAU, Z, WORK, LWORK,
14 * IWORK, LIWORK, RESULT, INFO )
15 *
16 * .. Scalar Arguments ..
17 * INTEGER INFO, LDA, LDU, LIWORK, LWORK, NOUNIT, NSIZES,
18 * $ NTYPES
19 * REAL THRESH
20 * ..
21 * .. Array Arguments ..
22 * LOGICAL DOTYPE( * )
23 * INTEGER ISEED( 4 ), IWORK( * ), NN( * )
24 * REAL A( LDA, * ), D1( * ), D2( * ), D3( * ),
25 * $ D4( * ), EVEIGS( * ), RESULT( * ), TAU( * ),
26 * $ U( LDU, * ), V( LDU, * ), WA1( * ), WA2( * ),
27 * $ WA3( * ), WORK( * ), Z( LDU, * )
28 * ..
29 *
30 *
31 *> \par Purpose:
32 * =============
33 *>
34 *> \verbatim
35 *>
36 *> SDRVST2STG checks the symmetric eigenvalue problem drivers.
37 *>
38 *> SSTEV computes all eigenvalues and, optionally,
39 *> eigenvectors of a real symmetric tridiagonal matrix.
40 *>
41 *> SSTEVX computes selected eigenvalues and, optionally,
42 *> eigenvectors of a real symmetric tridiagonal matrix.
43 *>
44 *> SSTEVR computes selected eigenvalues and, optionally,
45 *> eigenvectors of a real symmetric tridiagonal matrix
46 *> using the Relatively Robust Representation where it can.
47 *>
48 *> SSYEV computes all eigenvalues and, optionally,
49 *> eigenvectors of a real symmetric matrix.
50 *>
51 *> SSYEVX computes selected eigenvalues and, optionally,
52 *> eigenvectors of a real symmetric matrix.
53 *>
54 *> SSYEVR computes selected eigenvalues and, optionally,
55 *> eigenvectors of a real symmetric matrix
56 *> using the Relatively Robust Representation where it can.
57 *>
58 *> SSPEV computes all eigenvalues and, optionally,
59 *> eigenvectors of a real symmetric matrix in packed
60 *> storage.
61 *>
62 *> SSPEVX computes selected eigenvalues and, optionally,
63 *> eigenvectors of a real symmetric matrix in packed
64 *> storage.
65 *>
66 *> SSBEV computes all eigenvalues and, optionally,
67 *> eigenvectors of a real symmetric band matrix.
68 *>
69 *> SSBEVX computes selected eigenvalues and, optionally,
70 *> eigenvectors of a real symmetric band matrix.
71 *>
72 *> SSYEVD computes all eigenvalues and, optionally,
73 *> eigenvectors of a real symmetric matrix using
74 *> a divide and conquer algorithm.
75 *>
76 *> SSPEVD computes all eigenvalues and, optionally,
77 *> eigenvectors of a real symmetric matrix in packed
78 *> storage, using a divide and conquer algorithm.
79 *>
80 *> SSBEVD computes all eigenvalues and, optionally,
81 *> eigenvectors of a real symmetric band matrix,
82 *> using a divide and conquer algorithm.
83 *>
84 *> When SDRVST2STG is called, a number of matrix "sizes" ("n's") and a
85 *> number of matrix "types" are specified. For each size ("n")
86 *> and each type of matrix, one matrix will be generated and used
87 *> to test the appropriate drivers. For each matrix and each
88 *> driver routine called, the following tests will be performed:
89 *>
90 *> (1) | A - Z D Z' | / ( |A| n ulp )
91 *>
92 *> (2) | I - Z Z' | / ( n ulp )
93 *>
94 *> (3) | D1 - D2 | / ( |D1| ulp )
95 *>
96 *> where Z is the matrix of eigenvectors returned when the
97 *> eigenvector option is given and D1 and D2 are the eigenvalues
98 *> returned with and without the eigenvector option.
99 *>
100 *> The "sizes" are specified by an array NN(1:NSIZES); the value of
101 *> each element NN(j) specifies one size.
102 *> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
103 *> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
104 *> Currently, the list of possible types is:
105 *>
106 *> (1) The zero matrix.
107 *> (2) The identity matrix.
108 *>
109 *> (3) A diagonal matrix with evenly spaced eigenvalues
110 *> 1, ..., ULP and random signs.
111 *> (ULP = (first number larger than 1) - 1 )
112 *> (4) A diagonal matrix with geometrically spaced eigenvalues
113 *> 1, ..., ULP and random signs.
114 *> (5) A diagonal matrix with "clustered" eigenvalues
115 *> 1, ULP, ..., ULP and random signs.
116 *>
117 *> (6) Same as (4), but multiplied by SQRT( overflow threshold )
118 *> (7) Same as (4), but multiplied by SQRT( underflow threshold )
119 *>
120 *> (8) A matrix of the form U' D U, where U is orthogonal and
121 *> D has evenly spaced entries 1, ..., ULP with random signs
122 *> on the diagonal.
123 *>
124 *> (9) A matrix of the form U' D U, where U is orthogonal and
125 *> D has geometrically spaced entries 1, ..., ULP with random
126 *> signs on the diagonal.
127 *>
128 *> (10) A matrix of the form U' D U, where U is orthogonal and
129 *> D has "clustered" entries 1, ULP,..., ULP with random
130 *> signs on the diagonal.
131 *>
132 *> (11) Same as (8), but multiplied by SQRT( overflow threshold )
133 *> (12) Same as (8), but multiplied by SQRT( underflow threshold )
134 *>
135 *> (13) Symmetric matrix with random entries chosen from (-1,1).
136 *> (14) Same as (13), but multiplied by SQRT( overflow threshold )
137 *> (15) Same as (13), but multiplied by SQRT( underflow threshold )
138 *> (16) A band matrix with half bandwidth randomly chosen between
139 *> 0 and N-1, with evenly spaced eigenvalues 1, ..., ULP
140 *> with random signs.
141 *> (17) Same as (16), but multiplied by SQRT( overflow threshold )
142 *> (18) Same as (16), but multiplied by SQRT( underflow threshold )
143 *> \endverbatim
144 *
145 * Arguments:
146 * ==========
147 *
148 *> \verbatim
149 *> NSIZES INTEGER
150 *> The number of sizes of matrices to use. If it is zero,
151 *> SDRVST2STG does nothing. It must be at least zero.
152 *> Not modified.
153 *>
154 *> NN INTEGER array, dimension (NSIZES)
155 *> An array containing the sizes to be used for the matrices.
156 *> Zero values will be skipped. The values must be at least
157 *> zero.
158 *> Not modified.
159 *>
160 *> NTYPES INTEGER
161 *> The number of elements in DOTYPE. If it is zero, SDRVST2STG
162 *> does nothing. It must be at least zero. If it is MAXTYP+1
163 *> and NSIZES is 1, then an additional type, MAXTYP+1 is
164 *> defined, which is to use whatever matrix is in A. This
165 *> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
166 *> DOTYPE(MAXTYP+1) is .TRUE. .
167 *> Not modified.
168 *>
169 *> DOTYPE LOGICAL array, dimension (NTYPES)
170 *> If DOTYPE(j) is .TRUE., then for each size in NN a
171 *> matrix of that size and of type j will be generated.
172 *> If NTYPES is smaller than the maximum number of types
173 *> defined (PARAMETER MAXTYP), then types NTYPES+1 through
174 *> MAXTYP will not be generated. If NTYPES is larger
175 *> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
176 *> will be ignored.
177 *> Not modified.
178 *>
179 *> ISEED INTEGER array, dimension (4)
180 *> On entry ISEED specifies the seed of the random number
181 *> generator. The array elements should be between 0 and 4095;
182 *> if not they will be reduced mod 4096. Also, ISEED(4) must
183 *> be odd. The random number generator uses a linear
184 *> congruential sequence limited to small integers, and so
185 *> should produce machine independent random numbers. The
186 *> values of ISEED are changed on exit, and can be used in the
187 *> next call to SDRVST2STG to continue the same random number
188 *> sequence.
189 *> Modified.
190 *>
191 *> THRESH REAL
192 *> A test will count as "failed" if the "error", computed as
193 *> described above, exceeds THRESH. Note that the error
194 *> is scaled to be O(1), so THRESH should be a reasonably
195 *> small multiple of 1, e.g., 10 or 100. In particular,
196 *> it should not depend on the precision (single vs. double)
197 *> or the size of the matrix. It must be at least zero.
198 *> Not modified.
199 *>
200 *> NOUNIT INTEGER
201 *> The FORTRAN unit number for printing out error messages
202 *> (e.g., if a routine returns IINFO not equal to 0.)
203 *> Not modified.
204 *>
205 *> A REAL array, dimension (LDA , max(NN))
206 *> Used to hold the matrix whose eigenvalues are to be
207 *> computed. On exit, A contains the last matrix actually
208 *> used.
209 *> Modified.
210 *>
211 *> LDA INTEGER
212 *> The leading dimension of A. It must be at
213 *> least 1 and at least max( NN ).
214 *> Not modified.
215 *>
216 *> D1 REAL array, dimension (max(NN))
217 *> The eigenvalues of A, as computed by SSTEQR simlutaneously
218 *> with Z. On exit, the eigenvalues in D1 correspond with the
219 *> matrix in A.
220 *> Modified.
221 *>
222 *> D2 REAL array, dimension (max(NN))
223 *> The eigenvalues of A, as computed by SSTEQR if Z is not
224 *> computed. On exit, the eigenvalues in D2 correspond with
225 *> the matrix in A.
226 *> Modified.
227 *>
228 *> D3 REAL array, dimension (max(NN))
229 *> The eigenvalues of A, as computed by SSTERF. On exit, the
230 *> eigenvalues in D3 correspond with the matrix in A.
231 *> Modified.
232 *>
233 *> D4 REAL array, dimension
234 *>
235 *> EVEIGS REAL array, dimension (max(NN))
236 *> The eigenvalues as computed by SSTEV('N', ... )
237 *> (I reserve the right to change this to the output of
238 *> whichever algorithm computes the most accurate eigenvalues).
239 *>
240 *> WA1 REAL array, dimension
241 *>
242 *> WA2 REAL array, dimension
243 *>
244 *> WA3 REAL array, dimension
245 *>
246 *> U REAL array, dimension (LDU, max(NN))
247 *> The orthogonal matrix computed by SSYTRD + SORGTR.
248 *> Modified.
249 *>
250 *> LDU INTEGER
251 *> The leading dimension of U, Z, and V. It must be at
252 *> least 1 and at least max( NN ).
253 *> Not modified.
254 *>
255 *> V REAL array, dimension (LDU, max(NN))
256 *> The Housholder vectors computed by SSYTRD in reducing A to
257 *> tridiagonal form.
258 *> Modified.
259 *>
260 *> TAU REAL array, dimension (max(NN))
261 *> The Householder factors computed by SSYTRD in reducing A
262 *> to tridiagonal form.
263 *> Modified.
264 *>
265 *> Z REAL array, dimension (LDU, max(NN))
266 *> The orthogonal matrix of eigenvectors computed by SSTEQR,
267 *> SPTEQR, and SSTEIN.
268 *> Modified.
269 *>
270 *> WORK REAL array, dimension (LWORK)
271 *> Workspace.
272 *> Modified.
273 *>
274 *> LWORK INTEGER
275 *> The number of entries in WORK. This must be at least
276 *> 1 + 4 * Nmax + 2 * Nmax * lg Nmax + 4 * Nmax**2
277 *> where Nmax = max( NN(j), 2 ) and lg = log base 2.
278 *> Not modified.
279 *>
280 *> IWORK INTEGER array,
281 *> dimension (6 + 6*Nmax + 5 * Nmax * lg Nmax )
282 *> where Nmax = max( NN(j), 2 ) and lg = log base 2.
283 *> Workspace.
284 *> Modified.
285 *>
286 *> RESULT REAL array, dimension (105)
287 *> The values computed by the tests described above.
288 *> The values are currently limited to 1/ulp, to avoid
289 *> overflow.
290 *> Modified.
291 *>
292 *> INFO INTEGER
293 *> If 0, then everything ran OK.
294 *> -1: NSIZES < 0
295 *> -2: Some NN(j) < 0
296 *> -3: NTYPES < 0
297 *> -5: THRESH < 0
298 *> -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
299 *> -16: LDU < 1 or LDU < NMAX.
300 *> -21: LWORK too small.
301 *> If SLATMR, SLATMS, SSYTRD, SORGTR, SSTEQR, SSTERF,
302 *> or SORMTR returns an error code, the
303 *> absolute value of it is returned.
304 *> Modified.
305 *>
306 *>-----------------------------------------------------------------------
307 *>
308 *> Some Local Variables and Parameters:
309 *> ---- ----- --------- --- ----------
310 *> ZERO, ONE Real 0 and 1.
311 *> MAXTYP The number of types defined.
312 *> NTEST The number of tests performed, or which can
313 *> be performed so far, for the current matrix.
314 *> NTESTT The total number of tests performed so far.
315 *> NMAX Largest value in NN.
316 *> NMATS The number of matrices generated so far.
317 *> NERRS The number of tests which have exceeded THRESH
318 *> so far (computed by SLAFTS).
319 *> COND, IMODE Values to be passed to the matrix generators.
320 *> ANORM Norm of A; passed to matrix generators.
321 *>
322 *> OVFL, UNFL Overflow and underflow thresholds.
323 *> ULP, ULPINV Finest relative precision and its inverse.
324 *> RTOVFL, RTUNFL Square roots of the previous 2 values.
325 *> The following four arrays decode JTYPE:
326 *> KTYPE(j) The general type (1-10) for type "j".
327 *> KMODE(j) The MODE value to be passed to the matrix
328 *> generator for type "j".
329 *> KMAGN(j) The order of magnitude ( O(1),
330 *> O(overflow^(1/2) ), O(underflow^(1/2) )
331 *>
332 *> The tests performed are: Routine tested
333 *> 1= | A - U S U' | / ( |A| n ulp ) SSTEV('V', ... )
334 *> 2= | I - U U' | / ( n ulp ) SSTEV('V', ... )
335 *> 3= |D(with Z) - D(w/o Z)| / (|D| ulp) SSTEV('N', ... )
336 *> 4= | A - U S U' | / ( |A| n ulp ) SSTEVX('V','A', ... )
337 *> 5= | I - U U' | / ( n ulp ) SSTEVX('V','A', ... )
338 *> 6= |D(with Z) - EVEIGS| / (|D| ulp) SSTEVX('N','A', ... )
339 *> 7= | A - U S U' | / ( |A| n ulp ) SSTEVR('V','A', ... )
340 *> 8= | I - U U' | / ( n ulp ) SSTEVR('V','A', ... )
341 *> 9= |D(with Z) - EVEIGS| / (|D| ulp) SSTEVR('N','A', ... )
342 *> 10= | A - U S U' | / ( |A| n ulp ) SSTEVX('V','I', ... )
343 *> 11= | I - U U' | / ( n ulp ) SSTEVX('V','I', ... )
344 *> 12= |D(with Z) - D(w/o Z)| / (|D| ulp) SSTEVX('N','I', ... )
345 *> 13= | A - U S U' | / ( |A| n ulp ) SSTEVX('V','V', ... )
346 *> 14= | I - U U' | / ( n ulp ) SSTEVX('V','V', ... )
347 *> 15= |D(with Z) - D(w/o Z)| / (|D| ulp) SSTEVX('N','V', ... )
348 *> 16= | A - U S U' | / ( |A| n ulp ) SSTEVD('V', ... )
349 *> 17= | I - U U' | / ( n ulp ) SSTEVD('V', ... )
350 *> 18= |D(with Z) - EVEIGS| / (|D| ulp) SSTEVD('N', ... )
351 *> 19= | A - U S U' | / ( |A| n ulp ) SSTEVR('V','I', ... )
352 *> 20= | I - U U' | / ( n ulp ) SSTEVR('V','I', ... )
353 *> 21= |D(with Z) - D(w/o Z)| / (|D| ulp) SSTEVR('N','I', ... )
354 *> 22= | A - U S U' | / ( |A| n ulp ) SSTEVR('V','V', ... )
355 *> 23= | I - U U' | / ( n ulp ) SSTEVR('V','V', ... )
356 *> 24= |D(with Z) - D(w/o Z)| / (|D| ulp) SSTEVR('N','V', ... )
357 *>
358 *> 25= | A - U S U' | / ( |A| n ulp ) SSYEV('L','V', ... )
359 *> 26= | I - U U' | / ( n ulp ) SSYEV('L','V', ... )
360 *> 27= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEV_2STAGE('L','N', ... )
361 *> 28= | A - U S U' | / ( |A| n ulp ) SSYEVX('L','V','A', ... )
362 *> 29= | I - U U' | / ( n ulp ) SSYEVX('L','V','A', ... )
363 *> 30= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVX_2STAGE('L','N','A', ... )
364 *> 31= | A - U S U' | / ( |A| n ulp ) SSYEVX('L','V','I', ... )
365 *> 32= | I - U U' | / ( n ulp ) SSYEVX('L','V','I', ... )
366 *> 33= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVX_2STAGE('L','N','I', ... )
367 *> 34= | A - U S U' | / ( |A| n ulp ) SSYEVX('L','V','V', ... )
368 *> 35= | I - U U' | / ( n ulp ) SSYEVX('L','V','V', ... )
369 *> 36= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVX_2STAGE('L','N','V', ... )
370 *> 37= | A - U S U' | / ( |A| n ulp ) SSPEV('L','V', ... )
371 *> 38= | I - U U' | / ( n ulp ) SSPEV('L','V', ... )
372 *> 39= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEV('L','N', ... )
373 *> 40= | A - U S U' | / ( |A| n ulp ) SSPEVX('L','V','A', ... )
374 *> 41= | I - U U' | / ( n ulp ) SSPEVX('L','V','A', ... )
375 *> 42= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVX('L','N','A', ... )
376 *> 43= | A - U S U' | / ( |A| n ulp ) SSPEVX('L','V','I', ... )
377 *> 44= | I - U U' | / ( n ulp ) SSPEVX('L','V','I', ... )
378 *> 45= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVX('L','N','I', ... )
379 *> 46= | A - U S U' | / ( |A| n ulp ) SSPEVX('L','V','V', ... )
380 *> 47= | I - U U' | / ( n ulp ) SSPEVX('L','V','V', ... )
381 *> 48= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVX('L','N','V', ... )
382 *> 49= | A - U S U' | / ( |A| n ulp ) SSBEV('L','V', ... )
383 *> 50= | I - U U' | / ( n ulp ) SSBEV('L','V', ... )
384 *> 51= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEV_2STAGE('L','N', ... )
385 *> 52= | A - U S U' | / ( |A| n ulp ) SSBEVX('L','V','A', ... )
386 *> 53= | I - U U' | / ( n ulp ) SSBEVX('L','V','A', ... )
387 *> 54= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVX_2STAGE('L','N','A', ... )
388 *> 55= | A - U S U' | / ( |A| n ulp ) SSBEVX('L','V','I', ... )
389 *> 56= | I - U U' | / ( n ulp ) SSBEVX('L','V','I', ... )
390 *> 57= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVX_2STAGE('L','N','I', ... )
391 *> 58= | A - U S U' | / ( |A| n ulp ) SSBEVX('L','V','V', ... )
392 *> 59= | I - U U' | / ( n ulp ) SSBEVX('L','V','V', ... )
393 *> 60= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVX_2STAGE('L','N','V', ... )
394 *> 61= | A - U S U' | / ( |A| n ulp ) SSYEVD('L','V', ... )
395 *> 62= | I - U U' | / ( n ulp ) SSYEVD('L','V', ... )
396 *> 63= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVD_2STAGE('L','N', ... )
397 *> 64= | A - U S U' | / ( |A| n ulp ) SSPEVD('L','V', ... )
398 *> 65= | I - U U' | / ( n ulp ) SSPEVD('L','V', ... )
399 *> 66= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVD('L','N', ... )
400 *> 67= | A - U S U' | / ( |A| n ulp ) SSBEVD('L','V', ... )
401 *> 68= | I - U U' | / ( n ulp ) SSBEVD('L','V', ... )
402 *> 69= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVD_2STAGE('L','N', ... )
403 *> 70= | A - U S U' | / ( |A| n ulp ) SSYEVR('L','V','A', ... )
404 *> 71= | I - U U' | / ( n ulp ) SSYEVR('L','V','A', ... )
405 *> 72= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVR_2STAGE('L','N','A', ... )
406 *> 73= | A - U S U' | / ( |A| n ulp ) SSYEVR('L','V','I', ... )
407 *> 74= | I - U U' | / ( n ulp ) SSYEVR('L','V','I', ... )
408 *> 75= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVR_2STAGE('L','N','I', ... )
409 *> 76= | A - U S U' | / ( |A| n ulp ) SSYEVR('L','V','V', ... )
410 *> 77= | I - U U' | / ( n ulp ) SSYEVR('L','V','V', ... )
411 *> 78= |D(with Z) - D(w/o Z)| / (|D| ulp) SSYEVR_2STAGE('L','N','V', ... )
412 *>
413 *> Tests 25 through 78 are repeated (as tests 79 through 132)
414 *> with UPLO='U'
415 *>
416 *> To be added in 1999
417 *>
418 *> 79= | A - U S U' | / ( |A| n ulp ) SSPEVR('L','V','A', ... )
419 *> 80= | I - U U' | / ( n ulp ) SSPEVR('L','V','A', ... )
420 *> 81= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVR('L','N','A', ... )
421 *> 82= | A - U S U' | / ( |A| n ulp ) SSPEVR('L','V','I', ... )
422 *> 83= | I - U U' | / ( n ulp ) SSPEVR('L','V','I', ... )
423 *> 84= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVR('L','N','I', ... )
424 *> 85= | A - U S U' | / ( |A| n ulp ) SSPEVR('L','V','V', ... )
425 *> 86= | I - U U' | / ( n ulp ) SSPEVR('L','V','V', ... )
426 *> 87= |D(with Z) - D(w/o Z)| / (|D| ulp) SSPEVR('L','N','V', ... )
427 *> 88= | A - U S U' | / ( |A| n ulp ) SSBEVR('L','V','A', ... )
428 *> 89= | I - U U' | / ( n ulp ) SSBEVR('L','V','A', ... )
429 *> 90= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVR('L','N','A', ... )
430 *> 91= | A - U S U' | / ( |A| n ulp ) SSBEVR('L','V','I', ... )
431 *> 92= | I - U U' | / ( n ulp ) SSBEVR('L','V','I', ... )
432 *> 93= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVR('L','N','I', ... )
433 *> 94= | A - U S U' | / ( |A| n ulp ) SSBEVR('L','V','V', ... )
434 *> 95= | I - U U' | / ( n ulp ) SSBEVR('L','V','V', ... )
435 *> 96= |D(with Z) - D(w/o Z)| / (|D| ulp) SSBEVR('L','N','V', ... )
436 *> \endverbatim
437 *
438 * Authors:
439 * ========
440 *
441 *> \author Univ. of Tennessee
442 *> \author Univ. of California Berkeley
443 *> \author Univ. of Colorado Denver
444 *> \author NAG Ltd.
445 *
446 *> \date December 2016
447 *
448 *> \ingroup single_eig
449 *
450 * =====================================================================
451  SUBROUTINE sdrvst2stg( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
452  $ nounit, a, lda, d1, d2, d3, d4, eveigs, wa1,
453  $ wa2, wa3, u, ldu, v, tau, z, work, lwork,
454  $ iwork, liwork, result, info )
455 *
456 * -- LAPACK test routine (version 3.7.0) --
457 * -- LAPACK is a software package provided by Univ. of Tennessee, --
458 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
459 * December 2016
460 *
461 * .. Scalar Arguments ..
462  INTEGER INFO, LDA, LDU, LIWORK, LWORK, NOUNIT, NSIZES,
463  $ ntypes
464  REAL THRESH
465 * ..
466 * .. Array Arguments ..
467  LOGICAL DOTYPE( * )
468  INTEGER ISEED( 4 ), IWORK( * ), NN( * )
469  REAL A( lda, * ), D1( * ), D2( * ), D3( * ),
470  $ d4( * ), eveigs( * ), result( * ), tau( * ),
471  $ u( ldu, * ), v( ldu, * ), wa1( * ), wa2( * ),
472  $ wa3( * ), work( * ), z( ldu, * )
473 * ..
474 *
475 * =====================================================================
476 *
477 * .. Parameters ..
478  REAL ZERO, ONE, TWO, TEN
479  parameter ( zero = 0.0e0, one = 1.0e0, two = 2.0e0,
480  $ ten = 10.0e0 )
481  REAL HALF
482  parameter ( half = 0.5e+0 )
483  INTEGER MAXTYP
484  parameter ( maxtyp = 18 )
485 * ..
486 * .. Local Scalars ..
487  LOGICAL BADNN
488  CHARACTER UPLO
489  INTEGER I, IDIAG, IHBW, IINFO, IL, IMODE, INDX, IROW,
490  $ itemp, itype, iu, iuplo, j, j1, j2, jcol,
491  $ jsize, jtype, kd, lgn, liwedc, lwedc, m, m2,
492  $ m3, mtypes, n, nerrs, nmats, nmax, ntest,
493  $ ntestt
494  REAL ABSTOL, ANINV, ANORM, COND, OVFL, RTOVFL,
495  $ rtunfl, temp1, temp2, temp3, ulp, ulpinv, unfl,
496  $ vl, vu
497 * ..
498 * .. Local Arrays ..
499  INTEGER IDUMMA( 1 ), IOLDSD( 4 ), ISEED2( 4 ),
500  $ iseed3( 4 ), kmagn( maxtyp ), kmode( maxtyp ),
501  $ ktype( maxtyp )
502 * ..
503 * .. External Functions ..
504  REAL SLAMCH, SLARND, SSXT1
505  EXTERNAL slamch, slarnd, ssxt1
506 * ..
507 * .. External Subroutines ..
508  EXTERNAL alasvm, slabad, slacpy, slafts, slaset, slatmr,
515  $ ssytrd_sb2st, ssyt22, xerbla
516 * ..
517 * .. Scalars in Common ..
518  CHARACTER*32 SRNAMT
519 * ..
520 * .. Common blocks ..
521  COMMON / srnamc / srnamt
522 * ..
523 * .. Intrinsic Functions ..
524  INTRINSIC abs, REAL, INT, LOG, MAX, MIN, SQRT
525 * ..
526 * .. Data statements ..
527  DATA ktype / 1, 2, 5*4, 5*5, 3*8, 3*9 /
528  DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
529  $ 2, 3, 1, 2, 3 /
530  DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
531  $ 0, 0, 4, 4, 4 /
532 * ..
533 * .. Executable Statements ..
534 *
535 * Keep ftrnchek happy
536 *
537  vl = zero
538  vu = zero
539 *
540 * 1) Check for errors
541 *
542  ntestt = 0
543  info = 0
544 *
545  badnn = .false.
546  nmax = 1
547  DO 10 j = 1, nsizes
548  nmax = max( nmax, nn( j ) )
549  IF( nn( j ).LT.0 )
550  $ badnn = .true.
551  10 CONTINUE
552 *
553 * Check for errors
554 *
555  IF( nsizes.LT.0 ) THEN
556  info = -1
557  ELSE IF( badnn ) THEN
558  info = -2
559  ELSE IF( ntypes.LT.0 ) THEN
560  info = -3
561  ELSE IF( lda.LT.nmax ) THEN
562  info = -9
563  ELSE IF( ldu.LT.nmax ) THEN
564  info = -16
565  ELSE IF( 2*max( 2, nmax )**2.GT.lwork ) THEN
566  info = -21
567  END IF
568 *
569  IF( info.NE.0 ) THEN
570  CALL xerbla( 'SDRVST2STG', -info )
571  RETURN
572  END IF
573 *
574 * Quick return if nothing to do
575 *
576  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
577  $ RETURN
578 *
579 * More Important constants
580 *
581  unfl = slamch( 'Safe minimum' )
582  ovfl = slamch( 'Overflow' )
583  CALL slabad( unfl, ovfl )
584  ulp = slamch( 'Epsilon' )*slamch( 'Base' )
585  ulpinv = one / ulp
586  rtunfl = sqrt( unfl )
587  rtovfl = sqrt( ovfl )
588 *
589 * Loop over sizes, types
590 *
591  DO 20 i = 1, 4
592  iseed2( i ) = iseed( i )
593  iseed3( i ) = iseed( i )
594  20 CONTINUE
595 *
596  nerrs = 0
597  nmats = 0
598 *
599 *
600  DO 1740 jsize = 1, nsizes
601  n = nn( jsize )
602  IF( n.GT.0 ) THEN
603  lgn = int( log( REAL( N ) ) / log( TWO ) )
604  IF( 2**lgn.LT.n )
605  $ lgn = lgn + 1
606  IF( 2**lgn.LT.n )
607  $ lgn = lgn + 1
608  lwedc = 1 + 4*n + 2*n*lgn + 4*n**2
609 c LIWEDC = 6 + 6*N + 5*N*LGN
610  liwedc = 3 + 5*n
611  ELSE
612  lwedc = 9
613 c LIWEDC = 12
614  liwedc = 8
615  END IF
616  aninv = one / REAL( MAX( 1, N ) )
617 *
618  IF( nsizes.NE.1 ) THEN
619  mtypes = min( maxtyp, ntypes )
620  ELSE
621  mtypes = min( maxtyp+1, ntypes )
622  END IF
623 *
624  DO 1730 jtype = 1, mtypes
625 *
626  IF( .NOT.dotype( jtype ) )
627  $ GO TO 1730
628  nmats = nmats + 1
629  ntest = 0
630 *
631  DO 30 j = 1, 4
632  ioldsd( j ) = iseed( j )
633  30 CONTINUE
634 *
635 * 2) Compute "A"
636 *
637 * Control parameters:
638 *
639 * KMAGN KMODE KTYPE
640 * =1 O(1) clustered 1 zero
641 * =2 large clustered 2 identity
642 * =3 small exponential (none)
643 * =4 arithmetic diagonal, (w/ eigenvalues)
644 * =5 random log symmetric, w/ eigenvalues
645 * =6 random (none)
646 * =7 random diagonal
647 * =8 random symmetric
648 * =9 band symmetric, w/ eigenvalues
649 *
650  IF( mtypes.GT.maxtyp )
651  $ GO TO 110
652 *
653  itype = ktype( jtype )
654  imode = kmode( jtype )
655 *
656 * Compute norm
657 *
658  GO TO ( 40, 50, 60 )kmagn( jtype )
659 *
660  40 CONTINUE
661  anorm = one
662  GO TO 70
663 *
664  50 CONTINUE
665  anorm = ( rtovfl*ulp )*aninv
666  GO TO 70
667 *
668  60 CONTINUE
669  anorm = rtunfl*n*ulpinv
670  GO TO 70
671 *
672  70 CONTINUE
673 *
674  CALL slaset( 'Full', lda, n, zero, zero, a, lda )
675  iinfo = 0
676  cond = ulpinv
677 *
678 * Special Matrices -- Identity & Jordan block
679 *
680 * Zero
681 *
682  IF( itype.EQ.1 ) THEN
683  iinfo = 0
684 *
685  ELSE IF( itype.EQ.2 ) THEN
686 *
687 * Identity
688 *
689  DO 80 jcol = 1, n
690  a( jcol, jcol ) = anorm
691  80 CONTINUE
692 *
693  ELSE IF( itype.EQ.4 ) THEN
694 *
695 * Diagonal Matrix, [Eigen]values Specified
696 *
697  CALL slatms( n, n, 'S', iseed, 'S', work, imode, cond,
698  $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
699  $ iinfo )
700 *
701  ELSE IF( itype.EQ.5 ) THEN
702 *
703 * Symmetric, eigenvalues specified
704 *
705  CALL slatms( n, n, 'S', iseed, 'S', work, imode, cond,
706  $ anorm, n, n, 'N', a, lda, work( n+1 ),
707  $ iinfo )
708 *
709  ELSE IF( itype.EQ.7 ) THEN
710 *
711 * Diagonal, random eigenvalues
712 *
713  idumma( 1 ) = 1
714  CALL slatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
715  $ 'T', 'N', work( n+1 ), 1, one,
716  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
717  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
718 *
719  ELSE IF( itype.EQ.8 ) THEN
720 *
721 * Symmetric, random eigenvalues
722 *
723  idumma( 1 ) = 1
724  CALL slatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
725  $ 'T', 'N', work( n+1 ), 1, one,
726  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
727  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
728 *
729  ELSE IF( itype.EQ.9 ) THEN
730 *
731 * Symmetric banded, eigenvalues specified
732 *
733  ihbw = int( ( n-1 )*slarnd( 1, iseed3 ) )
734  CALL slatms( n, n, 'S', iseed, 'S', work, imode, cond,
735  $ anorm, ihbw, ihbw, 'Z', u, ldu, work( n+1 ),
736  $ iinfo )
737 *
738 * Store as dense matrix for most routines.
739 *
740  CALL slaset( 'Full', lda, n, zero, zero, a, lda )
741  DO 100 idiag = -ihbw, ihbw
742  irow = ihbw - idiag + 1
743  j1 = max( 1, idiag+1 )
744  j2 = min( n, n+idiag )
745  DO 90 j = j1, j2
746  i = j - idiag
747  a( i, j ) = u( irow, j )
748  90 CONTINUE
749  100 CONTINUE
750  ELSE
751  iinfo = 1
752  END IF
753 *
754  IF( iinfo.NE.0 ) THEN
755  WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
756  $ ioldsd
757  info = abs( iinfo )
758  RETURN
759  END IF
760 *
761  110 CONTINUE
762 *
763  abstol = unfl + unfl
764  IF( n.LE.1 ) THEN
765  il = 1
766  iu = n
767  ELSE
768  il = 1 + int( ( n-1 )*slarnd( 1, iseed2 ) )
769  iu = 1 + int( ( n-1 )*slarnd( 1, iseed2 ) )
770  IF( il.GT.iu ) THEN
771  itemp = il
772  il = iu
773  iu = itemp
774  END IF
775  END IF
776 *
777 * 3) If matrix is tridiagonal, call SSTEV and SSTEVX.
778 *
779  IF( jtype.LE.7 ) THEN
780  ntest = 1
781  DO 120 i = 1, n
782  d1( i ) = REAL( A( I, I ) )
783  120 CONTINUE
784  DO 130 i = 1, n - 1
785  d2( i ) = REAL( A( I+1, I ) )
786  130 CONTINUE
787  srnamt = 'SSTEV'
788  CALL sstev( 'V', n, d1, d2, z, ldu, work, iinfo )
789  IF( iinfo.NE.0 ) THEN
790  WRITE( nounit, fmt = 9999 )'SSTEV(V)', iinfo, n,
791  $ jtype, ioldsd
792  info = abs( iinfo )
793  IF( iinfo.LT.0 ) THEN
794  RETURN
795  ELSE
796  result( 1 ) = ulpinv
797  result( 2 ) = ulpinv
798  result( 3 ) = ulpinv
799  GO TO 180
800  END IF
801  END IF
802 *
803 * Do tests 1 and 2.
804 *
805  DO 140 i = 1, n
806  d3( i ) = REAL( A( I, I ) )
807  140 CONTINUE
808  DO 150 i = 1, n - 1
809  d4( i ) = REAL( A( I+1, I ) )
810  150 CONTINUE
811  CALL sstt21( n, 0, d3, d4, d1, d2, z, ldu, work,
812  $ result( 1 ) )
813 *
814  ntest = 3
815  DO 160 i = 1, n - 1
816  d4( i ) = REAL( A( I+1, I ) )
817  160 CONTINUE
818  srnamt = 'SSTEV'
819  CALL sstev( 'N', n, d3, d4, z, ldu, work, iinfo )
820  IF( iinfo.NE.0 ) THEN
821  WRITE( nounit, fmt = 9999 )'SSTEV(N)', iinfo, n,
822  $ jtype, ioldsd
823  info = abs( iinfo )
824  IF( iinfo.LT.0 ) THEN
825  RETURN
826  ELSE
827  result( 3 ) = ulpinv
828  GO TO 180
829  END IF
830  END IF
831 *
832 * Do test 3.
833 *
834  temp1 = zero
835  temp2 = zero
836  DO 170 j = 1, n
837  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
838  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
839  170 CONTINUE
840  result( 3 ) = temp2 / max( unfl,
841  $ ulp*max( temp1, temp2 ) )
842 *
843  180 CONTINUE
844 *
845  ntest = 4
846  DO 190 i = 1, n
847  eveigs( i ) = d3( i )
848  d1( i ) = REAL( A( I, I ) )
849  190 CONTINUE
850  DO 200 i = 1, n - 1
851  d2( i ) = REAL( A( I+1, I ) )
852  200 CONTINUE
853  srnamt = 'SSTEVX'
854  CALL sstevx( 'V', 'A', n, d1, d2, vl, vu, il, iu, abstol,
855  $ m, wa1, z, ldu, work, iwork, iwork( 5*n+1 ),
856  $ iinfo )
857  IF( iinfo.NE.0 ) THEN
858  WRITE( nounit, fmt = 9999 )'SSTEVX(V,A)', iinfo, n,
859  $ jtype, ioldsd
860  info = abs( iinfo )
861  IF( iinfo.LT.0 ) THEN
862  RETURN
863  ELSE
864  result( 4 ) = ulpinv
865  result( 5 ) = ulpinv
866  result( 6 ) = ulpinv
867  GO TO 250
868  END IF
869  END IF
870  IF( n.GT.0 ) THEN
871  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
872  ELSE
873  temp3 = zero
874  END IF
875 *
876 * Do tests 4 and 5.
877 *
878  DO 210 i = 1, n
879  d3( i ) = REAL( A( I, I ) )
880  210 CONTINUE
881  DO 220 i = 1, n - 1
882  d4( i ) = REAL( A( I+1, I ) )
883  220 CONTINUE
884  CALL sstt21( n, 0, d3, d4, wa1, d2, z, ldu, work,
885  $ result( 4 ) )
886 *
887  ntest = 6
888  DO 230 i = 1, n - 1
889  d4( i ) = REAL( A( I+1, I ) )
890  230 CONTINUE
891  srnamt = 'SSTEVX'
892  CALL sstevx( 'N', 'A', n, d3, d4, vl, vu, il, iu, abstol,
893  $ m2, wa2, z, ldu, work, iwork,
894  $ iwork( 5*n+1 ), iinfo )
895  IF( iinfo.NE.0 ) THEN
896  WRITE( nounit, fmt = 9999 )'SSTEVX(N,A)', iinfo, n,
897  $ jtype, ioldsd
898  info = abs( iinfo )
899  IF( iinfo.LT.0 ) THEN
900  RETURN
901  ELSE
902  result( 6 ) = ulpinv
903  GO TO 250
904  END IF
905  END IF
906 *
907 * Do test 6.
908 *
909  temp1 = zero
910  temp2 = zero
911  DO 240 j = 1, n
912  temp1 = max( temp1, abs( wa2( j ) ),
913  $ abs( eveigs( j ) ) )
914  temp2 = max( temp2, abs( wa2( j )-eveigs( j ) ) )
915  240 CONTINUE
916  result( 6 ) = temp2 / max( unfl,
917  $ ulp*max( temp1, temp2 ) )
918 *
919  250 CONTINUE
920 *
921  ntest = 7
922  DO 260 i = 1, n
923  d1( i ) = REAL( A( I, I ) )
924  260 CONTINUE
925  DO 270 i = 1, n - 1
926  d2( i ) = REAL( A( I+1, I ) )
927  270 CONTINUE
928  srnamt = 'SSTEVR'
929  CALL sstevr( 'V', 'A', n, d1, d2, vl, vu, il, iu, abstol,
930  $ m, wa1, z, ldu, iwork, work, lwork,
931  $ iwork(2*n+1), liwork-2*n, iinfo )
932  IF( iinfo.NE.0 ) THEN
933  WRITE( nounit, fmt = 9999 )'SSTEVR(V,A)', iinfo, n,
934  $ jtype, ioldsd
935  info = abs( iinfo )
936  IF( iinfo.LT.0 ) THEN
937  RETURN
938  ELSE
939  result( 7 ) = ulpinv
940  result( 8 ) = ulpinv
941  GO TO 320
942  END IF
943  END IF
944  IF( n.GT.0 ) THEN
945  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
946  ELSE
947  temp3 = zero
948  END IF
949 *
950 * Do tests 7 and 8.
951 *
952  DO 280 i = 1, n
953  d3( i ) = REAL( A( I, I ) )
954  280 CONTINUE
955  DO 290 i = 1, n - 1
956  d4( i ) = REAL( A( I+1, I ) )
957  290 CONTINUE
958  CALL sstt21( n, 0, d3, d4, wa1, d2, z, ldu, work,
959  $ result( 7 ) )
960 *
961  ntest = 9
962  DO 300 i = 1, n - 1
963  d4( i ) = REAL( A( I+1, I ) )
964  300 CONTINUE
965  srnamt = 'SSTEVR'
966  CALL sstevr( 'N', 'A', n, d3, d4, vl, vu, il, iu, abstol,
967  $ m2, wa2, z, ldu, iwork, work, lwork,
968  $ iwork(2*n+1), liwork-2*n, iinfo )
969  IF( iinfo.NE.0 ) THEN
970  WRITE( nounit, fmt = 9999 )'SSTEVR(N,A)', iinfo, n,
971  $ jtype, ioldsd
972  info = abs( iinfo )
973  IF( iinfo.LT.0 ) THEN
974  RETURN
975  ELSE
976  result( 9 ) = ulpinv
977  GO TO 320
978  END IF
979  END IF
980 *
981 * Do test 9.
982 *
983  temp1 = zero
984  temp2 = zero
985  DO 310 j = 1, n
986  temp1 = max( temp1, abs( wa2( j ) ),
987  $ abs( eveigs( j ) ) )
988  temp2 = max( temp2, abs( wa2( j )-eveigs( j ) ) )
989  310 CONTINUE
990  result( 9 ) = temp2 / max( unfl,
991  $ ulp*max( temp1, temp2 ) )
992 *
993  320 CONTINUE
994 *
995 *
996  ntest = 10
997  DO 330 i = 1, n
998  d1( i ) = REAL( A( I, I ) )
999  330 CONTINUE
1000  DO 340 i = 1, n - 1
1001  d2( i ) = REAL( A( I+1, I ) )
1002  340 CONTINUE
1003  srnamt = 'SSTEVX'
1004  CALL sstevx( 'V', 'I', n, d1, d2, vl, vu, il, iu, abstol,
1005  $ m2, wa2, z, ldu, work, iwork,
1006  $ iwork( 5*n+1 ), iinfo )
1007  IF( iinfo.NE.0 ) THEN
1008  WRITE( nounit, fmt = 9999 )'SSTEVX(V,I)', iinfo, n,
1009  $ jtype, ioldsd
1010  info = abs( iinfo )
1011  IF( iinfo.LT.0 ) THEN
1012  RETURN
1013  ELSE
1014  result( 10 ) = ulpinv
1015  result( 11 ) = ulpinv
1016  result( 12 ) = ulpinv
1017  GO TO 380
1018  END IF
1019  END IF
1020 *
1021 * Do tests 10 and 11.
1022 *
1023  DO 350 i = 1, n
1024  d3( i ) = REAL( A( I, I ) )
1025  350 CONTINUE
1026  DO 360 i = 1, n - 1
1027  d4( i ) = REAL( A( I+1, I ) )
1028  360 CONTINUE
1029  CALL sstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1030  $ max( 1, m2 ), result( 10 ) )
1031 *
1032 *
1033  ntest = 12
1034  DO 370 i = 1, n - 1
1035  d4( i ) = REAL( A( I+1, I ) )
1036  370 CONTINUE
1037  srnamt = 'SSTEVX'
1038  CALL sstevx( 'N', 'I', n, d3, d4, vl, vu, il, iu, abstol,
1039  $ m3, wa3, z, ldu, work, iwork,
1040  $ iwork( 5*n+1 ), iinfo )
1041  IF( iinfo.NE.0 ) THEN
1042  WRITE( nounit, fmt = 9999 )'SSTEVX(N,I)', iinfo, n,
1043  $ jtype, ioldsd
1044  info = abs( iinfo )
1045  IF( iinfo.LT.0 ) THEN
1046  RETURN
1047  ELSE
1048  result( 12 ) = ulpinv
1049  GO TO 380
1050  END IF
1051  END IF
1052 *
1053 * Do test 12.
1054 *
1055  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1056  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1057  result( 12 ) = ( temp1+temp2 ) / max( unfl, ulp*temp3 )
1058 *
1059  380 CONTINUE
1060 *
1061  ntest = 12
1062  IF( n.GT.0 ) THEN
1063  IF( il.NE.1 ) THEN
1064  vl = wa1( il ) - max( half*
1065  $ ( wa1( il )-wa1( il-1 ) ), ten*ulp*temp3,
1066  $ ten*rtunfl )
1067  ELSE
1068  vl = wa1( 1 ) - max( half*( wa1( n )-wa1( 1 ) ),
1069  $ ten*ulp*temp3, ten*rtunfl )
1070  END IF
1071  IF( iu.NE.n ) THEN
1072  vu = wa1( iu ) + max( half*
1073  $ ( wa1( iu+1 )-wa1( iu ) ), ten*ulp*temp3,
1074  $ ten*rtunfl )
1075  ELSE
1076  vu = wa1( n ) + max( half*( wa1( n )-wa1( 1 ) ),
1077  $ ten*ulp*temp3, ten*rtunfl )
1078  END IF
1079  ELSE
1080  vl = zero
1081  vu = one
1082  END IF
1083 *
1084  DO 390 i = 1, n
1085  d1( i ) = REAL( A( I, I ) )
1086  390 CONTINUE
1087  DO 400 i = 1, n - 1
1088  d2( i ) = REAL( A( I+1, I ) )
1089  400 CONTINUE
1090  srnamt = 'SSTEVX'
1091  CALL sstevx( 'V', 'V', n, d1, d2, vl, vu, il, iu, abstol,
1092  $ m2, wa2, z, ldu, work, iwork,
1093  $ iwork( 5*n+1 ), iinfo )
1094  IF( iinfo.NE.0 ) THEN
1095  WRITE( nounit, fmt = 9999 )'SSTEVX(V,V)', iinfo, n,
1096  $ jtype, ioldsd
1097  info = abs( iinfo )
1098  IF( iinfo.LT.0 ) THEN
1099  RETURN
1100  ELSE
1101  result( 13 ) = ulpinv
1102  result( 14 ) = ulpinv
1103  result( 15 ) = ulpinv
1104  GO TO 440
1105  END IF
1106  END IF
1107 *
1108  IF( m2.EQ.0 .AND. n.GT.0 ) THEN
1109  result( 13 ) = ulpinv
1110  result( 14 ) = ulpinv
1111  result( 15 ) = ulpinv
1112  GO TO 440
1113  END IF
1114 *
1115 * Do tests 13 and 14.
1116 *
1117  DO 410 i = 1, n
1118  d3( i ) = REAL( A( I, I ) )
1119  410 CONTINUE
1120  DO 420 i = 1, n - 1
1121  d4( i ) = REAL( A( I+1, I ) )
1122  420 CONTINUE
1123  CALL sstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1124  $ max( 1, m2 ), result( 13 ) )
1125 *
1126  ntest = 15
1127  DO 430 i = 1, n - 1
1128  d4( i ) = REAL( A( I+1, I ) )
1129  430 CONTINUE
1130  srnamt = 'SSTEVX'
1131  CALL sstevx( 'N', 'V', n, d3, d4, vl, vu, il, iu, abstol,
1132  $ m3, wa3, z, ldu, work, iwork,
1133  $ iwork( 5*n+1 ), iinfo )
1134  IF( iinfo.NE.0 ) THEN
1135  WRITE( nounit, fmt = 9999 )'SSTEVX(N,V)', iinfo, n,
1136  $ jtype, ioldsd
1137  info = abs( iinfo )
1138  IF( iinfo.LT.0 ) THEN
1139  RETURN
1140  ELSE
1141  result( 15 ) = ulpinv
1142  GO TO 440
1143  END IF
1144  END IF
1145 *
1146 * Do test 15.
1147 *
1148  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1149  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1150  result( 15 ) = ( temp1+temp2 ) / max( unfl, temp3*ulp )
1151 *
1152  440 CONTINUE
1153 *
1154  ntest = 16
1155  DO 450 i = 1, n
1156  d1( i ) = REAL( A( I, I ) )
1157  450 CONTINUE
1158  DO 460 i = 1, n - 1
1159  d2( i ) = REAL( A( I+1, I ) )
1160  460 CONTINUE
1161  srnamt = 'SSTEVD'
1162  CALL sstevd( 'V', n, d1, d2, z, ldu, work, lwedc, iwork,
1163  $ liwedc, iinfo )
1164  IF( iinfo.NE.0 ) THEN
1165  WRITE( nounit, fmt = 9999 )'SSTEVD(V)', iinfo, n,
1166  $ jtype, ioldsd
1167  info = abs( iinfo )
1168  IF( iinfo.LT.0 ) THEN
1169  RETURN
1170  ELSE
1171  result( 16 ) = ulpinv
1172  result( 17 ) = ulpinv
1173  result( 18 ) = ulpinv
1174  GO TO 510
1175  END IF
1176  END IF
1177 *
1178 * Do tests 16 and 17.
1179 *
1180  DO 470 i = 1, n
1181  d3( i ) = REAL( A( I, I ) )
1182  470 CONTINUE
1183  DO 480 i = 1, n - 1
1184  d4( i ) = REAL( A( I+1, I ) )
1185  480 CONTINUE
1186  CALL sstt21( n, 0, d3, d4, d1, d2, z, ldu, work,
1187  $ result( 16 ) )
1188 *
1189  ntest = 18
1190  DO 490 i = 1, n - 1
1191  d4( i ) = REAL( A( I+1, I ) )
1192  490 CONTINUE
1193  srnamt = 'SSTEVD'
1194  CALL sstevd( 'N', n, d3, d4, z, ldu, work, lwedc, iwork,
1195  $ liwedc, iinfo )
1196  IF( iinfo.NE.0 ) THEN
1197  WRITE( nounit, fmt = 9999 )'SSTEVD(N)', iinfo, n,
1198  $ jtype, ioldsd
1199  info = abs( iinfo )
1200  IF( iinfo.LT.0 ) THEN
1201  RETURN
1202  ELSE
1203  result( 18 ) = ulpinv
1204  GO TO 510
1205  END IF
1206  END IF
1207 *
1208 * Do test 18.
1209 *
1210  temp1 = zero
1211  temp2 = zero
1212  DO 500 j = 1, n
1213  temp1 = max( temp1, abs( eveigs( j ) ),
1214  $ abs( d3( j ) ) )
1215  temp2 = max( temp2, abs( eveigs( j )-d3( j ) ) )
1216  500 CONTINUE
1217  result( 18 ) = temp2 / max( unfl,
1218  $ ulp*max( temp1, temp2 ) )
1219 *
1220  510 CONTINUE
1221 *
1222  ntest = 19
1223  DO 520 i = 1, n
1224  d1( i ) = REAL( A( I, I ) )
1225  520 CONTINUE
1226  DO 530 i = 1, n - 1
1227  d2( i ) = REAL( A( I+1, I ) )
1228  530 CONTINUE
1229  srnamt = 'SSTEVR'
1230  CALL sstevr( 'V', 'I', n, d1, d2, vl, vu, il, iu, abstol,
1231  $ m2, wa2, z, ldu, iwork, work, lwork,
1232  $ iwork(2*n+1), liwork-2*n, iinfo )
1233  IF( iinfo.NE.0 ) THEN
1234  WRITE( nounit, fmt = 9999 )'SSTEVR(V,I)', iinfo, n,
1235  $ jtype, ioldsd
1236  info = abs( iinfo )
1237  IF( iinfo.LT.0 ) THEN
1238  RETURN
1239  ELSE
1240  result( 19 ) = ulpinv
1241  result( 20 ) = ulpinv
1242  result( 21 ) = ulpinv
1243  GO TO 570
1244  END IF
1245  END IF
1246 *
1247 * DO tests 19 and 20.
1248 *
1249  DO 540 i = 1, n
1250  d3( i ) = REAL( A( I, I ) )
1251  540 CONTINUE
1252  DO 550 i = 1, n - 1
1253  d4( i ) = REAL( A( I+1, I ) )
1254  550 CONTINUE
1255  CALL sstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1256  $ max( 1, m2 ), result( 19 ) )
1257 *
1258 *
1259  ntest = 21
1260  DO 560 i = 1, n - 1
1261  d4( i ) = REAL( A( I+1, I ) )
1262  560 CONTINUE
1263  srnamt = 'SSTEVR'
1264  CALL sstevr( 'N', 'I', n, d3, d4, vl, vu, il, iu, abstol,
1265  $ m3, wa3, z, ldu, iwork, work, lwork,
1266  $ iwork(2*n+1), liwork-2*n, iinfo )
1267  IF( iinfo.NE.0 ) THEN
1268  WRITE( nounit, fmt = 9999 )'SSTEVR(N,I)', iinfo, n,
1269  $ jtype, ioldsd
1270  info = abs( iinfo )
1271  IF( iinfo.LT.0 ) THEN
1272  RETURN
1273  ELSE
1274  result( 21 ) = ulpinv
1275  GO TO 570
1276  END IF
1277  END IF
1278 *
1279 * Do test 21.
1280 *
1281  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1282  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1283  result( 21 ) = ( temp1+temp2 ) / max( unfl, ulp*temp3 )
1284 *
1285  570 CONTINUE
1286 *
1287  ntest = 21
1288  IF( n.GT.0 ) THEN
1289  IF( il.NE.1 ) THEN
1290  vl = wa1( il ) - max( half*
1291  $ ( wa1( il )-wa1( il-1 ) ), ten*ulp*temp3,
1292  $ ten*rtunfl )
1293  ELSE
1294  vl = wa1( 1 ) - max( half*( wa1( n )-wa1( 1 ) ),
1295  $ ten*ulp*temp3, ten*rtunfl )
1296  END IF
1297  IF( iu.NE.n ) THEN
1298  vu = wa1( iu ) + max( half*
1299  $ ( wa1( iu+1 )-wa1( iu ) ), ten*ulp*temp3,
1300  $ ten*rtunfl )
1301  ELSE
1302  vu = wa1( n ) + max( half*( wa1( n )-wa1( 1 ) ),
1303  $ ten*ulp*temp3, ten*rtunfl )
1304  END IF
1305  ELSE
1306  vl = zero
1307  vu = one
1308  END IF
1309 *
1310  DO 580 i = 1, n
1311  d1( i ) = REAL( A( I, I ) )
1312  580 CONTINUE
1313  DO 590 i = 1, n - 1
1314  d2( i ) = REAL( A( I+1, I ) )
1315  590 CONTINUE
1316  srnamt = 'SSTEVR'
1317  CALL sstevr( 'V', 'V', n, d1, d2, vl, vu, il, iu, abstol,
1318  $ m2, wa2, z, ldu, iwork, work, lwork,
1319  $ iwork(2*n+1), liwork-2*n, iinfo )
1320  IF( iinfo.NE.0 ) THEN
1321  WRITE( nounit, fmt = 9999 )'SSTEVR(V,V)', iinfo, n,
1322  $ jtype, ioldsd
1323  info = abs( iinfo )
1324  IF( iinfo.LT.0 ) THEN
1325  RETURN
1326  ELSE
1327  result( 22 ) = ulpinv
1328  result( 23 ) = ulpinv
1329  result( 24 ) = ulpinv
1330  GO TO 630
1331  END IF
1332  END IF
1333 *
1334  IF( m2.EQ.0 .AND. n.GT.0 ) THEN
1335  result( 22 ) = ulpinv
1336  result( 23 ) = ulpinv
1337  result( 24 ) = ulpinv
1338  GO TO 630
1339  END IF
1340 *
1341 * Do tests 22 and 23.
1342 *
1343  DO 600 i = 1, n
1344  d3( i ) = REAL( A( I, I ) )
1345  600 CONTINUE
1346  DO 610 i = 1, n - 1
1347  d4( i ) = REAL( A( I+1, I ) )
1348  610 CONTINUE
1349  CALL sstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1350  $ max( 1, m2 ), result( 22 ) )
1351 *
1352  ntest = 24
1353  DO 620 i = 1, n - 1
1354  d4( i ) = REAL( A( I+1, I ) )
1355  620 CONTINUE
1356  srnamt = 'SSTEVR'
1357  CALL sstevr( 'N', 'V', n, d3, d4, vl, vu, il, iu, abstol,
1358  $ m3, wa3, z, ldu, iwork, work, lwork,
1359  $ iwork(2*n+1), liwork-2*n, iinfo )
1360  IF( iinfo.NE.0 ) THEN
1361  WRITE( nounit, fmt = 9999 )'SSTEVR(N,V)', iinfo, n,
1362  $ jtype, ioldsd
1363  info = abs( iinfo )
1364  IF( iinfo.LT.0 ) THEN
1365  RETURN
1366  ELSE
1367  result( 24 ) = ulpinv
1368  GO TO 630
1369  END IF
1370  END IF
1371 *
1372 * Do test 24.
1373 *
1374  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1375  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1376  result( 24 ) = ( temp1+temp2 ) / max( unfl, temp3*ulp )
1377 *
1378  630 CONTINUE
1379 *
1380 *
1381 *
1382  ELSE
1383 *
1384  DO 640 i = 1, 24
1385  result( i ) = zero
1386  640 CONTINUE
1387  ntest = 24
1388  END IF
1389 *
1390 * Perform remaining tests storing upper or lower triangular
1391 * part of matrix.
1392 *
1393  DO 1720 iuplo = 0, 1
1394  IF( iuplo.EQ.0 ) THEN
1395  uplo = 'L'
1396  ELSE
1397  uplo = 'U'
1398  END IF
1399 *
1400 * 4) Call SSYEV and SSYEVX.
1401 *
1402  CALL slacpy( ' ', n, n, a, lda, v, ldu )
1403 *
1404  ntest = ntest + 1
1405  srnamt = 'SSYEV'
1406  CALL ssyev( 'V', uplo, n, a, ldu, d1, work, lwork,
1407  $ iinfo )
1408  IF( iinfo.NE.0 ) THEN
1409  WRITE( nounit, fmt = 9999 )'SSYEV(V,' // uplo // ')',
1410  $ iinfo, n, jtype, ioldsd
1411  info = abs( iinfo )
1412  IF( iinfo.LT.0 ) THEN
1413  RETURN
1414  ELSE
1415  result( ntest ) = ulpinv
1416  result( ntest+1 ) = ulpinv
1417  result( ntest+2 ) = ulpinv
1418  GO TO 660
1419  END IF
1420  END IF
1421 *
1422 * Do tests 25 and 26 (or +54)
1423 *
1424  CALL ssyt21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
1425  $ ldu, tau, work, result( ntest ) )
1426 *
1427  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1428 *
1429  ntest = ntest + 2
1430  srnamt = 'SSYEV_2STAGE'
1431  CALL ssyev_2stage( 'N', uplo, n, a, ldu, d3, work, lwork,
1432  $ iinfo )
1433  IF( iinfo.NE.0 ) THEN
1434  WRITE( nounit, fmt = 9999 )
1435  $ 'SSYEV_2STAGE(N,' // uplo // ')',
1436  $ iinfo, n, jtype, ioldsd
1437  info = abs( iinfo )
1438  IF( iinfo.LT.0 ) THEN
1439  RETURN
1440  ELSE
1441  result( ntest ) = ulpinv
1442  GO TO 660
1443  END IF
1444  END IF
1445 *
1446 * Do test 27 (or +54)
1447 *
1448  temp1 = zero
1449  temp2 = zero
1450  DO 650 j = 1, n
1451  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1452  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1453  650 CONTINUE
1454  result( ntest ) = temp2 / max( unfl,
1455  $ ulp*max( temp1, temp2 ) )
1456 *
1457  660 CONTINUE
1458  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1459 *
1460  ntest = ntest + 1
1461 *
1462  IF( n.GT.0 ) THEN
1463  temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1464  IF( il.NE.1 ) THEN
1465  vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1466  $ ten*ulp*temp3, ten*rtunfl )
1467  ELSE IF( n.GT.0 ) THEN
1468  vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1469  $ ten*ulp*temp3, ten*rtunfl )
1470  END IF
1471  IF( iu.NE.n ) THEN
1472  vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1473  $ ten*ulp*temp3, ten*rtunfl )
1474  ELSE IF( n.GT.0 ) THEN
1475  vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1476  $ ten*ulp*temp3, ten*rtunfl )
1477  END IF
1478  ELSE
1479  temp3 = zero
1480  vl = zero
1481  vu = one
1482  END IF
1483 *
1484  srnamt = 'SSYEVX'
1485  CALL ssyevx( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
1486  $ abstol, m, wa1, z, ldu, work, lwork, iwork,
1487  $ iwork( 5*n+1 ), iinfo )
1488  IF( iinfo.NE.0 ) THEN
1489  WRITE( nounit, fmt = 9999 )'SSYEVX(V,A,' // uplo //
1490  $ ')', iinfo, n, jtype, ioldsd
1491  info = abs( iinfo )
1492  IF( iinfo.LT.0 ) THEN
1493  RETURN
1494  ELSE
1495  result( ntest ) = ulpinv
1496  result( ntest+1 ) = ulpinv
1497  result( ntest+2 ) = ulpinv
1498  GO TO 680
1499  END IF
1500  END IF
1501 *
1502 * Do tests 28 and 29 (or +54)
1503 *
1504  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1505 *
1506  CALL ssyt21( 1, uplo, n, 0, a, ldu, d1, d2, z, ldu, v,
1507  $ ldu, tau, work, result( ntest ) )
1508 *
1509  ntest = ntest + 2
1510  srnamt = 'SSYEVX_2STAGE'
1511  CALL ssyevx_2stage( 'N', 'A', uplo, n, a, ldu, vl, vu,
1512  $ il, iu, abstol, m2, wa2, z, ldu, work,
1513  $ lwork, iwork, iwork( 5*n+1 ), iinfo )
1514  IF( iinfo.NE.0 ) THEN
1515  WRITE( nounit, fmt = 9999 )
1516  $ 'SSYEVX_2STAGE(N,A,' // uplo //
1517  $ ')', iinfo, n, jtype, ioldsd
1518  info = abs( iinfo )
1519  IF( iinfo.LT.0 ) THEN
1520  RETURN
1521  ELSE
1522  result( ntest ) = ulpinv
1523  GO TO 680
1524  END IF
1525  END IF
1526 *
1527 * Do test 30 (or +54)
1528 *
1529  temp1 = zero
1530  temp2 = zero
1531  DO 670 j = 1, n
1532  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1533  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1534  670 CONTINUE
1535  result( ntest ) = temp2 / max( unfl,
1536  $ ulp*max( temp1, temp2 ) )
1537 *
1538  680 CONTINUE
1539 *
1540  ntest = ntest + 1
1541  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1542  srnamt = 'SSYEVX'
1543  CALL ssyevx( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
1544  $ abstol, m2, wa2, z, ldu, work, lwork, iwork,
1545  $ iwork( 5*n+1 ), iinfo )
1546  IF( iinfo.NE.0 ) THEN
1547  WRITE( nounit, fmt = 9999 )'SSYEVX(V,I,' // uplo //
1548  $ ')', iinfo, n, jtype, ioldsd
1549  info = abs( iinfo )
1550  IF( iinfo.LT.0 ) THEN
1551  RETURN
1552  ELSE
1553  result( ntest ) = ulpinv
1554  result( ntest+1 ) = ulpinv
1555  result( ntest+2 ) = ulpinv
1556  GO TO 690
1557  END IF
1558  END IF
1559 *
1560 * Do tests 31 and 32 (or +54)
1561 *
1562  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1563 *
1564  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1565  $ v, ldu, tau, work, result( ntest ) )
1566 *
1567  ntest = ntest + 2
1568  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1569  srnamt = 'SSYEVX_2STAGE'
1570  CALL ssyevx_2stage( 'N', 'I', uplo, n, a, ldu, vl, vu,
1571  $ il, iu, abstol, m3, wa3, z, ldu, work,
1572  $ lwork, iwork, iwork( 5*n+1 ), iinfo )
1573  IF( iinfo.NE.0 ) THEN
1574  WRITE( nounit, fmt = 9999 )
1575  $ 'SSYEVX_2STAGE(N,I,' // uplo //
1576  $ ')', iinfo, n, jtype, ioldsd
1577  info = abs( iinfo )
1578  IF( iinfo.LT.0 ) THEN
1579  RETURN
1580  ELSE
1581  result( ntest ) = ulpinv
1582  GO TO 690
1583  END IF
1584  END IF
1585 *
1586 * Do test 33 (or +54)
1587 *
1588  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1589  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1590  result( ntest ) = ( temp1+temp2 ) /
1591  $ max( unfl, ulp*temp3 )
1592  690 CONTINUE
1593 *
1594  ntest = ntest + 1
1595  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1596  srnamt = 'SSYEVX'
1597  CALL ssyevx( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
1598  $ abstol, m2, wa2, z, ldu, work, lwork, iwork,
1599  $ iwork( 5*n+1 ), iinfo )
1600  IF( iinfo.NE.0 ) THEN
1601  WRITE( nounit, fmt = 9999 )'SSYEVX(V,V,' // uplo //
1602  $ ')', iinfo, n, jtype, ioldsd
1603  info = abs( iinfo )
1604  IF( iinfo.LT.0 ) THEN
1605  RETURN
1606  ELSE
1607  result( ntest ) = ulpinv
1608  result( ntest+1 ) = ulpinv
1609  result( ntest+2 ) = ulpinv
1610  GO TO 700
1611  END IF
1612  END IF
1613 *
1614 * Do tests 34 and 35 (or +54)
1615 *
1616  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1617 *
1618  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1619  $ v, ldu, tau, work, result( ntest ) )
1620 *
1621  ntest = ntest + 2
1622  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1623  srnamt = 'SSYEVX_2STAGE'
1624  CALL ssyevx_2stage( 'N', 'V', uplo, n, a, ldu, vl, vu,
1625  $ il, iu, abstol, m3, wa3, z, ldu, work,
1626  $ lwork, iwork, iwork( 5*n+1 ), iinfo )
1627  IF( iinfo.NE.0 ) THEN
1628  WRITE( nounit, fmt = 9999 )
1629  $ 'SSYEVX_2STAGE(N,V,' // uplo //
1630  $ ')', iinfo, n, jtype, ioldsd
1631  info = abs( iinfo )
1632  IF( iinfo.LT.0 ) THEN
1633  RETURN
1634  ELSE
1635  result( ntest ) = ulpinv
1636  GO TO 700
1637  END IF
1638  END IF
1639 *
1640  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1641  result( ntest ) = ulpinv
1642  GO TO 700
1643  END IF
1644 *
1645 * Do test 36 (or +54)
1646 *
1647  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1648  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1649  IF( n.GT.0 ) THEN
1650  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1651  ELSE
1652  temp3 = zero
1653  END IF
1654  result( ntest ) = ( temp1+temp2 ) /
1655  $ max( unfl, temp3*ulp )
1656 *
1657  700 CONTINUE
1658 *
1659 * 5) Call SSPEV and SSPEVX.
1660 *
1661  CALL slacpy( ' ', n, n, v, ldu, a, lda )
1662 *
1663 * Load array WORK with the upper or lower triangular
1664 * part of the matrix in packed form.
1665 *
1666  IF( iuplo.EQ.1 ) THEN
1667  indx = 1
1668  DO 720 j = 1, n
1669  DO 710 i = 1, j
1670  work( indx ) = a( i, j )
1671  indx = indx + 1
1672  710 CONTINUE
1673  720 CONTINUE
1674  ELSE
1675  indx = 1
1676  DO 740 j = 1, n
1677  DO 730 i = j, n
1678  work( indx ) = a( i, j )
1679  indx = indx + 1
1680  730 CONTINUE
1681  740 CONTINUE
1682  END IF
1683 *
1684  ntest = ntest + 1
1685  srnamt = 'SSPEV'
1686  CALL sspev( 'V', uplo, n, work, d1, z, ldu, v, iinfo )
1687  IF( iinfo.NE.0 ) THEN
1688  WRITE( nounit, fmt = 9999 )'SSPEV(V,' // uplo // ')',
1689  $ iinfo, n, jtype, ioldsd
1690  info = abs( iinfo )
1691  IF( iinfo.LT.0 ) THEN
1692  RETURN
1693  ELSE
1694  result( ntest ) = ulpinv
1695  result( ntest+1 ) = ulpinv
1696  result( ntest+2 ) = ulpinv
1697  GO TO 800
1698  END IF
1699  END IF
1700 *
1701 * Do tests 37 and 38 (or +54)
1702 *
1703  CALL ssyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1704  $ ldu, tau, work, result( ntest ) )
1705 *
1706  IF( iuplo.EQ.1 ) THEN
1707  indx = 1
1708  DO 760 j = 1, n
1709  DO 750 i = 1, j
1710  work( indx ) = a( i, j )
1711  indx = indx + 1
1712  750 CONTINUE
1713  760 CONTINUE
1714  ELSE
1715  indx = 1
1716  DO 780 j = 1, n
1717  DO 770 i = j, n
1718  work( indx ) = a( i, j )
1719  indx = indx + 1
1720  770 CONTINUE
1721  780 CONTINUE
1722  END IF
1723 *
1724  ntest = ntest + 2
1725  srnamt = 'SSPEV'
1726  CALL sspev( 'N', uplo, n, work, d3, z, ldu, v, iinfo )
1727  IF( iinfo.NE.0 ) THEN
1728  WRITE( nounit, fmt = 9999 )'SSPEV(N,' // uplo // ')',
1729  $ iinfo, n, jtype, ioldsd
1730  info = abs( iinfo )
1731  IF( iinfo.LT.0 ) THEN
1732  RETURN
1733  ELSE
1734  result( ntest ) = ulpinv
1735  GO TO 800
1736  END IF
1737  END IF
1738 *
1739 * Do test 39 (or +54)
1740 *
1741  temp1 = zero
1742  temp2 = zero
1743  DO 790 j = 1, n
1744  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1745  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1746  790 CONTINUE
1747  result( ntest ) = temp2 / max( unfl,
1748  $ ulp*max( temp1, temp2 ) )
1749 *
1750 * Load array WORK with the upper or lower triangular part
1751 * of the matrix in packed form.
1752 *
1753  800 CONTINUE
1754  IF( iuplo.EQ.1 ) THEN
1755  indx = 1
1756  DO 820 j = 1, n
1757  DO 810 i = 1, j
1758  work( indx ) = a( i, j )
1759  indx = indx + 1
1760  810 CONTINUE
1761  820 CONTINUE
1762  ELSE
1763  indx = 1
1764  DO 840 j = 1, n
1765  DO 830 i = j, n
1766  work( indx ) = a( i, j )
1767  indx = indx + 1
1768  830 CONTINUE
1769  840 CONTINUE
1770  END IF
1771 *
1772  ntest = ntest + 1
1773 *
1774  IF( n.GT.0 ) THEN
1775  temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1776  IF( il.NE.1 ) THEN
1777  vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1778  $ ten*ulp*temp3, ten*rtunfl )
1779  ELSE IF( n.GT.0 ) THEN
1780  vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1781  $ ten*ulp*temp3, ten*rtunfl )
1782  END IF
1783  IF( iu.NE.n ) THEN
1784  vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1785  $ ten*ulp*temp3, ten*rtunfl )
1786  ELSE IF( n.GT.0 ) THEN
1787  vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1788  $ ten*ulp*temp3, ten*rtunfl )
1789  END IF
1790  ELSE
1791  temp3 = zero
1792  vl = zero
1793  vu = one
1794  END IF
1795 *
1796  srnamt = 'SSPEVX'
1797  CALL sspevx( 'V', 'A', uplo, n, work, vl, vu, il, iu,
1798  $ abstol, m, wa1, z, ldu, v, iwork,
1799  $ iwork( 5*n+1 ), iinfo )
1800  IF( iinfo.NE.0 ) THEN
1801  WRITE( nounit, fmt = 9999 )'SSPEVX(V,A,' // uplo //
1802  $ ')', iinfo, n, jtype, ioldsd
1803  info = abs( iinfo )
1804  IF( iinfo.LT.0 ) THEN
1805  RETURN
1806  ELSE
1807  result( ntest ) = ulpinv
1808  result( ntest+1 ) = ulpinv
1809  result( ntest+2 ) = ulpinv
1810  GO TO 900
1811  END IF
1812  END IF
1813 *
1814 * Do tests 40 and 41 (or +54)
1815 *
1816  CALL ssyt21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1817  $ ldu, tau, work, result( ntest ) )
1818 *
1819  ntest = ntest + 2
1820 *
1821  IF( iuplo.EQ.1 ) THEN
1822  indx = 1
1823  DO 860 j = 1, n
1824  DO 850 i = 1, j
1825  work( indx ) = a( i, j )
1826  indx = indx + 1
1827  850 CONTINUE
1828  860 CONTINUE
1829  ELSE
1830  indx = 1
1831  DO 880 j = 1, n
1832  DO 870 i = j, n
1833  work( indx ) = a( i, j )
1834  indx = indx + 1
1835  870 CONTINUE
1836  880 CONTINUE
1837  END IF
1838 *
1839  srnamt = 'SSPEVX'
1840  CALL sspevx( 'N', 'A', uplo, n, work, vl, vu, il, iu,
1841  $ abstol, m2, wa2, z, ldu, v, iwork,
1842  $ iwork( 5*n+1 ), iinfo )
1843  IF( iinfo.NE.0 ) THEN
1844  WRITE( nounit, fmt = 9999 )'SSPEVX(N,A,' // uplo //
1845  $ ')', iinfo, n, jtype, ioldsd
1846  info = abs( iinfo )
1847  IF( iinfo.LT.0 ) THEN
1848  RETURN
1849  ELSE
1850  result( ntest ) = ulpinv
1851  GO TO 900
1852  END IF
1853  END IF
1854 *
1855 * Do test 42 (or +54)
1856 *
1857  temp1 = zero
1858  temp2 = zero
1859  DO 890 j = 1, n
1860  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1861  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1862  890 CONTINUE
1863  result( ntest ) = temp2 / max( unfl,
1864  $ ulp*max( temp1, temp2 ) )
1865 *
1866  900 CONTINUE
1867  IF( iuplo.EQ.1 ) THEN
1868  indx = 1
1869  DO 920 j = 1, n
1870  DO 910 i = 1, j
1871  work( indx ) = a( i, j )
1872  indx = indx + 1
1873  910 CONTINUE
1874  920 CONTINUE
1875  ELSE
1876  indx = 1
1877  DO 940 j = 1, n
1878  DO 930 i = j, n
1879  work( indx ) = a( i, j )
1880  indx = indx + 1
1881  930 CONTINUE
1882  940 CONTINUE
1883  END IF
1884 *
1885  ntest = ntest + 1
1886 *
1887  srnamt = 'SSPEVX'
1888  CALL sspevx( 'V', 'I', uplo, n, work, vl, vu, il, iu,
1889  $ abstol, m2, wa2, z, ldu, v, iwork,
1890  $ iwork( 5*n+1 ), iinfo )
1891  IF( iinfo.NE.0 ) THEN
1892  WRITE( nounit, fmt = 9999 )'SSPEVX(V,I,' // uplo //
1893  $ ')', iinfo, n, jtype, ioldsd
1894  info = abs( iinfo )
1895  IF( iinfo.LT.0 ) THEN
1896  RETURN
1897  ELSE
1898  result( ntest ) = ulpinv
1899  result( ntest+1 ) = ulpinv
1900  result( ntest+2 ) = ulpinv
1901  GO TO 990
1902  END IF
1903  END IF
1904 *
1905 * Do tests 43 and 44 (or +54)
1906 *
1907  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1908  $ v, ldu, tau, work, result( ntest ) )
1909 *
1910  ntest = ntest + 2
1911 *
1912  IF( iuplo.EQ.1 ) THEN
1913  indx = 1
1914  DO 960 j = 1, n
1915  DO 950 i = 1, j
1916  work( indx ) = a( i, j )
1917  indx = indx + 1
1918  950 CONTINUE
1919  960 CONTINUE
1920  ELSE
1921  indx = 1
1922  DO 980 j = 1, n
1923  DO 970 i = j, n
1924  work( indx ) = a( i, j )
1925  indx = indx + 1
1926  970 CONTINUE
1927  980 CONTINUE
1928  END IF
1929 *
1930  srnamt = 'SSPEVX'
1931  CALL sspevx( 'N', 'I', uplo, n, work, vl, vu, il, iu,
1932  $ abstol, m3, wa3, z, ldu, v, iwork,
1933  $ iwork( 5*n+1 ), iinfo )
1934  IF( iinfo.NE.0 ) THEN
1935  WRITE( nounit, fmt = 9999 )'SSPEVX(N,I,' // uplo //
1936  $ ')', iinfo, n, jtype, ioldsd
1937  info = abs( iinfo )
1938  IF( iinfo.LT.0 ) THEN
1939  RETURN
1940  ELSE
1941  result( ntest ) = ulpinv
1942  GO TO 990
1943  END IF
1944  END IF
1945 *
1946  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1947  result( ntest ) = ulpinv
1948  GO TO 990
1949  END IF
1950 *
1951 * Do test 45 (or +54)
1952 *
1953  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1954  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1955  IF( n.GT.0 ) THEN
1956  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1957  ELSE
1958  temp3 = zero
1959  END IF
1960  result( ntest ) = ( temp1+temp2 ) /
1961  $ max( unfl, temp3*ulp )
1962 *
1963  990 CONTINUE
1964  IF( iuplo.EQ.1 ) THEN
1965  indx = 1
1966  DO 1010 j = 1, n
1967  DO 1000 i = 1, j
1968  work( indx ) = a( i, j )
1969  indx = indx + 1
1970  1000 CONTINUE
1971  1010 CONTINUE
1972  ELSE
1973  indx = 1
1974  DO 1030 j = 1, n
1975  DO 1020 i = j, n
1976  work( indx ) = a( i, j )
1977  indx = indx + 1
1978  1020 CONTINUE
1979  1030 CONTINUE
1980  END IF
1981 *
1982  ntest = ntest + 1
1983 *
1984  srnamt = 'SSPEVX'
1985  CALL sspevx( 'V', 'V', uplo, n, work, vl, vu, il, iu,
1986  $ abstol, m2, wa2, z, ldu, v, iwork,
1987  $ iwork( 5*n+1 ), iinfo )
1988  IF( iinfo.NE.0 ) THEN
1989  WRITE( nounit, fmt = 9999 )'SSPEVX(V,V,' // uplo //
1990  $ ')', iinfo, n, jtype, ioldsd
1991  info = abs( iinfo )
1992  IF( iinfo.LT.0 ) THEN
1993  RETURN
1994  ELSE
1995  result( ntest ) = ulpinv
1996  result( ntest+1 ) = ulpinv
1997  result( ntest+2 ) = ulpinv
1998  GO TO 1080
1999  END IF
2000  END IF
2001 *
2002 * Do tests 46 and 47 (or +54)
2003 *
2004  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2005  $ v, ldu, tau, work, result( ntest ) )
2006 *
2007  ntest = ntest + 2
2008 *
2009  IF( iuplo.EQ.1 ) THEN
2010  indx = 1
2011  DO 1050 j = 1, n
2012  DO 1040 i = 1, j
2013  work( indx ) = a( i, j )
2014  indx = indx + 1
2015  1040 CONTINUE
2016  1050 CONTINUE
2017  ELSE
2018  indx = 1
2019  DO 1070 j = 1, n
2020  DO 1060 i = j, n
2021  work( indx ) = a( i, j )
2022  indx = indx + 1
2023  1060 CONTINUE
2024  1070 CONTINUE
2025  END IF
2026 *
2027  srnamt = 'SSPEVX'
2028  CALL sspevx( 'N', 'V', uplo, n, work, vl, vu, il, iu,
2029  $ abstol, m3, wa3, z, ldu, v, iwork,
2030  $ iwork( 5*n+1 ), iinfo )
2031  IF( iinfo.NE.0 ) THEN
2032  WRITE( nounit, fmt = 9999 )'SSPEVX(N,V,' // uplo //
2033  $ ')', iinfo, n, jtype, ioldsd
2034  info = abs( iinfo )
2035  IF( iinfo.LT.0 ) THEN
2036  RETURN
2037  ELSE
2038  result( ntest ) = ulpinv
2039  GO TO 1080
2040  END IF
2041  END IF
2042 *
2043  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2044  result( ntest ) = ulpinv
2045  GO TO 1080
2046  END IF
2047 *
2048 * Do test 48 (or +54)
2049 *
2050  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2051  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2052  IF( n.GT.0 ) THEN
2053  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2054  ELSE
2055  temp3 = zero
2056  END IF
2057  result( ntest ) = ( temp1+temp2 ) /
2058  $ max( unfl, temp3*ulp )
2059 *
2060  1080 CONTINUE
2061 *
2062 * 6) Call SSBEV and SSBEVX.
2063 *
2064  IF( jtype.LE.7 ) THEN
2065  kd = 1
2066  ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
2067  kd = max( n-1, 0 )
2068  ELSE
2069  kd = ihbw
2070  END IF
2071 *
2072 * Load array V with the upper or lower triangular part
2073 * of the matrix in band form.
2074 *
2075  IF( iuplo.EQ.1 ) THEN
2076  DO 1100 j = 1, n
2077  DO 1090 i = max( 1, j-kd ), j
2078  v( kd+1+i-j, j ) = a( i, j )
2079  1090 CONTINUE
2080  1100 CONTINUE
2081  ELSE
2082  DO 1120 j = 1, n
2083  DO 1110 i = j, min( n, j+kd )
2084  v( 1+i-j, j ) = a( i, j )
2085  1110 CONTINUE
2086  1120 CONTINUE
2087  END IF
2088 *
2089  ntest = ntest + 1
2090  srnamt = 'SSBEV'
2091  CALL ssbev( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
2092  $ iinfo )
2093  IF( iinfo.NE.0 ) THEN
2094  WRITE( nounit, fmt = 9999 )'SSBEV(V,' // uplo // ')',
2095  $ iinfo, n, jtype, ioldsd
2096  info = abs( iinfo )
2097  IF( iinfo.LT.0 ) THEN
2098  RETURN
2099  ELSE
2100  result( ntest ) = ulpinv
2101  result( ntest+1 ) = ulpinv
2102  result( ntest+2 ) = ulpinv
2103  GO TO 1180
2104  END IF
2105  END IF
2106 *
2107 * Do tests 49 and 50 (or ... )
2108 *
2109  CALL ssyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
2110  $ ldu, tau, work, result( ntest ) )
2111 *
2112  IF( iuplo.EQ.1 ) THEN
2113  DO 1140 j = 1, n
2114  DO 1130 i = max( 1, j-kd ), j
2115  v( kd+1+i-j, j ) = a( i, j )
2116  1130 CONTINUE
2117  1140 CONTINUE
2118  ELSE
2119  DO 1160 j = 1, n
2120  DO 1150 i = j, min( n, j+kd )
2121  v( 1+i-j, j ) = a( i, j )
2122  1150 CONTINUE
2123  1160 CONTINUE
2124  END IF
2125 *
2126  ntest = ntest + 2
2127  srnamt = 'SSBEV_2STAGE'
2128  CALL ssbev_2stage( 'N', uplo, n, kd, v, ldu, d3, z, ldu,
2129  $ work, lwork, iinfo )
2130  IF( iinfo.NE.0 ) THEN
2131  WRITE( nounit, fmt = 9999 )
2132  $ 'SSBEV_2STAGE(N,' // uplo // ')',
2133  $ iinfo, n, jtype, ioldsd
2134  info = abs( iinfo )
2135  IF( iinfo.LT.0 ) THEN
2136  RETURN
2137  ELSE
2138  result( ntest ) = ulpinv
2139  GO TO 1180
2140  END IF
2141  END IF
2142 *
2143 * Do test 51 (or +54)
2144 *
2145  temp1 = zero
2146  temp2 = zero
2147  DO 1170 j = 1, n
2148  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2149  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2150  1170 CONTINUE
2151  result( ntest ) = temp2 / max( unfl,
2152  $ ulp*max( temp1, temp2 ) )
2153 *
2154 * Load array V with the upper or lower triangular part
2155 * of the matrix in band form.
2156 *
2157  1180 CONTINUE
2158  IF( iuplo.EQ.1 ) THEN
2159  DO 1200 j = 1, n
2160  DO 1190 i = max( 1, j-kd ), j
2161  v( kd+1+i-j, j ) = a( i, j )
2162  1190 CONTINUE
2163  1200 CONTINUE
2164  ELSE
2165  DO 1220 j = 1, n
2166  DO 1210 i = j, min( n, j+kd )
2167  v( 1+i-j, j ) = a( i, j )
2168  1210 CONTINUE
2169  1220 CONTINUE
2170  END IF
2171 *
2172  ntest = ntest + 1
2173  srnamt = 'SSBEVX'
2174  CALL ssbevx( 'V', 'A', uplo, n, kd, v, ldu, u, ldu, vl,
2175  $ vu, il, iu, abstol, m, wa2, z, ldu, work,
2176  $ iwork, iwork( 5*n+1 ), iinfo )
2177  IF( iinfo.NE.0 ) THEN
2178  WRITE( nounit, fmt = 9999 )'SSBEVX(V,A,' // uplo //
2179  $ ')', iinfo, n, jtype, ioldsd
2180  info = abs( iinfo )
2181  IF( iinfo.LT.0 ) THEN
2182  RETURN
2183  ELSE
2184  result( ntest ) = ulpinv
2185  result( ntest+1 ) = ulpinv
2186  result( ntest+2 ) = ulpinv
2187  GO TO 1280
2188  END IF
2189  END IF
2190 *
2191 * Do tests 52 and 53 (or +54)
2192 *
2193  CALL ssyt21( 1, uplo, n, 0, a, ldu, wa2, d2, z, ldu, v,
2194  $ ldu, tau, work, result( ntest ) )
2195 *
2196  ntest = ntest + 2
2197 *
2198  IF( iuplo.EQ.1 ) THEN
2199  DO 1240 j = 1, n
2200  DO 1230 i = max( 1, j-kd ), j
2201  v( kd+1+i-j, j ) = a( i, j )
2202  1230 CONTINUE
2203  1240 CONTINUE
2204  ELSE
2205  DO 1260 j = 1, n
2206  DO 1250 i = j, min( n, j+kd )
2207  v( 1+i-j, j ) = a( i, j )
2208  1250 CONTINUE
2209  1260 CONTINUE
2210  END IF
2211 *
2212  srnamt = 'SSBEVX_2STAGE'
2213  CALL ssbevx_2stage( 'N', 'A', uplo, n, kd, v, ldu,
2214  $ u, ldu, vl, vu, il, iu, abstol, m3, wa3,
2215  $ z, ldu, work, lwork, iwork, iwork( 5*n+1 ),
2216  $ iinfo )
2217  IF( iinfo.NE.0 ) THEN
2218  WRITE( nounit, fmt = 9999 )
2219  $ 'SSBEVX_2STAGE(N,A,' // uplo //
2220  $ ')', iinfo, n, jtype, ioldsd
2221  info = abs( iinfo )
2222  IF( iinfo.LT.0 ) THEN
2223  RETURN
2224  ELSE
2225  result( ntest ) = ulpinv
2226  GO TO 1280
2227  END IF
2228  END IF
2229 *
2230 * Do test 54 (or +54)
2231 *
2232  temp1 = zero
2233  temp2 = zero
2234  DO 1270 j = 1, n
2235  temp1 = max( temp1, abs( wa2( j ) ), abs( wa3( j ) ) )
2236  temp2 = max( temp2, abs( wa2( j )-wa3( j ) ) )
2237  1270 CONTINUE
2238  result( ntest ) = temp2 / max( unfl,
2239  $ ulp*max( temp1, temp2 ) )
2240 *
2241  1280 CONTINUE
2242  ntest = ntest + 1
2243  IF( iuplo.EQ.1 ) THEN
2244  DO 1300 j = 1, n
2245  DO 1290 i = max( 1, j-kd ), j
2246  v( kd+1+i-j, j ) = a( i, j )
2247  1290 CONTINUE
2248  1300 CONTINUE
2249  ELSE
2250  DO 1320 j = 1, n
2251  DO 1310 i = j, min( n, j+kd )
2252  v( 1+i-j, j ) = a( i, j )
2253  1310 CONTINUE
2254  1320 CONTINUE
2255  END IF
2256 *
2257  srnamt = 'SSBEVX'
2258  CALL ssbevx( 'V', 'I', uplo, n, kd, v, ldu, u, ldu, vl,
2259  $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
2260  $ iwork, iwork( 5*n+1 ), iinfo )
2261  IF( iinfo.NE.0 ) THEN
2262  WRITE( nounit, fmt = 9999 )'SSBEVX(V,I,' // uplo //
2263  $ ')', iinfo, n, jtype, ioldsd
2264  info = abs( iinfo )
2265  IF( iinfo.LT.0 ) THEN
2266  RETURN
2267  ELSE
2268  result( ntest ) = ulpinv
2269  result( ntest+1 ) = ulpinv
2270  result( ntest+2 ) = ulpinv
2271  GO TO 1370
2272  END IF
2273  END IF
2274 *
2275 * Do tests 55 and 56 (or +54)
2276 *
2277  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2278  $ v, ldu, tau, work, result( ntest ) )
2279 *
2280  ntest = ntest + 2
2281 *
2282  IF( iuplo.EQ.1 ) THEN
2283  DO 1340 j = 1, n
2284  DO 1330 i = max( 1, j-kd ), j
2285  v( kd+1+i-j, j ) = a( i, j )
2286  1330 CONTINUE
2287  1340 CONTINUE
2288  ELSE
2289  DO 1360 j = 1, n
2290  DO 1350 i = j, min( n, j+kd )
2291  v( 1+i-j, j ) = a( i, j )
2292  1350 CONTINUE
2293  1360 CONTINUE
2294  END IF
2295 *
2296  srnamt = 'SSBEVX_2STAGE'
2297  CALL ssbevx_2stage( 'N', 'I', uplo, n, kd, v, ldu,
2298  $ u, ldu, vl, vu, il, iu, abstol, m3, wa3,
2299  $ z, ldu, work, lwork, iwork, iwork( 5*n+1 ),
2300  $ iinfo )
2301  IF( iinfo.NE.0 ) THEN
2302  WRITE( nounit, fmt = 9999 )
2303  $ 'SSBEVX_2STAGE(N,I,' // uplo //
2304  $ ')', iinfo, n, jtype, ioldsd
2305  info = abs( iinfo )
2306  IF( iinfo.LT.0 ) THEN
2307  RETURN
2308  ELSE
2309  result( ntest ) = ulpinv
2310  GO TO 1370
2311  END IF
2312  END IF
2313 *
2314 * Do test 57 (or +54)
2315 *
2316  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2317  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2318  IF( n.GT.0 ) THEN
2319  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2320  ELSE
2321  temp3 = zero
2322  END IF
2323  result( ntest ) = ( temp1+temp2 ) /
2324  $ max( unfl, temp3*ulp )
2325 *
2326  1370 CONTINUE
2327  ntest = ntest + 1
2328  IF( iuplo.EQ.1 ) THEN
2329  DO 1390 j = 1, n
2330  DO 1380 i = max( 1, j-kd ), j
2331  v( kd+1+i-j, j ) = a( i, j )
2332  1380 CONTINUE
2333  1390 CONTINUE
2334  ELSE
2335  DO 1410 j = 1, n
2336  DO 1400 i = j, min( n, j+kd )
2337  v( 1+i-j, j ) = a( i, j )
2338  1400 CONTINUE
2339  1410 CONTINUE
2340  END IF
2341 *
2342  srnamt = 'SSBEVX'
2343  CALL ssbevx( 'V', 'V', uplo, n, kd, v, ldu, u, ldu, vl,
2344  $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
2345  $ iwork, iwork( 5*n+1 ), iinfo )
2346  IF( iinfo.NE.0 ) THEN
2347  WRITE( nounit, fmt = 9999 )'SSBEVX(V,V,' // uplo //
2348  $ ')', iinfo, n, jtype, ioldsd
2349  info = abs( iinfo )
2350  IF( iinfo.LT.0 ) THEN
2351  RETURN
2352  ELSE
2353  result( ntest ) = ulpinv
2354  result( ntest+1 ) = ulpinv
2355  result( ntest+2 ) = ulpinv
2356  GO TO 1460
2357  END IF
2358  END IF
2359 *
2360 * Do tests 58 and 59 (or +54)
2361 *
2362  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2363  $ v, ldu, tau, work, result( ntest ) )
2364 *
2365  ntest = ntest + 2
2366 *
2367  IF( iuplo.EQ.1 ) THEN
2368  DO 1430 j = 1, n
2369  DO 1420 i = max( 1, j-kd ), j
2370  v( kd+1+i-j, j ) = a( i, j )
2371  1420 CONTINUE
2372  1430 CONTINUE
2373  ELSE
2374  DO 1450 j = 1, n
2375  DO 1440 i = j, min( n, j+kd )
2376  v( 1+i-j, j ) = a( i, j )
2377  1440 CONTINUE
2378  1450 CONTINUE
2379  END IF
2380 *
2381  srnamt = 'SSBEVX_2STAGE'
2382  CALL ssbevx_2stage( 'N', 'V', uplo, n, kd, v, ldu,
2383  $ u, ldu, vl, vu, il, iu, abstol, m3, wa3,
2384  $ z, ldu, work, lwork, iwork, iwork( 5*n+1 ),
2385  $ iinfo )
2386  IF( iinfo.NE.0 ) THEN
2387  WRITE( nounit, fmt = 9999 )
2388  $ 'SSBEVX_2STAGE(N,V,' // uplo //
2389  $ ')', iinfo, n, jtype, ioldsd
2390  info = abs( iinfo )
2391  IF( iinfo.LT.0 ) THEN
2392  RETURN
2393  ELSE
2394  result( ntest ) = ulpinv
2395  GO TO 1460
2396  END IF
2397  END IF
2398 *
2399  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2400  result( ntest ) = ulpinv
2401  GO TO 1460
2402  END IF
2403 *
2404 * Do test 60 (or +54)
2405 *
2406  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2407  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2408  IF( n.GT.0 ) THEN
2409  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2410  ELSE
2411  temp3 = zero
2412  END IF
2413  result( ntest ) = ( temp1+temp2 ) /
2414  $ max( unfl, temp3*ulp )
2415 *
2416  1460 CONTINUE
2417 *
2418 * 7) Call SSYEVD
2419 *
2420  CALL slacpy( ' ', n, n, a, lda, v, ldu )
2421 *
2422  ntest = ntest + 1
2423  srnamt = 'SSYEVD'
2424  CALL ssyevd( 'V', uplo, n, a, ldu, d1, work, lwedc,
2425  $ iwork, liwedc, iinfo )
2426  IF( iinfo.NE.0 ) THEN
2427  WRITE( nounit, fmt = 9999 )'SSYEVD(V,' // uplo //
2428  $ ')', iinfo, n, jtype, ioldsd
2429  info = abs( iinfo )
2430  IF( iinfo.LT.0 ) THEN
2431  RETURN
2432  ELSE
2433  result( ntest ) = ulpinv
2434  result( ntest+1 ) = ulpinv
2435  result( ntest+2 ) = ulpinv
2436  GO TO 1480
2437  END IF
2438  END IF
2439 *
2440 * Do tests 61 and 62 (or +54)
2441 *
2442  CALL ssyt21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
2443  $ ldu, tau, work, result( ntest ) )
2444 *
2445  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2446 *
2447  ntest = ntest + 2
2448  srnamt = 'SSYEVD_2STAGE'
2449  CALL ssyevd_2stage( 'N', uplo, n, a, ldu, d3, work,
2450  $ lwork, iwork, liwedc, iinfo )
2451  IF( iinfo.NE.0 ) THEN
2452  WRITE( nounit, fmt = 9999 )
2453  $ 'SSYEVD_2STAGE(N,' // uplo //
2454  $ ')', iinfo, n, jtype, ioldsd
2455  info = abs( iinfo )
2456  IF( iinfo.LT.0 ) THEN
2457  RETURN
2458  ELSE
2459  result( ntest ) = ulpinv
2460  GO TO 1480
2461  END IF
2462  END IF
2463 *
2464 * Do test 63 (or +54)
2465 *
2466  temp1 = zero
2467  temp2 = zero
2468  DO 1470 j = 1, n
2469  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2470  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2471  1470 CONTINUE
2472  result( ntest ) = temp2 / max( unfl,
2473  $ ulp*max( temp1, temp2 ) )
2474 *
2475  1480 CONTINUE
2476 *
2477 * 8) Call SSPEVD.
2478 *
2479  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2480 *
2481 * Load array WORK with the upper or lower triangular
2482 * part of the matrix in packed form.
2483 *
2484  IF( iuplo.EQ.1 ) THEN
2485  indx = 1
2486  DO 1500 j = 1, n
2487  DO 1490 i = 1, j
2488  work( indx ) = a( i, j )
2489  indx = indx + 1
2490  1490 CONTINUE
2491  1500 CONTINUE
2492  ELSE
2493  indx = 1
2494  DO 1520 j = 1, n
2495  DO 1510 i = j, n
2496  work( indx ) = a( i, j )
2497  indx = indx + 1
2498  1510 CONTINUE
2499  1520 CONTINUE
2500  END IF
2501 *
2502  ntest = ntest + 1
2503  srnamt = 'SSPEVD'
2504  CALL sspevd( 'V', uplo, n, work, d1, z, ldu,
2505  $ work( indx ), lwedc-indx+1, iwork, liwedc,
2506  $ iinfo )
2507  IF( iinfo.NE.0 ) THEN
2508  WRITE( nounit, fmt = 9999 )'SSPEVD(V,' // uplo //
2509  $ ')', iinfo, n, jtype, ioldsd
2510  info = abs( iinfo )
2511  IF( iinfo.LT.0 ) THEN
2512  RETURN
2513  ELSE
2514  result( ntest ) = ulpinv
2515  result( ntest+1 ) = ulpinv
2516  result( ntest+2 ) = ulpinv
2517  GO TO 1580
2518  END IF
2519  END IF
2520 *
2521 * Do tests 64 and 65 (or +54)
2522 *
2523  CALL ssyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
2524  $ ldu, tau, work, result( ntest ) )
2525 *
2526  IF( iuplo.EQ.1 ) THEN
2527  indx = 1
2528  DO 1540 j = 1, n
2529  DO 1530 i = 1, j
2530 *
2531  work( indx ) = a( i, j )
2532  indx = indx + 1
2533  1530 CONTINUE
2534  1540 CONTINUE
2535  ELSE
2536  indx = 1
2537  DO 1560 j = 1, n
2538  DO 1550 i = j, n
2539  work( indx ) = a( i, j )
2540  indx = indx + 1
2541  1550 CONTINUE
2542  1560 CONTINUE
2543  END IF
2544 *
2545  ntest = ntest + 2
2546  srnamt = 'SSPEVD'
2547  CALL sspevd( 'N', uplo, n, work, d3, z, ldu,
2548  $ work( indx ), lwedc-indx+1, iwork, liwedc,
2549  $ iinfo )
2550  IF( iinfo.NE.0 ) THEN
2551  WRITE( nounit, fmt = 9999 )'SSPEVD(N,' // uplo //
2552  $ ')', iinfo, n, jtype, ioldsd
2553  info = abs( iinfo )
2554  IF( iinfo.LT.0 ) THEN
2555  RETURN
2556  ELSE
2557  result( ntest ) = ulpinv
2558  GO TO 1580
2559  END IF
2560  END IF
2561 *
2562 * Do test 66 (or +54)
2563 *
2564  temp1 = zero
2565  temp2 = zero
2566  DO 1570 j = 1, n
2567  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2568  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2569  1570 CONTINUE
2570  result( ntest ) = temp2 / max( unfl,
2571  $ ulp*max( temp1, temp2 ) )
2572  1580 CONTINUE
2573 *
2574 * 9) Call SSBEVD.
2575 *
2576  IF( jtype.LE.7 ) THEN
2577  kd = 1
2578  ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
2579  kd = max( n-1, 0 )
2580  ELSE
2581  kd = ihbw
2582  END IF
2583 *
2584 * Load array V with the upper or lower triangular part
2585 * of the matrix in band form.
2586 *
2587  IF( iuplo.EQ.1 ) THEN
2588  DO 1600 j = 1, n
2589  DO 1590 i = max( 1, j-kd ), j
2590  v( kd+1+i-j, j ) = a( i, j )
2591  1590 CONTINUE
2592  1600 CONTINUE
2593  ELSE
2594  DO 1620 j = 1, n
2595  DO 1610 i = j, min( n, j+kd )
2596  v( 1+i-j, j ) = a( i, j )
2597  1610 CONTINUE
2598  1620 CONTINUE
2599  END IF
2600 *
2601  ntest = ntest + 1
2602  srnamt = 'SSBEVD'
2603  CALL ssbevd( 'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
2604  $ lwedc, iwork, liwedc, iinfo )
2605  IF( iinfo.NE.0 ) THEN
2606  WRITE( nounit, fmt = 9999 )'SSBEVD(V,' // uplo //
2607  $ ')', iinfo, n, jtype, ioldsd
2608  info = abs( iinfo )
2609  IF( iinfo.LT.0 ) THEN
2610  RETURN
2611  ELSE
2612  result( ntest ) = ulpinv
2613  result( ntest+1 ) = ulpinv
2614  result( ntest+2 ) = ulpinv
2615  GO TO 1680
2616  END IF
2617  END IF
2618 *
2619 * Do tests 67 and 68 (or +54)
2620 *
2621  CALL ssyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
2622  $ ldu, tau, work, result( ntest ) )
2623 *
2624  IF( iuplo.EQ.1 ) THEN
2625  DO 1640 j = 1, n
2626  DO 1630 i = max( 1, j-kd ), j
2627  v( kd+1+i-j, j ) = a( i, j )
2628  1630 CONTINUE
2629  1640 CONTINUE
2630  ELSE
2631  DO 1660 j = 1, n
2632  DO 1650 i = j, min( n, j+kd )
2633  v( 1+i-j, j ) = a( i, j )
2634  1650 CONTINUE
2635  1660 CONTINUE
2636  END IF
2637 *
2638  ntest = ntest + 2
2639  srnamt = 'SSBEVD_2STAGE'
2640  CALL ssbevd_2stage( 'N', uplo, n, kd, v, ldu, d3, z, ldu,
2641  $ work, lwork, iwork, liwedc, iinfo )
2642  IF( iinfo.NE.0 ) THEN
2643  WRITE( nounit, fmt = 9999 )
2644  $ 'SSBEVD_2STAGE(N,' // uplo //
2645  $ ')', iinfo, n, jtype, ioldsd
2646  info = abs( iinfo )
2647  IF( iinfo.LT.0 ) THEN
2648  RETURN
2649  ELSE
2650  result( ntest ) = ulpinv
2651  GO TO 1680
2652  END IF
2653  END IF
2654 *
2655 * Do test 69 (or +54)
2656 *
2657  temp1 = zero
2658  temp2 = zero
2659  DO 1670 j = 1, n
2660  temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2661  temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2662  1670 CONTINUE
2663  result( ntest ) = temp2 / max( unfl,
2664  $ ulp*max( temp1, temp2 ) )
2665 *
2666  1680 CONTINUE
2667 *
2668 *
2669  CALL slacpy( ' ', n, n, a, lda, v, ldu )
2670  ntest = ntest + 1
2671  srnamt = 'SSYEVR'
2672  CALL ssyevr( 'V', 'A', uplo, n, a, ldu, vl, vu, il, iu,
2673  $ abstol, m, wa1, z, ldu, iwork, work, lwork,
2674  $ iwork(2*n+1), liwork-2*n, iinfo )
2675  IF( iinfo.NE.0 ) THEN
2676  WRITE( nounit, fmt = 9999 )'SSYEVR(V,A,' // uplo //
2677  $ ')', iinfo, n, jtype, ioldsd
2678  info = abs( iinfo )
2679  IF( iinfo.LT.0 ) THEN
2680  RETURN
2681  ELSE
2682  result( ntest ) = ulpinv
2683  result( ntest+1 ) = ulpinv
2684  result( ntest+2 ) = ulpinv
2685  GO TO 1700
2686  END IF
2687  END IF
2688 *
2689 * Do tests 70 and 71 (or ... )
2690 *
2691  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2692 *
2693  CALL ssyt21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
2694  $ ldu, tau, work, result( ntest ) )
2695 *
2696  ntest = ntest + 2
2697  srnamt = 'SSYEVR_2STAGE'
2698  CALL ssyevr_2stage( 'N', 'A', uplo, n, a, ldu, vl, vu,
2699  $ il, iu, abstol, m2, wa2, z, ldu, iwork,
2700  $ work, lwork, iwork(2*n+1), liwork-2*n,
2701  $ iinfo )
2702  IF( iinfo.NE.0 ) THEN
2703  WRITE( nounit, fmt = 9999 )
2704  $ 'SSYEVR_2STAGE(N,A,' // uplo //
2705  $ ')', iinfo, n, jtype, ioldsd
2706  info = abs( iinfo )
2707  IF( iinfo.LT.0 ) THEN
2708  RETURN
2709  ELSE
2710  result( ntest ) = ulpinv
2711  GO TO 1700
2712  END IF
2713  END IF
2714 *
2715 * Do test 72 (or ... )
2716 *
2717  temp1 = zero
2718  temp2 = zero
2719  DO 1690 j = 1, n
2720  temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
2721  temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
2722  1690 CONTINUE
2723  result( ntest ) = temp2 / max( unfl,
2724  $ ulp*max( temp1, temp2 ) )
2725 *
2726  1700 CONTINUE
2727 *
2728  ntest = ntest + 1
2729  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2730  srnamt = 'SSYEVR'
2731  CALL ssyevr( 'V', 'I', uplo, n, a, ldu, vl, vu, il, iu,
2732  $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
2733  $ iwork(2*n+1), liwork-2*n, iinfo )
2734  IF( iinfo.NE.0 ) THEN
2735  WRITE( nounit, fmt = 9999 )'SSYEVR(V,I,' // uplo //
2736  $ ')', iinfo, n, jtype, ioldsd
2737  info = abs( iinfo )
2738  IF( iinfo.LT.0 ) THEN
2739  RETURN
2740  ELSE
2741  result( ntest ) = ulpinv
2742  result( ntest+1 ) = ulpinv
2743  result( ntest+2 ) = ulpinv
2744  GO TO 1710
2745  END IF
2746  END IF
2747 *
2748 * Do tests 73 and 74 (or +54)
2749 *
2750  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2751 *
2752  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2753  $ v, ldu, tau, work, result( ntest ) )
2754 *
2755  ntest = ntest + 2
2756  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2757  srnamt = 'SSYEVR_2STAGE'
2758  CALL ssyevr_2stage( 'N', 'I', uplo, n, a, ldu, vl, vu,
2759  $ il, iu, abstol, m3, wa3, z, ldu, iwork,
2760  $ work, lwork, iwork(2*n+1), liwork-2*n,
2761  $ iinfo )
2762  IF( iinfo.NE.0 ) THEN
2763  WRITE( nounit, fmt = 9999 )
2764  $ 'SSYEVR_2STAGE(N,I,' // uplo //
2765  $ ')', iinfo, n, jtype, ioldsd
2766  info = abs( iinfo )
2767  IF( iinfo.LT.0 ) THEN
2768  RETURN
2769  ELSE
2770  result( ntest ) = ulpinv
2771  GO TO 1710
2772  END IF
2773  END IF
2774 *
2775 * Do test 75 (or +54)
2776 *
2777  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2778  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2779  result( ntest ) = ( temp1+temp2 ) /
2780  $ max( unfl, ulp*temp3 )
2781  1710 CONTINUE
2782 *
2783  ntest = ntest + 1
2784  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2785  srnamt = 'SSYEVR'
2786  CALL ssyevr( 'V', 'V', uplo, n, a, ldu, vl, vu, il, iu,
2787  $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
2788  $ iwork(2*n+1), liwork-2*n, iinfo )
2789  IF( iinfo.NE.0 ) THEN
2790  WRITE( nounit, fmt = 9999 )'SSYEVR(V,V,' // uplo //
2791  $ ')', iinfo, n, jtype, ioldsd
2792  info = abs( iinfo )
2793  IF( iinfo.LT.0 ) THEN
2794  RETURN
2795  ELSE
2796  result( ntest ) = ulpinv
2797  result( ntest+1 ) = ulpinv
2798  result( ntest+2 ) = ulpinv
2799  GO TO 700
2800  END IF
2801  END IF
2802 *
2803 * Do tests 76 and 77 (or +54)
2804 *
2805  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2806 *
2807  CALL ssyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2808  $ v, ldu, tau, work, result( ntest ) )
2809 *
2810  ntest = ntest + 2
2811  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2812  srnamt = 'SSYEVR_2STAGE'
2813  CALL ssyevr_2stage( 'N', 'V', uplo, n, a, ldu, vl, vu,
2814  $ il, iu, abstol, m3, wa3, z, ldu, iwork,
2815  $ work, lwork, iwork(2*n+1), liwork-2*n,
2816  $ iinfo )
2817  IF( iinfo.NE.0 ) THEN
2818  WRITE( nounit, fmt = 9999 )
2819  $ 'SSYEVR_2STAGE(N,V,' // uplo //
2820  $ ')', iinfo, n, jtype, ioldsd
2821  info = abs( iinfo )
2822  IF( iinfo.LT.0 ) THEN
2823  RETURN
2824  ELSE
2825  result( ntest ) = ulpinv
2826  GO TO 700
2827  END IF
2828  END IF
2829 *
2830  IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2831  result( ntest ) = ulpinv
2832  GO TO 700
2833  END IF
2834 *
2835 * Do test 78 (or +54)
2836 *
2837  temp1 = ssxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2838  temp2 = ssxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2839  IF( n.GT.0 ) THEN
2840  temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2841  ELSE
2842  temp3 = zero
2843  END IF
2844  result( ntest ) = ( temp1+temp2 ) /
2845  $ max( unfl, temp3*ulp )
2846 *
2847  CALL slacpy( ' ', n, n, v, ldu, a, lda )
2848 *
2849  1720 CONTINUE
2850 *
2851 * End of Loop -- Check for RESULT(j) > THRESH
2852 *
2853  ntestt = ntestt + ntest
2854 *
2855  CALL slafts( 'SST', n, n, jtype, ntest, result, ioldsd,
2856  $ thresh, nounit, nerrs )
2857 *
2858  1730 CONTINUE
2859  1740 CONTINUE
2860 *
2861 * Summary
2862 *
2863  CALL alasvm( 'SST', nounit, nerrs, ntestt, 0 )
2864 *
2865  9999 FORMAT( ' SDRVST2STG: ', a, ' returned INFO=', i6, '.', / 9x,
2866  $ 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
2867 *
2868  RETURN
2869 *
2870 * End of SDRVST2STG
2871 *
2872  END
subroutine ssyevd_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO)
SSYEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY ma...
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine ssyevr(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSYEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: ssyevr.f:338
subroutine sspevd(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: sspevd.f:180
subroutine sspev(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, INFO)
SSPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: sspev.f:132
subroutine ssyev_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO)
SSYEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...
Definition: ssyev_2stage.f:185
subroutine sstt22(N, M, KBAND, AD, AE, SD, SE, U, LDU, WORK, LDWORK, RESULT)
SSTT22
Definition: sstt22.f:141
subroutine sdrvst2stg(NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH, NOUNIT, A, LDA, D1, D2, D3, D4, EVEIGS, WA1, WA2, WA3, U, LDU, V, TAU, Z, WORK, LWORK, IWORK, LIWORK, RESULT, INFO)
SDRVST2STG
Definition: sdrvst2stg.f:455
subroutine sstevr(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSTEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: sstevr.f:308
subroutine ssyevx(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO)
SSYEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: ssyevx.f:255
subroutine ssbevd(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSBEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: ssbevd.f:195
subroutine sstevx(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
SSTEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: sstevx.f:229
subroutine slatmr(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER, CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM, PACK, A, LDA, IWORK, INFO)
SLATMR
Definition: slatmr.f:473
subroutine slafts(TYPE, M, N, IMAT, NTESTS, RESULT, ISEED, THRESH, IOUNIT, IE)
SLAFTS
Definition: slafts.f:101
subroutine ssyevr_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSYEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY ma...
subroutine sstev(JOBZ, N, D, E, Z, LDZ, WORK, INFO)
SSTEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: sstev.f:118
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine ssytrd_sy2sb(UPLO, N, KD, A, LDA, AB, LDAB, TAU, WORK, LWORK, INFO)
SSYTRD_SY2SB
Definition: ssytrd_sy2sb.f:245
subroutine sspevx(JOBZ, RANGE, UPLO, N, AP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
SSPEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: sspevx.f:236
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
subroutine ssyevx_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO)
SSYEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY ma...
subroutine ssyevd(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO)
SSYEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: ssyevd.f:185
subroutine ssyt21(ITYPE, UPLO, N, KBAND, A, LDA, D, E, U, LDU, V, LDV, TAU, WORK, RESULT)
SSYT21
Definition: ssyt21.f:207
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:76
subroutine ssytrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
SSYTRD_2STAGE
subroutine ssbevx(JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
SSBEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: ssbevx.f:267
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine ssyev(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO)
SSYEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: ssyev.f:134
subroutine sstevd(JOBZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSTEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: sstevd.f:165
subroutine ssbevd_2stage(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER...
subroutine ssyt22(ITYPE, UPLO, N, M, KBAND, A, LDA, D, E, U, LDU, V, LDV, TAU, WORK, RESULT)
SSYT22
Definition: ssyt22.f:157
subroutine ssbevx_2stage(JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO)
SSBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER...
subroutine ssbev(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, INFO)
SSBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: ssbev.f:148
subroutine ssbev_2stage(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, INFO)
SSBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
Definition: ssbev_2stage.f:206
subroutine sstt21(N, KBAND, AD, AE, SD, SE, U, LDU, WORK, RESULT)
SSTT21
Definition: sstt21.f:129