LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine derrst ( character*3  PATH,
integer  NUNIT 
)

DERRST

Purpose:
 DERRST tests the error exits for DSYTRD, DORGTR, DORMTR, DSPTRD,
 DOPGTR, DOPMTR, DSTEQR, SSTERF, SSTEBZ, SSTEIN, DPTEQR, DSBTRD,
 DSYEV, SSYEVX, SSYEVD, DSBEV, SSBEVX, SSBEVD,
 DSPEV, SSPEVX, SSPEVD, DSTEV, SSTEVX, SSTEVD, and SSTEDC.
Parameters
[in]PATH
          PATH is CHARACTER*3
          The LAPACK path name for the routines to be tested.
[in]NUNIT
          NUNIT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 59 of file derrst.f.

59 *
60 * -- LAPACK test routine (version 3.4.0) --
61 * -- LAPACK is a software package provided by Univ. of Tennessee, --
62 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
63 * November 2011
64 *
65 * .. Scalar Arguments ..
66  CHARACTER*3 path
67  INTEGER nunit
68 * ..
69 *
70 * =====================================================================
71 *
72 * NMAX has to be at least 3 or LIW may be too small
73 * .. Parameters ..
74  INTEGER nmax, liw, lw
75  parameter ( nmax = 3, liw = 12*nmax, lw = 20*nmax )
76 * ..
77 * .. Local Scalars ..
78  CHARACTER*2 c2
79  INTEGER i, info, j, m, n, nsplit, nt
80 * ..
81 * .. Local Arrays ..
82  INTEGER i1( nmax ), i2( nmax ), i3( nmax ), iw( liw )
83  DOUBLE PRECISION a( nmax, nmax ), c( nmax, nmax ), d( nmax ),
84  $ e( nmax ), q( nmax, nmax ), r( nmax ),
85  $ tau( nmax ), w( lw ), x( nmax ),
86  $ z( nmax, nmax )
87 * ..
88 * .. External Functions ..
89  LOGICAL lsamen
90  EXTERNAL lsamen
91 * ..
92 * .. External Subroutines ..
93  EXTERNAL chkxer, dopgtr, dopmtr, dorgtr, dormtr, dpteqr,
98 * ..
99 * .. Scalars in Common ..
100  LOGICAL lerr, ok
101  CHARACTER*32 srnamt
102  INTEGER infot, nout
103 * ..
104 * .. Common blocks ..
105  COMMON / infoc / infot, nout, ok, lerr
106  COMMON / srnamc / srnamt
107 * ..
108 * .. Intrinsic Functions ..
109  INTRINSIC dble
110 * ..
111 * .. Executable Statements ..
112 *
113  nout = nunit
114  WRITE( nout, fmt = * )
115  c2 = path( 2: 3 )
116 *
117 * Set the variables to innocuous values.
118 *
119  DO 20 j = 1, nmax
120  DO 10 i = 1, nmax
121  a( i, j ) = 1.d0 / dble( i+j )
122  10 CONTINUE
123  20 CONTINUE
124  DO 30 j = 1, nmax
125  d( j ) = dble( j )
126  e( j ) = 0.0d0
127  i1( j ) = j
128  i2( j ) = j
129  tau( j ) = 1.d0
130  30 CONTINUE
131  ok = .true.
132  nt = 0
133 *
134 * Test error exits for the ST path.
135 *
136  IF( lsamen( 2, c2, 'ST' ) ) THEN
137 *
138 * DSYTRD
139 *
140  srnamt = 'DSYTRD'
141  infot = 1
142  CALL dsytrd( '/', 0, a, 1, d, e, tau, w, 1, info )
143  CALL chkxer( 'DSYTRD', infot, nout, lerr, ok )
144  infot = 2
145  CALL dsytrd( 'U', -1, a, 1, d, e, tau, w, 1, info )
146  CALL chkxer( 'DSYTRD', infot, nout, lerr, ok )
147  infot = 4
148  CALL dsytrd( 'U', 2, a, 1, d, e, tau, w, 1, info )
149  CALL chkxer( 'DSYTRD', infot, nout, lerr, ok )
150  infot = 9
151  CALL dsytrd( 'U', 0, a, 1, d, e, tau, w, 0, info )
152  CALL chkxer( 'DSYTRD', infot, nout, lerr, ok )
153  nt = nt + 4
154 *
155 * DORGTR
156 *
157  srnamt = 'DORGTR'
158  infot = 1
159  CALL dorgtr( '/', 0, a, 1, tau, w, 1, info )
160  CALL chkxer( 'DORGTR', infot, nout, lerr, ok )
161  infot = 2
162  CALL dorgtr( 'U', -1, a, 1, tau, w, 1, info )
163  CALL chkxer( 'DORGTR', infot, nout, lerr, ok )
164  infot = 4
165  CALL dorgtr( 'U', 2, a, 1, tau, w, 1, info )
166  CALL chkxer( 'DORGTR', infot, nout, lerr, ok )
167  infot = 7
168  CALL dorgtr( 'U', 3, a, 3, tau, w, 1, info )
169  CALL chkxer( 'DORGTR', infot, nout, lerr, ok )
170  nt = nt + 4
171 *
172 * DORMTR
173 *
174  srnamt = 'DORMTR'
175  infot = 1
176  CALL dormtr( '/', 'U', 'N', 0, 0, a, 1, tau, c, 1, w, 1, info )
177  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
178  infot = 2
179  CALL dormtr( 'L', '/', 'N', 0, 0, a, 1, tau, c, 1, w, 1, info )
180  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
181  infot = 3
182  CALL dormtr( 'L', 'U', '/', 0, 0, a, 1, tau, c, 1, w, 1, info )
183  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
184  infot = 4
185  CALL dormtr( 'L', 'U', 'N', -1, 0, a, 1, tau, c, 1, w, 1,
186  $ info )
187  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
188  infot = 5
189  CALL dormtr( 'L', 'U', 'N', 0, -1, a, 1, tau, c, 1, w, 1,
190  $ info )
191  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
192  infot = 7
193  CALL dormtr( 'L', 'U', 'N', 2, 0, a, 1, tau, c, 2, w, 1, info )
194  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
195  infot = 7
196  CALL dormtr( 'R', 'U', 'N', 0, 2, a, 1, tau, c, 1, w, 1, info )
197  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
198  infot = 10
199  CALL dormtr( 'L', 'U', 'N', 2, 0, a, 2, tau, c, 1, w, 1, info )
200  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
201  infot = 12
202  CALL dormtr( 'L', 'U', 'N', 0, 2, a, 1, tau, c, 1, w, 1, info )
203  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
204  infot = 12
205  CALL dormtr( 'R', 'U', 'N', 2, 0, a, 1, tau, c, 2, w, 1, info )
206  CALL chkxer( 'DORMTR', infot, nout, lerr, ok )
207  nt = nt + 10
208 *
209 * DSPTRD
210 *
211  srnamt = 'DSPTRD'
212  infot = 1
213  CALL dsptrd( '/', 0, a, d, e, tau, info )
214  CALL chkxer( 'DSPTRD', infot, nout, lerr, ok )
215  infot = 2
216  CALL dsptrd( 'U', -1, a, d, e, tau, info )
217  CALL chkxer( 'DSPTRD', infot, nout, lerr, ok )
218  nt = nt + 2
219 *
220 * DOPGTR
221 *
222  srnamt = 'DOPGTR'
223  infot = 1
224  CALL dopgtr( '/', 0, a, tau, z, 1, w, info )
225  CALL chkxer( 'DOPGTR', infot, nout, lerr, ok )
226  infot = 2
227  CALL dopgtr( 'U', -1, a, tau, z, 1, w, info )
228  CALL chkxer( 'DOPGTR', infot, nout, lerr, ok )
229  infot = 6
230  CALL dopgtr( 'U', 2, a, tau, z, 1, w, info )
231  CALL chkxer( 'DOPGTR', infot, nout, lerr, ok )
232  nt = nt + 3
233 *
234 * DOPMTR
235 *
236  srnamt = 'DOPMTR'
237  infot = 1
238  CALL dopmtr( '/', 'U', 'N', 0, 0, a, tau, c, 1, w, info )
239  CALL chkxer( 'DOPMTR', infot, nout, lerr, ok )
240  infot = 2
241  CALL dopmtr( 'L', '/', 'N', 0, 0, a, tau, c, 1, w, info )
242  CALL chkxer( 'DOPMTR', infot, nout, lerr, ok )
243  infot = 3
244  CALL dopmtr( 'L', 'U', '/', 0, 0, a, tau, c, 1, w, info )
245  CALL chkxer( 'DOPMTR', infot, nout, lerr, ok )
246  infot = 4
247  CALL dopmtr( 'L', 'U', 'N', -1, 0, a, tau, c, 1, w, info )
248  CALL chkxer( 'DOPMTR', infot, nout, lerr, ok )
249  infot = 5
250  CALL dopmtr( 'L', 'U', 'N', 0, -1, a, tau, c, 1, w, info )
251  CALL chkxer( 'DOPMTR', infot, nout, lerr, ok )
252  infot = 9
253  CALL dopmtr( 'L', 'U', 'N', 2, 0, a, tau, c, 1, w, info )
254  CALL chkxer( 'DOPMTR', infot, nout, lerr, ok )
255  nt = nt + 6
256 *
257 * DPTEQR
258 *
259  srnamt = 'DPTEQR'
260  infot = 1
261  CALL dpteqr( '/', 0, d, e, z, 1, w, info )
262  CALL chkxer( 'DPTEQR', infot, nout, lerr, ok )
263  infot = 2
264  CALL dpteqr( 'N', -1, d, e, z, 1, w, info )
265  CALL chkxer( 'DPTEQR', infot, nout, lerr, ok )
266  infot = 6
267  CALL dpteqr( 'V', 2, d, e, z, 1, w, info )
268  CALL chkxer( 'DPTEQR', infot, nout, lerr, ok )
269  nt = nt + 3
270 *
271 * DSTEBZ
272 *
273  srnamt = 'DSTEBZ'
274  infot = 1
275  CALL dstebz( '/', 'E', 0, 0.0d0, 1.0d0, 1, 0, 0.0d0, d, e, m,
276  $ nsplit, x, i1, i2, w, iw, info )
277  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
278  infot = 2
279  CALL dstebz( 'A', '/', 0, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
280  $ nsplit, x, i1, i2, w, iw, info )
281  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
282  infot = 3
283  CALL dstebz( 'A', 'E', -1, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
284  $ nsplit, x, i1, i2, w, iw, info )
285  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
286  infot = 5
287  CALL dstebz( 'V', 'E', 0, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
288  $ nsplit, x, i1, i2, w, iw, info )
289  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
290  infot = 6
291  CALL dstebz( 'I', 'E', 0, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
292  $ nsplit, x, i1, i2, w, iw, info )
293  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
294  infot = 6
295  CALL dstebz( 'I', 'E', 1, 0.0d0, 0.0d0, 2, 1, 0.0d0, d, e, m,
296  $ nsplit, x, i1, i2, w, iw, info )
297  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
298  infot = 7
299  CALL dstebz( 'I', 'E', 1, 0.0d0, 0.0d0, 1, 0, 0.0d0, d, e, m,
300  $ nsplit, x, i1, i2, w, iw, info )
301  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
302  infot = 7
303  CALL dstebz( 'I', 'E', 1, 0.0d0, 0.0d0, 1, 2, 0.0d0, d, e, m,
304  $ nsplit, x, i1, i2, w, iw, info )
305  CALL chkxer( 'DSTEBZ', infot, nout, lerr, ok )
306  nt = nt + 8
307 *
308 * DSTEIN
309 *
310  srnamt = 'DSTEIN'
311  infot = 1
312  CALL dstein( -1, d, e, 0, x, i1, i2, z, 1, w, iw, i3, info )
313  CALL chkxer( 'DSTEIN', infot, nout, lerr, ok )
314  infot = 4
315  CALL dstein( 0, d, e, -1, x, i1, i2, z, 1, w, iw, i3, info )
316  CALL chkxer( 'DSTEIN', infot, nout, lerr, ok )
317  infot = 4
318  CALL dstein( 0, d, e, 1, x, i1, i2, z, 1, w, iw, i3, info )
319  CALL chkxer( 'DSTEIN', infot, nout, lerr, ok )
320  infot = 9
321  CALL dstein( 2, d, e, 0, x, i1, i2, z, 1, w, iw, i3, info )
322  CALL chkxer( 'DSTEIN', infot, nout, lerr, ok )
323  nt = nt + 4
324 *
325 * DSTEQR
326 *
327  srnamt = 'DSTEQR'
328  infot = 1
329  CALL dsteqr( '/', 0, d, e, z, 1, w, info )
330  CALL chkxer( 'DSTEQR', infot, nout, lerr, ok )
331  infot = 2
332  CALL dsteqr( 'N', -1, d, e, z, 1, w, info )
333  CALL chkxer( 'DSTEQR', infot, nout, lerr, ok )
334  infot = 6
335  CALL dsteqr( 'V', 2, d, e, z, 1, w, info )
336  CALL chkxer( 'DSTEQR', infot, nout, lerr, ok )
337  nt = nt + 3
338 *
339 * DSTERF
340 *
341  srnamt = 'DSTERF'
342  infot = 1
343  CALL dsterf( -1, d, e, info )
344  CALL chkxer( 'DSTERF', infot, nout, lerr, ok )
345  nt = nt + 1
346 *
347 * DSTEDC
348 *
349  srnamt = 'DSTEDC'
350  infot = 1
351  CALL dstedc( '/', 0, d, e, z, 1, w, 1, iw, 1, info )
352  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
353  infot = 2
354  CALL dstedc( 'N', -1, d, e, z, 1, w, 1, iw, 1, info )
355  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
356  infot = 6
357  CALL dstedc( 'V', 2, d, e, z, 1, w, 23, iw, 28, info )
358  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
359  infot = 8
360  CALL dstedc( 'N', 1, d, e, z, 1, w, 0, iw, 1, info )
361  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
362  infot = 8
363  CALL dstedc( 'I', 2, d, e, z, 2, w, 0, iw, 12, info )
364  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
365  infot = 8
366  CALL dstedc( 'V', 2, d, e, z, 2, w, 0, iw, 28, info )
367  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
368  infot = 10
369  CALL dstedc( 'N', 1, d, e, z, 1, w, 1, iw, 0, info )
370  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
371  infot = 10
372  CALL dstedc( 'I', 2, d, e, z, 2, w, 19, iw, 0, info )
373  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
374  infot = 10
375  CALL dstedc( 'V', 2, d, e, z, 2, w, 23, iw, 0, info )
376  CALL chkxer( 'DSTEDC', infot, nout, lerr, ok )
377  nt = nt + 9
378 *
379 * DSTEVD
380 *
381  srnamt = 'DSTEVD'
382  infot = 1
383  CALL dstevd( '/', 0, d, e, z, 1, w, 1, iw, 1, info )
384  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
385  infot = 2
386  CALL dstevd( 'N', -1, d, e, z, 1, w, 1, iw, 1, info )
387  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
388  infot = 6
389  CALL dstevd( 'V', 2, d, e, z, 1, w, 19, iw, 12, info )
390  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
391  infot = 8
392  CALL dstevd( 'N', 1, d, e, z, 1, w, 0, iw, 1, info )
393  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
394  infot = 8
395  CALL dstevd( 'V', 2, d, e, z, 2, w, 12, iw, 12, info )
396  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
397  infot = 10
398  CALL dstevd( 'N', 0, d, e, z, 1, w, 1, iw, 0, info )
399  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
400  infot = 10
401  CALL dstevd( 'V', 2, d, e, z, 2, w, 19, iw, 11, info )
402  CALL chkxer( 'DSTEVD', infot, nout, lerr, ok )
403  nt = nt + 7
404 *
405 * DSTEV
406 *
407  srnamt = 'DSTEV '
408  infot = 1
409  CALL dstev( '/', 0, d, e, z, 1, w, info )
410  CALL chkxer( 'DSTEV ', infot, nout, lerr, ok )
411  infot = 2
412  CALL dstev( 'N', -1, d, e, z, 1, w, info )
413  CALL chkxer( 'DSTEV ', infot, nout, lerr, ok )
414  infot = 6
415  CALL dstev( 'V', 2, d, e, z, 1, w, info )
416  CALL chkxer( 'DSTEV ', infot, nout, lerr, ok )
417  nt = nt + 3
418 *
419 * DSTEVX
420 *
421  srnamt = 'DSTEVX'
422  infot = 1
423  CALL dstevx( '/', 'A', 0, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
424  $ x, z, 1, w, iw, i3, info )
425  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
426  infot = 2
427  CALL dstevx( 'N', '/', 0, d, e, 0.0d0, 1.0d0, 1, 0, 0.0d0, m,
428  $ x, z, 1, w, iw, i3, info )
429  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
430  infot = 3
431  CALL dstevx( 'N', 'A', -1, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
432  $ x, z, 1, w, iw, i3, info )
433  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
434  infot = 7
435  CALL dstevx( 'N', 'V', 1, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
436  $ x, z, 1, w, iw, i3, info )
437  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
438  infot = 8
439  CALL dstevx( 'N', 'I', 1, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
440  $ x, z, 1, w, iw, i3, info )
441  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
442  infot = 8
443  CALL dstevx( 'N', 'I', 1, d, e, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
444  $ x, z, 1, w, iw, i3, info )
445  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
446  infot = 9
447  CALL dstevx( 'N', 'I', 2, d, e, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
448  $ x, z, 1, w, iw, i3, info )
449  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
450  infot = 9
451  CALL dstevx( 'N', 'I', 1, d, e, 0.0d0, 0.0d0, 1, 2, 0.0d0, m,
452  $ x, z, 1, w, iw, i3, info )
453  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
454  infot = 14
455  CALL dstevx( 'V', 'A', 2, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
456  $ x, z, 1, w, iw, i3, info )
457  CALL chkxer( 'DSTEVX', infot, nout, lerr, ok )
458  nt = nt + 9
459 *
460 * DSTEVR
461 *
462  n = 1
463  srnamt = 'DSTEVR'
464  infot = 1
465  CALL dstevr( '/', 'A', 0, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
466  $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
467  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
468  infot = 2
469  CALL dstevr( 'V', '/', 0, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
470  $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
471  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
472  infot = 3
473  CALL dstevr( 'V', 'A', -1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
474  $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
475  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
476  infot = 7
477  CALL dstevr( 'V', 'V', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
478  $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
479  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
480  infot = 8
481  CALL dstevr( 'V', 'I', 1, d, e, 0.0d0, 0.0d0, 0, 1, 0.0d0, m,
482  $ w, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
483  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
484  infot = 9
485  n = 2
486  CALL dstevr( 'V', 'I', 2, d, e, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
487  $ w, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
488  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
489  infot = 14
490  n = 1
491  CALL dstevr( 'V', 'I', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
492  $ w, z, 0, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
493  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
494  infot = 17
495  CALL dstevr( 'V', 'I', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
496  $ w, z, 1, iw, x, 20*n-1, iw( 2*n+1 ), 10*n, info )
497  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
498  infot = 19
499  CALL dstevr( 'V', 'I', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
500  $ w, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n-1, info )
501  CALL chkxer( 'DSTEVR', infot, nout, lerr, ok )
502  nt = nt + 9
503 *
504 * DSYEVD
505 *
506  srnamt = 'DSYEVD'
507  infot = 1
508  CALL dsyevd( '/', 'U', 0, a, 1, x, w, 1, iw, 1, info )
509  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
510  infot = 2
511  CALL dsyevd( 'N', '/', 0, a, 1, x, w, 1, iw, 1, info )
512  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
513  infot = 3
514  CALL dsyevd( 'N', 'U', -1, a, 1, x, w, 1, iw, 1, info )
515  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
516  infot = 5
517  CALL dsyevd( 'N', 'U', 2, a, 1, x, w, 3, iw, 1, info )
518  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
519  infot = 8
520  CALL dsyevd( 'N', 'U', 1, a, 1, x, w, 0, iw, 1, info )
521  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
522  infot = 8
523  CALL dsyevd( 'N', 'U', 2, a, 2, x, w, 4, iw, 1, info )
524  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
525  infot = 8
526  CALL dsyevd( 'V', 'U', 2, a, 2, x, w, 20, iw, 12, info )
527  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
528  infot = 10
529  CALL dsyevd( 'N', 'U', 1, a, 1, x, w, 1, iw, 0, info )
530  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
531  infot = 10
532  CALL dsyevd( 'N', 'U', 2, a, 2, x, w, 5, iw, 0, info )
533  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
534  infot = 10
535  CALL dsyevd( 'V', 'U', 2, a, 2, x, w, 27, iw, 11, info )
536  CALL chkxer( 'DSYEVD', infot, nout, lerr, ok )
537  nt = nt + 10
538 *
539 * DSYEVR
540 *
541  srnamt = 'DSYEVR'
542  n = 1
543  infot = 1
544  CALL dsyevr( '/', 'A', 'U', 0, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
545  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
546  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
547  infot = 2
548  CALL dsyevr( 'V', '/', 'U', 0, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
549  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
550  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
551  infot = 3
552  CALL dsyevr( 'V', 'A', '/', -1, a, 1, 0.0d0, 0.0d0, 1, 1,
553  $ 0.0d0, m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n,
554  $ info )
555  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
556  infot = 4
557  CALL dsyevr( 'V', 'A', 'U', -1, a, 1, 0.0d0, 0.0d0, 1, 1,
558  $ 0.0d0, m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n,
559  $ info )
560  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
561  infot = 6
562  CALL dsyevr( 'V', 'A', 'U', 2, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
563  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
564  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
565  infot = 8
566  CALL dsyevr( 'V', 'V', 'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
567  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
568  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
569  infot = 9
570  CALL dsyevr( 'V', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 0, 1, 0.0d0,
571  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
572  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
573  infot = 10
574 *
575  CALL dsyevr( 'V', 'I', 'U', 2, a, 2, 0.0d0, 0.0d0, 2, 1, 0.0d0,
576  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
577  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
578  infot = 15
579  CALL dsyevr( 'V', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
580  $ m, r, z, 0, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
581  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
582  infot = 18
583  CALL dsyevr( 'V', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
584  $ m, r, z, 1, iw, q, 26*n-1, iw( 2*n+1 ), 10*n,
585  $ info )
586  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
587  infot = 20
588  CALL dsyevr( 'V', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
589  $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n-1,
590  $ info )
591  CALL chkxer( 'DSYEVR', infot, nout, lerr, ok )
592  nt = nt + 11
593 *
594 * DSYEV
595 *
596  srnamt = 'DSYEV '
597  infot = 1
598  CALL dsyev( '/', 'U', 0, a, 1, x, w, 1, info )
599  CALL chkxer( 'DSYEV ', infot, nout, lerr, ok )
600  infot = 2
601  CALL dsyev( 'N', '/', 0, a, 1, x, w, 1, info )
602  CALL chkxer( 'DSYEV ', infot, nout, lerr, ok )
603  infot = 3
604  CALL dsyev( 'N', 'U', -1, a, 1, x, w, 1, info )
605  CALL chkxer( 'DSYEV ', infot, nout, lerr, ok )
606  infot = 5
607  CALL dsyev( 'N', 'U', 2, a, 1, x, w, 3, info )
608  CALL chkxer( 'DSYEV ', infot, nout, lerr, ok )
609  infot = 8
610  CALL dsyev( 'N', 'U', 1, a, 1, x, w, 1, info )
611  CALL chkxer( 'DSYEV ', infot, nout, lerr, ok )
612  nt = nt + 5
613 *
614 * DSYEVX
615 *
616  srnamt = 'DSYEVX'
617  infot = 1
618  CALL dsyevx( '/', 'A', 'U', 0, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
619  $ m, x, z, 1, w, 1, iw, i3, info )
620  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
621  infot = 2
622  CALL dsyevx( 'N', '/', 'U', 0, a, 1, 0.0d0, 1.0d0, 1, 0, 0.0d0,
623  $ m, x, z, 1, w, 1, iw, i3, info )
624  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
625  infot = 3
626  CALL dsyevx( 'N', 'A', '/', 0, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
627  $ m, x, z, 1, w, 1, iw, i3, info )
628  infot = 4
629  CALL dsyevx( 'N', 'A', 'U', -1, a, 1, 0.0d0, 0.0d0, 0, 0,
630  $ 0.0d0, m, x, z, 1, w, 1, iw, i3, info )
631  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
632  infot = 6
633  CALL dsyevx( 'N', 'A', 'U', 2, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
634  $ m, x, z, 1, w, 16, iw, i3, info )
635  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
636  infot = 8
637  CALL dsyevx( 'N', 'V', 'U', 1, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
638  $ m, x, z, 1, w, 8, iw, i3, info )
639  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
640  infot = 9
641  CALL dsyevx( 'N', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
642  $ m, x, z, 1, w, 8, iw, i3, info )
643  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
644  infot = 9
645  CALL dsyevx( 'N', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 2, 1, 0.0d0,
646  $ m, x, z, 1, w, 8, iw, i3, info )
647  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
648  infot = 10
649  CALL dsyevx( 'N', 'I', 'U', 2, a, 2, 0.0d0, 0.0d0, 2, 1, 0.0d0,
650  $ m, x, z, 1, w, 16, iw, i3, info )
651  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
652  infot = 10
653  CALL dsyevx( 'N', 'I', 'U', 1, a, 1, 0.0d0, 0.0d0, 1, 2, 0.0d0,
654  $ m, x, z, 1, w, 8, iw, i3, info )
655  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
656  infot = 15
657  CALL dsyevx( 'V', 'A', 'U', 2, a, 2, 0.0d0, 0.0d0, 0, 0, 0.0d0,
658  $ m, x, z, 1, w, 16, iw, i3, info )
659  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
660  infot = 17
661  CALL dsyevx( 'V', 'A', 'U', 1, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
662  $ m, x, z, 1, w, 0, iw, i3, info )
663  CALL chkxer( 'DSYEVX', infot, nout, lerr, ok )
664  nt = nt + 12
665 *
666 * DSPEVD
667 *
668  srnamt = 'DSPEVD'
669  infot = 1
670  CALL dspevd( '/', 'U', 0, a, x, z, 1, w, 1, iw, 1, info )
671  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
672  infot = 2
673  CALL dspevd( 'N', '/', 0, a, x, z, 1, w, 1, iw, 1, info )
674  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
675  infot = 3
676  CALL dspevd( 'N', 'U', -1, a, x, z, 1, w, 1, iw, 1, info )
677  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
678  infot = 7
679  CALL dspevd( 'V', 'U', 2, a, x, z, 1, w, 23, iw, 12, info )
680  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
681  infot = 9
682  CALL dspevd( 'N', 'U', 1, a, x, z, 1, w, 0, iw, 1, info )
683  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
684  infot = 9
685  CALL dspevd( 'N', 'U', 2, a, x, z, 1, w, 3, iw, 1, info )
686  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
687  infot = 9
688  CALL dspevd( 'V', 'U', 2, a, x, z, 2, w, 16, iw, 12, info )
689  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
690  infot = 11
691  CALL dspevd( 'N', 'U', 1, a, x, z, 1, w, 1, iw, 0, info )
692  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
693  infot = 11
694  CALL dspevd( 'N', 'U', 2, a, x, z, 1, w, 4, iw, 0, info )
695  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
696  infot = 11
697  CALL dspevd( 'V', 'U', 2, a, x, z, 2, w, 23, iw, 11, info )
698  CALL chkxer( 'DSPEVD', infot, nout, lerr, ok )
699  nt = nt + 10
700 *
701 * DSPEV
702 *
703  srnamt = 'DSPEV '
704  infot = 1
705  CALL dspev( '/', 'U', 0, a, w, z, 1, x, info )
706  CALL chkxer( 'DSPEV ', infot, nout, lerr, ok )
707  infot = 2
708  CALL dspev( 'N', '/', 0, a, w, z, 1, x, info )
709  CALL chkxer( 'DSPEV ', infot, nout, lerr, ok )
710  infot = 3
711  CALL dspev( 'N', 'U', -1, a, w, z, 1, x, info )
712  CALL chkxer( 'DSPEV ', infot, nout, lerr, ok )
713  infot = 7
714  CALL dspev( 'V', 'U', 2, a, w, z, 1, x, info )
715  CALL chkxer( 'DSPEV ', infot, nout, lerr, ok )
716  nt = nt + 4
717 *
718 * DSPEVX
719 *
720  srnamt = 'DSPEVX'
721  infot = 1
722  CALL dspevx( '/', 'A', 'U', 0, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
723  $ x, z, 1, w, iw, i3, info )
724  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
725  infot = 2
726  CALL dspevx( 'N', '/', 'U', 0, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
727  $ x, z, 1, w, iw, i3, info )
728  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
729  infot = 3
730  CALL dspevx( 'N', 'A', '/', 0, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
731  $ x, z, 1, w, iw, i3, info )
732  infot = 4
733  CALL dspevx( 'N', 'A', 'U', -1, a, 0.0d0, 0.0d0, 0, 0, 0.0d0,
734  $ m, x, z, 1, w, iw, i3, info )
735  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
736  infot = 7
737  CALL dspevx( 'N', 'V', 'U', 1, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
738  $ x, z, 1, w, iw, i3, info )
739  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
740  infot = 8
741  CALL dspevx( 'N', 'I', 'U', 1, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
742  $ x, z, 1, w, iw, i3, info )
743  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
744  infot = 8
745  CALL dspevx( 'N', 'I', 'U', 1, a, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
746  $ x, z, 1, w, iw, i3, info )
747  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
748  infot = 9
749  CALL dspevx( 'N', 'I', 'U', 2, a, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
750  $ x, z, 1, w, iw, i3, info )
751  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
752  infot = 9
753  CALL dspevx( 'N', 'I', 'U', 1, a, 0.0d0, 0.0d0, 1, 2, 0.0d0, m,
754  $ x, z, 1, w, iw, i3, info )
755  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
756  infot = 14
757  CALL dspevx( 'V', 'A', 'U', 2, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
758  $ x, z, 1, w, iw, i3, info )
759  CALL chkxer( 'DSPEVX', infot, nout, lerr, ok )
760  nt = nt + 10
761 *
762 * Test error exits for the SB path.
763 *
764  ELSE IF( lsamen( 2, c2, 'SB' ) ) THEN
765 *
766 * DSBTRD
767 *
768  srnamt = 'DSBTRD'
769  infot = 1
770  CALL dsbtrd( '/', 'U', 0, 0, a, 1, d, e, z, 1, w, info )
771  CALL chkxer( 'DSBTRD', infot, nout, lerr, ok )
772  infot = 2
773  CALL dsbtrd( 'N', '/', 0, 0, a, 1, d, e, z, 1, w, info )
774  CALL chkxer( 'DSBTRD', infot, nout, lerr, ok )
775  infot = 3
776  CALL dsbtrd( 'N', 'U', -1, 0, a, 1, d, e, z, 1, w, info )
777  CALL chkxer( 'DSBTRD', infot, nout, lerr, ok )
778  infot = 4
779  CALL dsbtrd( 'N', 'U', 0, -1, a, 1, d, e, z, 1, w, info )
780  CALL chkxer( 'DSBTRD', infot, nout, lerr, ok )
781  infot = 6
782  CALL dsbtrd( 'N', 'U', 1, 1, a, 1, d, e, z, 1, w, info )
783  CALL chkxer( 'DSBTRD', infot, nout, lerr, ok )
784  infot = 10
785  CALL dsbtrd( 'V', 'U', 2, 0, a, 1, d, e, z, 1, w, info )
786  CALL chkxer( 'DSBTRD', infot, nout, lerr, ok )
787  nt = nt + 6
788 *
789 * DSBEVD
790 *
791  srnamt = 'DSBEVD'
792  infot = 1
793  CALL dsbevd( '/', 'U', 0, 0, a, 1, x, z, 1, w, 1, iw, 1, info )
794  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
795  infot = 2
796  CALL dsbevd( 'N', '/', 0, 0, a, 1, x, z, 1, w, 1, iw, 1, info )
797  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
798  infot = 3
799  CALL dsbevd( 'N', 'U', -1, 0, a, 1, x, z, 1, w, 1, iw, 1,
800  $ info )
801  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
802  infot = 4
803  CALL dsbevd( 'N', 'U', 0, -1, a, 1, x, z, 1, w, 1, iw, 1,
804  $ info )
805  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
806  infot = 6
807  CALL dsbevd( 'N', 'U', 2, 1, a, 1, x, z, 1, w, 4, iw, 1, info )
808  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
809  infot = 9
810  CALL dsbevd( 'V', 'U', 2, 1, a, 2, x, z, 1, w, 25, iw, 12,
811  $ info )
812  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
813  infot = 11
814  CALL dsbevd( 'N', 'U', 1, 0, a, 1, x, z, 1, w, 0, iw, 1, info )
815  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
816  infot = 11
817  CALL dsbevd( 'N', 'U', 2, 0, a, 1, x, z, 1, w, 3, iw, 1, info )
818  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
819  infot = 11
820  CALL dsbevd( 'V', 'U', 2, 0, a, 1, x, z, 2, w, 18, iw, 12,
821  $ info )
822  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
823  infot = 13
824  CALL dsbevd( 'N', 'U', 1, 0, a, 1, x, z, 1, w, 1, iw, 0, info )
825  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
826  infot = 13
827  CALL dsbevd( 'V', 'U', 2, 0, a, 1, x, z, 2, w, 25, iw, 11,
828  $ info )
829  CALL chkxer( 'DSBEVD', infot, nout, lerr, ok )
830  nt = nt + 11
831 *
832 * DSBEV
833 *
834  srnamt = 'DSBEV '
835  infot = 1
836  CALL dsbev( '/', 'U', 0, 0, a, 1, x, z, 1, w, info )
837  CALL chkxer( 'DSBEV ', infot, nout, lerr, ok )
838  infot = 2
839  CALL dsbev( 'N', '/', 0, 0, a, 1, x, z, 1, w, info )
840  CALL chkxer( 'DSBEV ', infot, nout, lerr, ok )
841  infot = 3
842  CALL dsbev( 'N', 'U', -1, 0, a, 1, x, z, 1, w, info )
843  CALL chkxer( 'DSBEV ', infot, nout, lerr, ok )
844  infot = 4
845  CALL dsbev( 'N', 'U', 0, -1, a, 1, x, z, 1, w, info )
846  CALL chkxer( 'DSBEV ', infot, nout, lerr, ok )
847  infot = 6
848  CALL dsbev( 'N', 'U', 2, 1, a, 1, x, z, 1, w, info )
849  CALL chkxer( 'DSBEV ', infot, nout, lerr, ok )
850  infot = 9
851  CALL dsbev( 'V', 'U', 2, 0, a, 1, x, z, 1, w, info )
852  CALL chkxer( 'DSBEV ', infot, nout, lerr, ok )
853  nt = nt + 6
854 *
855 * DSBEVX
856 *
857  srnamt = 'DSBEVX'
858  infot = 1
859  CALL dsbevx( '/', 'A', 'U', 0, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
860  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
861  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
862  infot = 2
863  CALL dsbevx( 'N', '/', 'U', 0, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
864  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
865  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
866  infot = 3
867  CALL dsbevx( 'N', 'A', '/', 0, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
868  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
869  infot = 4
870  CALL dsbevx( 'N', 'A', 'U', -1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
871  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
872  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
873  infot = 5
874  CALL dsbevx( 'N', 'A', 'U', 0, -1, a, 1, q, 1, 0.0d0, 0.0d0, 0,
875  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
876  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
877  infot = 7
878  CALL dsbevx( 'N', 'A', 'U', 2, 1, a, 1, q, 1, 0.0d0, 0.0d0, 0,
879  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
880  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
881  infot = 9
882  CALL dsbevx( 'V', 'A', 'U', 2, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
883  $ 0, 0.0d0, m, x, z, 2, w, iw, i3, info )
884  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
885  infot = 11
886  CALL dsbevx( 'N', 'V', 'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
887  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
888  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
889  infot = 12
890  CALL dsbevx( 'N', 'I', 'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
891  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
892  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
893  infot = 12
894  CALL dsbevx( 'N', 'I', 'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 2,
895  $ 1, 0.0d0, m, x, z, 1, w, iw, i3, info )
896  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
897  infot = 13
898  CALL dsbevx( 'N', 'I', 'U', 2, 0, a, 1, q, 1, 0.0d0, 0.0d0, 2,
899  $ 1, 0.0d0, m, x, z, 1, w, iw, i3, info )
900  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
901  infot = 13
902  CALL dsbevx( 'N', 'I', 'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 1,
903  $ 2, 0.0d0, m, x, z, 1, w, iw, i3, info )
904  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
905  infot = 18
906  CALL dsbevx( 'V', 'A', 'U', 2, 0, a, 1, q, 2, 0.0d0, 0.0d0, 0,
907  $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
908  CALL chkxer( 'DSBEVX', infot, nout, lerr, ok )
909  nt = nt + 13
910  END IF
911 *
912 * Print a summary line.
913 *
914  IF( ok ) THEN
915  WRITE( nout, fmt = 9999 )path, nt
916  ELSE
917  WRITE( nout, fmt = 9998 )path
918  END IF
919 *
920  9999 FORMAT( 1x, a3, ' routines passed the tests of the error exits',
921  $ ' (', i3, ' tests done)' )
922  9998 FORMAT( ' *** ', a3, ' routines failed the tests of the error ',
923  $ 'exits ***' )
924 *
925  RETURN
926 *
927 * End of DERRST
928 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine dstevr(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dstevr.f:306
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine dsbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
DSBTRD
Definition: dsbtrd.f:165
subroutine dstev(JOBZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: dstev.f:118
subroutine dorgtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
DORGTR
Definition: dorgtr.f:125
subroutine dsyev(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO)
DSYEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: dsyev.f:134
subroutine dsytrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
DSYTRD
Definition: dsytrd.f:194
subroutine dopmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
DOPMTR
Definition: dopmtr.f:152
subroutine dsbevx(JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSBEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dsbevx.f:267
logical function lsamen(N, CA, CB)
LSAMEN
Definition: lsamen.f:76
subroutine dspevd(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dspevd.f:181
subroutine dsbev(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, INFO)
DSBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: dsbev.f:148
subroutine dpteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DPTEQR
Definition: dpteqr.f:147
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:133
subroutine dspev(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, INFO)
DSPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: dspev.f:132
subroutine dsbevd(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSBEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dsbevd.f:195
subroutine dsyevr(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSYEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: dsyevr.f:334
subroutine dsyevx(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO)
DSYEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: dsyevx.f:255
subroutine dsyevd(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO)
DSYEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: dsyevd.f:187
subroutine chkxer(SRNAMT, INFOT, NOUT, LERR, OK)
Definition: cblat2.f:3199
subroutine dspevx(JOBZ, RANGE, UPLO, N, AP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSPEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dspevx.f:236
subroutine dormtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMTR
Definition: dormtr.f:173
subroutine dsptrd(UPLO, N, AP, D, E, TAU, INFO)
DSPTRD
Definition: dsptrd.f:152
subroutine dstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEDC
Definition: dstedc.f:191
subroutine dstevd(JOBZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dstevd.f:165
subroutine dopgtr(UPLO, N, AP, TAU, Q, LDQ, WORK, INFO)
DOPGTR
Definition: dopgtr.f:116
subroutine dstevx(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: dstevx.f:229
subroutine dstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEIN
Definition: dstein.f:176

Here is the call graph for this function:

Here is the caller graph for this function: