76 parameter( nmax = 3, liw = 12*nmax, lw = 20*nmax )
80 INTEGER I, INFO, J, M, N, NSPLIT, NT
83 INTEGER I1( NMAX ), I2( NMAX ), I3( NMAX ), IW( LIW )
84 DOUBLE PRECISION A( NMAX, NMAX ), C( NMAX, NMAX ), D( NMAX ),
85 $ E( NMAX ), Q( NMAX, NMAX ), R( NMAX ),
86 $ TAU( NMAX ), W( LW ), X( NMAX ),
110 COMMON / infoc / infot, nout, ok, lerr
111 COMMON / srnamc / srnamt
119 WRITE( nout, fmt = * )
126 a( i, j ) = 1.d0 / dble( i+j )
141 IF(
lsamen( 2, c2,
'ST' ) )
THEN
147 CALL dsytrd(
'/', 0, a, 1, d, e, tau, w, 1, info )
148 CALL chkxer(
'DSYTRD', infot, nout, lerr, ok )
150 CALL dsytrd(
'U', -1, a, 1, d, e, tau, w, 1, info )
151 CALL chkxer(
'DSYTRD', infot, nout, lerr, ok )
153 CALL dsytrd(
'U', 2, a, 1, d, e, tau, w, 1, info )
154 CALL chkxer(
'DSYTRD', infot, nout, lerr, ok )
156 CALL dsytrd(
'U', 0, a, 1, d, e, tau, w, 0, info )
157 CALL chkxer(
'DSYTRD', infot, nout, lerr, ok )
162 srnamt =
'DSYTRD_2STAGE'
166 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
170 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
174 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
178 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
182 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
186 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
190 CALL chkxer(
'DSYTRD_2STAGE', infot, nout, lerr, ok )
195 srnamt =
'DSYTRD_SY2SB'
197 CALL dsytrd_sy2sb(
'/', 0, 0, a, 1, c, 1, tau, w, 1, info )
198 CALL chkxer(
'DSYTRD_SY2SB', infot, nout, lerr, ok )
200 CALL dsytrd_sy2sb(
'U', -1, 0, a, 1, c, 1, tau, w, 1, info )
201 CALL chkxer(
'DSYTRD_SY2SB', infot, nout, lerr, ok )
203 CALL dsytrd_sy2sb(
'U', 0, -1, a, 1, c, 1, tau, w, 1, info )
204 CALL chkxer(
'DSYTRD_SY2SB', infot, nout, lerr, ok )
206 CALL dsytrd_sy2sb(
'U', 2, 0, a, 1, c, 1, tau, w, 1, info )
207 CALL chkxer(
'DSYTRD_SY2SB', infot, nout, lerr, ok )
209 CALL dsytrd_sy2sb(
'U', 0, 2, a, 1, c, 1, tau, w, 1, info )
210 CALL chkxer(
'DSYTRD_SY2SB', infot, nout, lerr, ok )
212 CALL dsytrd_sy2sb(
'U', 0, 0, a, 1, c, 1, tau, w, 0, info )
213 CALL chkxer(
'DSYTRD_SY2SB', infot, nout, lerr, ok )
218 srnamt =
'DSYTRD_SB2ST'
222 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
226 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
230 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
234 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
238 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
242 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
246 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
250 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
254 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
261 CALL dorgtr(
'/', 0, a, 1, tau, w, 1, info )
262 CALL chkxer(
'DORGTR', infot, nout, lerr, ok )
264 CALL dorgtr(
'U', -1, a, 1, tau, w, 1, info )
265 CALL chkxer(
'DORGTR', infot, nout, lerr, ok )
267 CALL dorgtr(
'U', 2, a, 1, tau, w, 1, info )
268 CALL chkxer(
'DORGTR', infot, nout, lerr, ok )
270 CALL dorgtr(
'U', 3, a, 3, tau, w, 1, info )
271 CALL chkxer(
'DORGTR', infot, nout, lerr, ok )
278 CALL dormtr(
'/',
'U',
'N', 0, 0, a, 1, tau, c, 1, w, 1, info )
279 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
281 CALL dormtr(
'L',
'/',
'N', 0, 0, a, 1, tau, c, 1, w, 1, info )
282 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
284 CALL dormtr(
'L',
'U',
'/', 0, 0, a, 1, tau, c, 1, w, 1, info )
285 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
287 CALL dormtr(
'L',
'U',
'N', -1, 0, a, 1, tau, c, 1, w, 1,
289 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
291 CALL dormtr(
'L',
'U',
'N', 0, -1, a, 1, tau, c, 1, w, 1,
293 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
295 CALL dormtr(
'L',
'U',
'N', 2, 0, a, 1, tau, c, 2, w, 1, info )
296 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
298 CALL dormtr(
'R',
'U',
'N', 0, 2, a, 1, tau, c, 1, w, 1, info )
299 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
301 CALL dormtr(
'L',
'U',
'N', 2, 0, a, 2, tau, c, 1, w, 1, info )
302 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
304 CALL dormtr(
'L',
'U',
'N', 0, 2, a, 1, tau, c, 1, w, 1, info )
305 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
307 CALL dormtr(
'R',
'U',
'N', 2, 0, a, 1, tau, c, 2, w, 1, info )
308 CALL chkxer(
'DORMTR', infot, nout, lerr, ok )
315 CALL dsptrd(
'/', 0, a, d, e, tau, info )
316 CALL chkxer(
'DSPTRD', infot, nout, lerr, ok )
318 CALL dsptrd(
'U', -1, a, d, e, tau, info )
319 CALL chkxer(
'DSPTRD', infot, nout, lerr, ok )
326 CALL dopgtr(
'/', 0, a, tau, z, 1, w, info )
327 CALL chkxer(
'DOPGTR', infot, nout, lerr, ok )
329 CALL dopgtr(
'U', -1, a, tau, z, 1, w, info )
330 CALL chkxer(
'DOPGTR', infot, nout, lerr, ok )
332 CALL dopgtr(
'U', 2, a, tau, z, 1, w, info )
333 CALL chkxer(
'DOPGTR', infot, nout, lerr, ok )
340 CALL dopmtr(
'/',
'U',
'N', 0, 0, a, tau, c, 1, w, info )
341 CALL chkxer(
'DOPMTR', infot, nout, lerr, ok )
343 CALL dopmtr(
'L',
'/',
'N', 0, 0, a, tau, c, 1, w, info )
344 CALL chkxer(
'DOPMTR', infot, nout, lerr, ok )
346 CALL dopmtr(
'L',
'U',
'/', 0, 0, a, tau, c, 1, w, info )
347 CALL chkxer(
'DOPMTR', infot, nout, lerr, ok )
349 CALL dopmtr(
'L',
'U',
'N', -1, 0, a, tau, c, 1, w, info )
350 CALL chkxer(
'DOPMTR', infot, nout, lerr, ok )
352 CALL dopmtr(
'L',
'U',
'N', 0, -1, a, tau, c, 1, w, info )
353 CALL chkxer(
'DOPMTR', infot, nout, lerr, ok )
355 CALL dopmtr(
'L',
'U',
'N', 2, 0, a, tau, c, 1, w, info )
356 CALL chkxer(
'DOPMTR', infot, nout, lerr, ok )
363 CALL dpteqr(
'/', 0, d, e, z, 1, w, info )
364 CALL chkxer(
'DPTEQR', infot, nout, lerr, ok )
366 CALL dpteqr(
'N', -1, d, e, z, 1, w, info )
367 CALL chkxer(
'DPTEQR', infot, nout, lerr, ok )
369 CALL dpteqr(
'V', 2, d, e, z, 1, w, info )
370 CALL chkxer(
'DPTEQR', infot, nout, lerr, ok )
377 CALL dstebz(
'/',
'E', 0, 0.0d0, 1.0d0, 1, 0, 0.0d0, d, e, m,
378 $ nsplit, x, i1, i2, w, iw, info )
379 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
381 CALL dstebz(
'A',
'/', 0, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
382 $ nsplit, x, i1, i2, w, iw, info )
383 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
385 CALL dstebz(
'A',
'E', -1, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
386 $ nsplit, x, i1, i2, w, iw, info )
387 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
389 CALL dstebz(
'V',
'E', 0, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
390 $ nsplit, x, i1, i2, w, iw, info )
391 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
393 CALL dstebz(
'I',
'E', 0, 0.0d0, 0.0d0, 0, 0, 0.0d0, d, e, m,
394 $ nsplit, x, i1, i2, w, iw, info )
395 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
397 CALL dstebz(
'I',
'E', 1, 0.0d0, 0.0d0, 2, 1, 0.0d0, d, e, m,
398 $ nsplit, x, i1, i2, w, iw, info )
399 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
401 CALL dstebz(
'I',
'E', 1, 0.0d0, 0.0d0, 1, 0, 0.0d0, d, e, m,
402 $ nsplit, x, i1, i2, w, iw, info )
403 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
405 CALL dstebz(
'I',
'E', 1, 0.0d0, 0.0d0, 1, 2, 0.0d0, d, e, m,
406 $ nsplit, x, i1, i2, w, iw, info )
407 CALL chkxer(
'DSTEBZ', infot, nout, lerr, ok )
414 CALL dstein( -1, d, e, 0, x, i1, i2, z, 1, w, iw, i3, info )
415 CALL chkxer(
'DSTEIN', infot, nout, lerr, ok )
417 CALL dstein( 0, d, e, -1, x, i1, i2, z, 1, w, iw, i3, info )
418 CALL chkxer(
'DSTEIN', infot, nout, lerr, ok )
420 CALL dstein( 0, d, e, 1, x, i1, i2, z, 1, w, iw, i3, info )
421 CALL chkxer(
'DSTEIN', infot, nout, lerr, ok )
423 CALL dstein( 2, d, e, 0, x, i1, i2, z, 1, w, iw, i3, info )
424 CALL chkxer(
'DSTEIN', infot, nout, lerr, ok )
431 CALL dsteqr(
'/', 0, d, e, z, 1, w, info )
432 CALL chkxer(
'DSTEQR', infot, nout, lerr, ok )
434 CALL dsteqr(
'N', -1, d, e, z, 1, w, info )
435 CALL chkxer(
'DSTEQR', infot, nout, lerr, ok )
437 CALL dsteqr(
'V', 2, d, e, z, 1, w, info )
438 CALL chkxer(
'DSTEQR', infot, nout, lerr, ok )
445 CALL dsterf( -1, d, e, info )
446 CALL chkxer(
'DSTERF', infot, nout, lerr, ok )
453 CALL dstedc(
'/', 0, d, e, z, 1, w, 1, iw, 1, info )
454 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
456 CALL dstedc(
'N', -1, d, e, z, 1, w, 1, iw, 1, info )
457 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
459 CALL dstedc(
'V', 2, d, e, z, 1, w, 23, iw, 28, info )
460 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
462 CALL dstedc(
'N', 1, d, e, z, 1, w, 0, iw, 1, info )
463 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
465 CALL dstedc(
'I', 2, d, e, z, 2, w, 0, iw, 12, info )
466 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
468 CALL dstedc(
'V', 2, d, e, z, 2, w, 0, iw, 28, info )
469 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
471 CALL dstedc(
'N', 1, d, e, z, 1, w, 1, iw, 0, info )
472 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
474 CALL dstedc(
'I', 2, d, e, z, 2, w, 19, iw, 0, info )
475 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
477 CALL dstedc(
'V', 2, d, e, z, 2, w, 23, iw, 0, info )
478 CALL chkxer(
'DSTEDC', infot, nout, lerr, ok )
485 CALL dstevd(
'/', 0, d, e, z, 1, w, 1, iw, 1, info )
486 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
488 CALL dstevd(
'N', -1, d, e, z, 1, w, 1, iw, 1, info )
489 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
491 CALL dstevd(
'V', 2, d, e, z, 1, w, 19, iw, 12, info )
492 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
494 CALL dstevd(
'N', 1, d, e, z, 1, w, 0, iw, 1, info )
495 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
497 CALL dstevd(
'V', 2, d, e, z, 2, w, 12, iw, 12, info )
498 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
500 CALL dstevd(
'N', 0, d, e, z, 1, w, 1, iw, 0, info )
501 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
503 CALL dstevd(
'V', 2, d, e, z, 2, w, 19, iw, 11, info )
504 CALL chkxer(
'DSTEVD', infot, nout, lerr, ok )
511 CALL dstev(
'/', 0, d, e, z, 1, w, info )
512 CALL chkxer(
'DSTEV ', infot, nout, lerr, ok )
514 CALL dstev(
'N', -1, d, e, z, 1, w, info )
515 CALL chkxer(
'DSTEV ', infot, nout, lerr, ok )
517 CALL dstev(
'V', 2, d, e, z, 1, w, info )
518 CALL chkxer(
'DSTEV ', infot, nout, lerr, ok )
525 CALL dstevx(
'/',
'A', 0, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
526 $ x, z, 1, w, iw, i3, info )
527 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
529 CALL dstevx(
'N',
'/', 0, d, e, 0.0d0, 1.0d0, 1, 0, 0.0d0, m,
530 $ x, z, 1, w, iw, i3, info )
531 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
533 CALL dstevx(
'N',
'A', -1, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
534 $ x, z, 1, w, iw, i3, info )
535 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
537 CALL dstevx(
'N',
'V', 1, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
538 $ x, z, 1, w, iw, i3, info )
539 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
541 CALL dstevx(
'N',
'I', 1, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
542 $ x, z, 1, w, iw, i3, info )
543 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
545 CALL dstevx(
'N',
'I', 1, d, e, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
546 $ x, z, 1, w, iw, i3, info )
547 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
549 CALL dstevx(
'N',
'I', 2, d, e, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
550 $ x, z, 1, w, iw, i3, info )
551 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
553 CALL dstevx(
'N',
'I', 1, d, e, 0.0d0, 0.0d0, 1, 2, 0.0d0, m,
554 $ x, z, 1, w, iw, i3, info )
555 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
557 CALL dstevx(
'V',
'A', 2, d, e, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
558 $ x, z, 1, w, iw, i3, info )
559 CALL chkxer(
'DSTEVX', infot, nout, lerr, ok )
567 CALL dstevr(
'/',
'A', 0, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
568 $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
569 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
571 CALL dstevr(
'V',
'/', 0, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
572 $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
573 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
575 CALL dstevr(
'V',
'A', -1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
576 $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
577 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
579 CALL dstevr(
'V',
'V', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
580 $ r, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
581 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
583 CALL dstevr(
'V',
'I', 1, d, e, 0.0d0, 0.0d0, 0, 1, 0.0d0, m,
584 $ w, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
585 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
588 CALL dstevr(
'V',
'I', 2, d, e, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
589 $ w, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
590 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
593 CALL dstevr(
'V',
'I', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
594 $ w, z, 0, iw, x, 20*n, iw( 2*n+1 ), 10*n, info )
595 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
597 CALL dstevr(
'V',
'I', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
598 $ w, z, 1, iw, x, 20*n-1, iw( 2*n+1 ), 10*n, info )
599 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
601 CALL dstevr(
'V',
'I', 1, d, e, 0.0d0, 0.0d0, 1, 1, 0.0d0, m,
602 $ w, z, 1, iw, x, 20*n, iw( 2*n+1 ), 10*n-1, info )
603 CALL chkxer(
'DSTEVR', infot, nout, lerr, ok )
610 CALL dsyevd(
'/',
'U', 0, a, 1, x, w, 1, iw, 1, info )
611 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
613 CALL dsyevd(
'N',
'/', 0, a, 1, x, w, 1, iw, 1, info )
614 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
616 CALL dsyevd(
'N',
'U', -1, a, 1, x, w, 1, iw, 1, info )
617 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
619 CALL dsyevd(
'N',
'U', 2, a, 1, x, w, 3, iw, 1, info )
620 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
622 CALL dsyevd(
'N',
'U', 1, a, 1, x, w, 0, iw, 1, info )
623 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
625 CALL dsyevd(
'N',
'U', 2, a, 2, x, w, 4, iw, 1, info )
626 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
628 CALL dsyevd(
'V',
'U', 2, a, 2, x, w, 20, iw, 12, info )
629 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
631 CALL dsyevd(
'N',
'U', 1, a, 1, x, w, 1, iw, 0, info )
632 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
634 CALL dsyevd(
'N',
'U', 2, a, 2, x, w, 5, iw, 0, info )
635 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
637 CALL dsyevd(
'V',
'U', 2, a, 2, x, w, 27, iw, 11, info )
638 CALL chkxer(
'DSYEVD', infot, nout, lerr, ok )
643 srnamt =
'DSYEVD_2STAGE'
645 CALL dsyevd_2stage(
'/',
'U', 0, a, 1, x, w, 1, iw, 1, info )
646 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
648 CALL dsyevd_2stage(
'V',
'U', 0, a, 1, x, w, 1, iw, 1, info )
649 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
651 CALL dsyevd_2stage(
'N',
'/', 0, a, 1, x, w, 1, iw, 1, info )
652 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
654 CALL dsyevd_2stage(
'N',
'U', -1, a, 1, x, w, 1, iw, 1, info )
655 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
657 CALL dsyevd_2stage(
'N',
'U', 2, a, 1, x, w, 3, iw, 1, info )
658 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
660 CALL dsyevd_2stage(
'N',
'U', 1, a, 1, x, w, 0, iw, 1, info )
661 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
663 CALL dsyevd_2stage(
'N',
'U', 2, a, 2, x, w, 4, iw, 1, info )
664 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
669 CALL dsyevd_2stage(
'N',
'U', 1, a, 1, x, w, 1, iw, 0, info )
670 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
672 CALL dsyevd_2stage(
'N',
'U', 2, a, 2, x, w, 25, iw, 0, info )
673 CALL chkxer(
'DSYEVD_2STAGE', infot, nout, lerr, ok )
684 CALL dsyevr(
'/',
'A',
'U', 0, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
685 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
686 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
688 CALL dsyevr(
'V',
'/',
'U', 0, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
689 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
690 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
692 CALL dsyevr(
'V',
'A',
'/', -1, a, 1, 0.0d0, 0.0d0, 1, 1,
693 $ 0.0d0, m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n,
695 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
697 CALL dsyevr(
'V',
'A',
'U', -1, a, 1, 0.0d0, 0.0d0, 1, 1,
698 $ 0.0d0, m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n,
700 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
702 CALL dsyevr(
'V',
'A',
'U', 2, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
703 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
704 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
706 CALL dsyevr(
'V',
'V',
'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
707 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
708 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
710 CALL dsyevr(
'V',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 0, 1, 0.0d0,
711 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
712 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
715 CALL dsyevr(
'V',
'I',
'U', 2, a, 2, 0.0d0, 0.0d0, 2, 1, 0.0d0,
716 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
717 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
719 CALL dsyevr(
'V',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
720 $ m, r, z, 0, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
721 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
723 CALL dsyevr(
'V',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
724 $ m, r, z, 1, iw, q, 26*n-1, iw( 2*n+1 ), 10*n,
726 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
728 CALL dsyevr(
'V',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 1, 1, 0.0d0,
729 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n-1,
731 CALL chkxer(
'DSYEVR', infot, nout, lerr, ok )
736 srnamt =
'DSYEVR_2STAGE'
740 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
741 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
742 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
745 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
746 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
747 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
750 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
751 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
752 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
755 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
756 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
757 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
760 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
761 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
762 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
765 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
766 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
767 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
770 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
771 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
772 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
775 $ 0.0d0, 0.0d0, 0, 1, 0.0d0,
776 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
777 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
780 $ 0.0d0, 0.0d0, 2, 1, 0.0d0,
781 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
782 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
785 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
786 $ m, r, z, 0, iw, q, 26*n, iw( 2*n+1 ), 10*n, info )
787 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
790 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
791 $ m, r, z, 1, iw, q, 0, iw( 2*n+1 ), 10*n,
793 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
796 $ 0.0d0, 0.0d0, 1, 1, 0.0d0,
797 $ m, r, z, 1, iw, q, 26*n, iw( 2*n+1 ), 0,
799 CALL chkxer(
'DSYEVR_2STAGE', infot, nout, lerr, ok )
806 CALL dsyev(
'/',
'U', 0, a, 1, x, w, 1, info )
807 CALL chkxer(
'DSYEV ', infot, nout, lerr, ok )
809 CALL dsyev(
'N',
'/', 0, a, 1, x, w, 1, info )
810 CALL chkxer(
'DSYEV ', infot, nout, lerr, ok )
812 CALL dsyev(
'N',
'U', -1, a, 1, x, w, 1, info )
813 CALL chkxer(
'DSYEV ', infot, nout, lerr, ok )
815 CALL dsyev(
'N',
'U', 2, a, 1, x, w, 3, info )
816 CALL chkxer(
'DSYEV ', infot, nout, lerr, ok )
818 CALL dsyev(
'N',
'U', 1, a, 1, x, w, 1, info )
819 CALL chkxer(
'DSYEV ', infot, nout, lerr, ok )
824 srnamt =
'DSYEV_2STAGE '
827 CALL chkxer(
'DSYEV_2STAGE ', infot, nout, lerr, ok )
830 CALL chkxer(
'DSYEV_2STAGE ', infot, nout, lerr, ok )
833 CALL chkxer(
'DSYEV_2STAGE ', infot, nout, lerr, ok )
836 CALL chkxer(
'DSYEV_2STAGE ', infot, nout, lerr, ok )
839 CALL chkxer(
'DSYEV_2STAGE ', infot, nout, lerr, ok )
842 CALL chkxer(
'DSYEV_2STAGE ', infot, nout, lerr, ok )
849 CALL dsyevx(
'/',
'A',
'U', 0, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
850 $ m, x, z, 1, w, 1, iw, i3, info )
851 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
853 CALL dsyevx(
'N',
'/',
'U', 0, a, 1, 0.0d0, 1.0d0, 1, 0, 0.0d0,
854 $ m, x, z, 1, w, 1, iw, i3, info )
855 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
857 CALL dsyevx(
'N',
'A',
'/', 0, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
858 $ m, x, z, 1, w, 1, iw, i3, info )
860 CALL dsyevx(
'N',
'A',
'U', -1, a, 1, 0.0d0, 0.0d0, 0, 0,
861 $ 0.0d0, m, x, z, 1, w, 1, iw, i3, info )
862 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
864 CALL dsyevx(
'N',
'A',
'U', 2, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
865 $ m, x, z, 1, w, 16, iw, i3, info )
866 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
868 CALL dsyevx(
'N',
'V',
'U', 1, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
869 $ m, x, z, 1, w, 8, iw, i3, info )
870 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
872 CALL dsyevx(
'N',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
873 $ m, x, z, 1, w, 8, iw, i3, info )
874 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
876 CALL dsyevx(
'N',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 2, 1, 0.0d0,
877 $ m, x, z, 1, w, 8, iw, i3, info )
878 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
880 CALL dsyevx(
'N',
'I',
'U', 2, a, 2, 0.0d0, 0.0d0, 2, 1, 0.0d0,
881 $ m, x, z, 1, w, 16, iw, i3, info )
882 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
884 CALL dsyevx(
'N',
'I',
'U', 1, a, 1, 0.0d0, 0.0d0, 1, 2, 0.0d0,
885 $ m, x, z, 1, w, 8, iw, i3, info )
886 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
888 CALL dsyevx(
'V',
'A',
'U', 2, a, 2, 0.0d0, 0.0d0, 0, 0, 0.0d0,
889 $ m, x, z, 1, w, 16, iw, i3, info )
890 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
892 CALL dsyevx(
'V',
'A',
'U', 1, a, 1, 0.0d0, 0.0d0, 0, 0, 0.0d0,
893 $ m, x, z, 1, w, 0, iw, i3, info )
894 CALL chkxer(
'DSYEVX', infot, nout, lerr, ok )
899 srnamt =
'DSYEVX_2STAGE'
902 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
903 $ m, x, z, 1, w, 1, iw, i3, info )
904 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
907 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
908 $ m, x, z, 1, w, 1, iw, i3, info )
909 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
912 $ 0.0d0, 1.0d0, 1, 0, 0.0d0,
913 $ m, x, z, 1, w, 1, iw, i3, info )
914 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
917 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
918 $ m, x, z, 1, w, 1, iw, i3, info )
921 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
922 $ m, x, z, 1, w, 1, iw, i3, info )
923 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
926 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
927 $ m, x, z, 1, w, 16, iw, i3, info )
928 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
931 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
932 $ m, x, z, 1, w, 8, iw, i3, info )
933 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
936 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
937 $ m, x, z, 1, w, 8, iw, i3, info )
938 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
941 $ 0.0d0, 0.0d0, 2, 1, 0.0d0,
942 $ m, x, z, 1, w, 8, iw, i3, info )
943 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
946 $ 0.0d0, 0.0d0, 2, 1, 0.0d0,
947 $ m, x, z, 1, w, 16, iw, i3, info )
948 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
951 $ 0.0d0, 0.0d0, 1, 2, 0.0d0,
952 $ m, x, z, 1, w, 8, iw, i3, info )
953 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
956 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
957 $ m, x, z, 0, w, 16, iw, i3, info )
958 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
961 $ 0.0d0, 0.0d0, 0, 0, 0.0d0,
962 $ m, x, z, 1, w, 0, iw, i3, info )
963 CALL chkxer(
'DSYEVX_2STAGE', infot, nout, lerr, ok )
970 CALL dspevd(
'/',
'U', 0, a, x, z, 1, w, 1, iw, 1, info )
971 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
973 CALL dspevd(
'N',
'/', 0, a, x, z, 1, w, 1, iw, 1, info )
974 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
976 CALL dspevd(
'N',
'U', -1, a, x, z, 1, w, 1, iw, 1, info )
977 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
979 CALL dspevd(
'V',
'U', 2, a, x, z, 1, w, 23, iw, 12, info )
980 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
982 CALL dspevd(
'N',
'U', 1, a, x, z, 1, w, 0, iw, 1, info )
983 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
985 CALL dspevd(
'N',
'U', 2, a, x, z, 1, w, 3, iw, 1, info )
986 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
988 CALL dspevd(
'V',
'U', 2, a, x, z, 2, w, 16, iw, 12, info )
989 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
991 CALL dspevd(
'N',
'U', 1, a, x, z, 1, w, 1, iw, 0, info )
992 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
994 CALL dspevd(
'N',
'U', 2, a, x, z, 1, w, 4, iw, 0, info )
995 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
997 CALL dspevd(
'V',
'U', 2, a, x, z, 2, w, 23, iw, 11, info )
998 CALL chkxer(
'DSPEVD', infot, nout, lerr, ok )
1005 CALL dspev(
'/',
'U', 0, a, w, z, 1, x, info )
1006 CALL chkxer(
'DSPEV ', infot, nout, lerr, ok )
1008 CALL dspev(
'N',
'/', 0, a, w, z, 1, x, info )
1009 CALL chkxer(
'DSPEV ', infot, nout, lerr, ok )
1011 CALL dspev(
'N',
'U', -1, a, w, z, 1, x, info )
1012 CALL chkxer(
'DSPEV ', infot, nout, lerr, ok )
1014 CALL dspev(
'V',
'U', 2, a, w, z, 1, x, info )
1015 CALL chkxer(
'DSPEV ', infot, nout, lerr, ok )
1022 CALL dspevx(
'/',
'A',
'U', 0, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
1023 $ x, z, 1, w, iw, i3, info )
1024 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1026 CALL dspevx(
'N',
'/',
'U', 0, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
1027 $ x, z, 1, w, iw, i3, info )
1028 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1030 CALL dspevx(
'N',
'A',
'/', 0, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
1031 $ x, z, 1, w, iw, i3, info )
1033 CALL dspevx(
'N',
'A',
'U', -1, a, 0.0d0, 0.0d0, 0, 0, 0.0d0,
1034 $ m, x, z, 1, w, iw, i3, info )
1035 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1037 CALL dspevx(
'N',
'V',
'U', 1, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
1038 $ x, z, 1, w, iw, i3, info )
1039 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1041 CALL dspevx(
'N',
'I',
'U', 1, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
1042 $ x, z, 1, w, iw, i3, info )
1043 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1045 CALL dspevx(
'N',
'I',
'U', 1, a, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
1046 $ x, z, 1, w, iw, i3, info )
1047 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1049 CALL dspevx(
'N',
'I',
'U', 2, a, 0.0d0, 0.0d0, 2, 1, 0.0d0, m,
1050 $ x, z, 1, w, iw, i3, info )
1051 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1053 CALL dspevx(
'N',
'I',
'U', 1, a, 0.0d0, 0.0d0, 1, 2, 0.0d0, m,
1054 $ x, z, 1, w, iw, i3, info )
1055 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1057 CALL dspevx(
'V',
'A',
'U', 2, a, 0.0d0, 0.0d0, 0, 0, 0.0d0, m,
1058 $ x, z, 1, w, iw, i3, info )
1059 CALL chkxer(
'DSPEVX', infot, nout, lerr, ok )
1064 ELSE IF(
lsamen( 2, c2,
'SB' ) )
THEN
1070 CALL dsbtrd(
'/',
'U', 0, 0, a, 1, d, e, z, 1, w, info )
1071 CALL chkxer(
'DSBTRD', infot, nout, lerr, ok )
1073 CALL dsbtrd(
'N',
'/', 0, 0, a, 1, d, e, z, 1, w, info )
1074 CALL chkxer(
'DSBTRD', infot, nout, lerr, ok )
1076 CALL dsbtrd(
'N',
'U', -1, 0, a, 1, d, e, z, 1, w, info )
1077 CALL chkxer(
'DSBTRD', infot, nout, lerr, ok )
1079 CALL dsbtrd(
'N',
'U', 0, -1, a, 1, d, e, z, 1, w, info )
1080 CALL chkxer(
'DSBTRD', infot, nout, lerr, ok )
1082 CALL dsbtrd(
'N',
'U', 1, 1, a, 1, d, e, z, 1, w, info )
1083 CALL chkxer(
'DSBTRD', infot, nout, lerr, ok )
1085 CALL dsbtrd(
'V',
'U', 2, 0, a, 1, d, e, z, 1, w, info )
1086 CALL chkxer(
'DSBTRD', infot, nout, lerr, ok )
1091 srnamt =
'DSYTRD_SB2ST'
1094 $ c, 1, w, 1, info )
1095 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1098 $ c, 1, w, 1, info )
1099 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1102 $ c, 1, w, 1, info )
1103 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1106 $ c, 1, w, 1, info )
1107 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1110 $ c, 1, w, 1, info )
1111 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1114 $ c, 1, w, 1, info )
1115 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1118 $ c, 1, w, 1, info )
1119 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1122 $ c, 0, w, 1, info )
1123 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1126 $ c, 1, w, 0, info )
1127 CALL chkxer(
'DSYTRD_SB2ST', infot, nout, lerr, ok )
1134 CALL dsbevd(
'/',
'U', 0, 0, a, 1, x, z, 1, w, 1, iw, 1, info )
1135 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1137 CALL dsbevd(
'N',
'/', 0, 0, a, 1, x, z, 1, w, 1, iw, 1, info )
1138 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1140 CALL dsbevd(
'N',
'U', -1, 0, a, 1, x, z, 1, w, 1, iw, 1,
1142 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1144 CALL dsbevd(
'N',
'U', 0, -1, a, 1, x, z, 1, w, 1, iw, 1,
1146 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1148 CALL dsbevd(
'N',
'U', 2, 1, a, 1, x, z, 1, w, 4, iw, 1, info )
1149 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1151 CALL dsbevd(
'V',
'U', 2, 1, a, 2, x, z, 1, w, 25, iw, 12,
1153 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1155 CALL dsbevd(
'N',
'U', 1, 0, a, 1, x, z, 1, w, 0, iw, 1, info )
1156 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1158 CALL dsbevd(
'N',
'U', 2, 0, a, 1, x, z, 1, w, 3, iw, 1, info )
1159 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1161 CALL dsbevd(
'V',
'U', 2, 0, a, 1, x, z, 2, w, 18, iw, 12,
1163 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1165 CALL dsbevd(
'N',
'U', 1, 0, a, 1, x, z, 1, w, 1, iw, 0, info )
1166 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1168 CALL dsbevd(
'V',
'U', 2, 0, a, 1, x, z, 2, w, 25, iw, 11,
1170 CALL chkxer(
'DSBEVD', infot, nout, lerr, ok )
1175 srnamt =
'DSBEVD_2STAGE'
1177 CALL dsbevd_2stage(
'/',
'U', 0, 0, a, 1, x, z, 1, w,
1179 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1181 CALL dsbevd_2stage(
'V',
'U', 0, 0, a, 1, x, z, 1, w,
1183 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1185 CALL dsbevd_2stage(
'N',
'/', 0, 0, a, 1, x, z, 1, w,
1187 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1189 CALL dsbevd_2stage(
'N',
'U', -1, 0, a, 1, x, z, 1, w,
1191 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1193 CALL dsbevd_2stage(
'N',
'U', 0, -1, a, 1, x, z, 1, w,
1195 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1197 CALL dsbevd_2stage(
'N',
'U', 2, 1, a, 1, x, z, 1, w,
1199 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1205 CALL dsbevd_2stage(
'N',
'U', 1, 0, a, 1, x, z, 1, w,
1207 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1209 CALL dsbevd_2stage(
'N',
'U', 2, 0, a, 1, x, z, 1, w,
1211 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1217 CALL dsbevd_2stage(
'N',
'U', 1, 0, a, 1, x, z, 1, w,
1219 CALL chkxer(
'DSBEVD_2STAGE', infot, nout, lerr, ok )
1231 CALL dsbev(
'/',
'U', 0, 0, a, 1, x, z, 1, w, info )
1232 CALL chkxer(
'DSBEV ', infot, nout, lerr, ok )
1234 CALL dsbev(
'N',
'/', 0, 0, a, 1, x, z, 1, w, info )
1235 CALL chkxer(
'DSBEV ', infot, nout, lerr, ok )
1237 CALL dsbev(
'N',
'U', -1, 0, a, 1, x, z, 1, w, info )
1238 CALL chkxer(
'DSBEV ', infot, nout, lerr, ok )
1240 CALL dsbev(
'N',
'U', 0, -1, a, 1, x, z, 1, w, info )
1241 CALL chkxer(
'DSBEV ', infot, nout, lerr, ok )
1243 CALL dsbev(
'N',
'U', 2, 1, a, 1, x, z, 1, w, info )
1244 CALL chkxer(
'DSBEV ', infot, nout, lerr, ok )
1246 CALL dsbev(
'V',
'U', 2, 0, a, 1, x, z, 1, w, info )
1247 CALL chkxer(
'DSBEV ', infot, nout, lerr, ok )
1252 srnamt =
'DSBEV_2STAGE '
1254 CALL dsbev_2stage(
'/',
'U', 0, 0, a, 1, x, z, 1, w, 0, info )
1255 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1257 CALL dsbev_2stage(
'V',
'U', 0, 0, a, 1, x, z, 1, w, 0, info )
1258 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1260 CALL dsbev_2stage(
'N',
'/', 0, 0, a, 1, x, z, 1, w, 0, info )
1261 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1263 CALL dsbev_2stage(
'N',
'U', -1, 0, a, 1, x, z, 1, w, 0, info )
1264 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1266 CALL dsbev_2stage(
'N',
'U', 0, -1, a, 1, x, z, 1, w, 0, info )
1267 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1269 CALL dsbev_2stage(
'N',
'U', 2, 1, a, 1, x, z, 1, w, 0, info )
1270 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1272 CALL dsbev_2stage(
'N',
'U', 2, 0, a, 1, x, z, 0, w, 0, info )
1273 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1275 CALL dsbev_2stage(
'N',
'U', 0, 0, a, 1, x, z, 1, w, 0, info )
1276 CALL chkxer(
'DSBEV_2STAGE ', infot, nout, lerr, ok )
1283 CALL dsbevx(
'/',
'A',
'U', 0, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1284 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1285 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1287 CALL dsbevx(
'N',
'/',
'U', 0, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1288 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1289 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1291 CALL dsbevx(
'N',
'A',
'/', 0, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1292 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1293 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1295 CALL dsbevx(
'N',
'A',
'U', -1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1296 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1297 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1299 CALL dsbevx(
'N',
'A',
'U', 0, -1, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1300 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1301 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1303 CALL dsbevx(
'N',
'A',
'U', 2, 1, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1304 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1305 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1307 CALL dsbevx(
'V',
'A',
'U', 2, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1308 $ 0, 0.0d0, m, x, z, 2, w, iw, i3, info )
1309 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1311 CALL dsbevx(
'N',
'V',
'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1312 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1313 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1315 CALL dsbevx(
'N',
'I',
'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 0,
1316 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1317 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1319 CALL dsbevx(
'N',
'I',
'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 2,
1320 $ 1, 0.0d0, m, x, z, 1, w, iw, i3, info )
1321 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1323 CALL dsbevx(
'N',
'I',
'U', 2, 0, a, 1, q, 1, 0.0d0, 0.0d0, 2,
1324 $ 1, 0.0d0, m, x, z, 1, w, iw, i3, info )
1325 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1327 CALL dsbevx(
'N',
'I',
'U', 1, 0, a, 1, q, 1, 0.0d0, 0.0d0, 1,
1328 $ 2, 0.0d0, m, x, z, 1, w, iw, i3, info )
1329 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1331 CALL dsbevx(
'V',
'A',
'U', 2, 0, a, 1, q, 2, 0.0d0, 0.0d0, 0,
1332 $ 0, 0.0d0, m, x, z, 1, w, iw, i3, info )
1333 CALL chkxer(
'DSBEVX', infot, nout, lerr, ok )
1338 srnamt =
'DSBEVX_2STAGE'
1340 CALL dsbevx_2stage(
'/',
'A',
'U', 0, 0, a, 1, q, 1, 0.0d0,
1341 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1342 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1344 CALL dsbevx_2stage(
'V',
'A',
'U', 0, 0, a, 1, q, 1, 0.0d0,
1345 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1346 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1348 CALL dsbevx_2stage(
'N',
'/',
'U', 0, 0, a, 1, q, 1, 0.0d0,
1349 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1350 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1352 CALL dsbevx_2stage(
'N',
'A',
'/', 0, 0, a, 1, q, 1, 0.0d0,
1353 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1354 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1356 CALL dsbevx_2stage(
'N',
'A',
'U', -1, 0, a, 1, q, 1, 0.0d0,
1357 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1358 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1360 CALL dsbevx_2stage(
'N',
'A',
'U', 0, -1, a, 1, q, 1, 0.0d0,
1361 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1362 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1364 CALL dsbevx_2stage(
'N',
'A',
'U', 2, 1, a, 1, q, 1, 0.0d0,
1365 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1366 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1372 CALL dsbevx_2stage(
'N',
'V',
'U', 1, 0, a, 1, q, 1, 0.0d0,
1373 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1374 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1376 CALL dsbevx_2stage(
'N',
'I',
'U', 1, 0, a, 1, q, 1, 0.0d0,
1377 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1378 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1380 CALL dsbevx_2stage(
'N',
'I',
'U', 1, 0, a, 1, q, 1, 0.0d0,
1381 $ 0.0d0, 2, 1, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1382 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1384 CALL dsbevx_2stage(
'N',
'I',
'U', 2, 0, a, 1, q, 1, 0.0d0,
1385 $ 0.0d0, 2, 1, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1386 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1388 CALL dsbevx_2stage(
'N',
'I',
'U', 1, 0, a, 1, q, 1, 0.0d0,
1389 $ 0.0d0, 1, 2, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1390 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1396 CALL dsbevx_2stage(
'N',
'A',
'U', 0, 0, a, 1, q, 1, 0.0d0,
1397 $ 0.0d0, 0, 0, 0.0d0, m, x, z, 1, w, 0, iw, i3, info )
1398 CALL chkxer(
'DSBEVX_2STAGE', infot, nout, lerr, ok )
1406 WRITE( nout, fmt = 9999 )path, nt
1408 WRITE( nout, fmt = 9998 )path
1411 9999
FORMAT( 1x, a3,
' routines passed the tests of the error exits',
1412 $
' (', i3,
' tests done)' )
1413 9998
FORMAT(
' *** ', a3,
' routines failed the tests of the error ',
subroutine chkxer(SRNAMT, INFOT, NOUT, LERR, OK)
subroutine dsytrd_sb2st(STAGE1, VECT, UPLO, N, KD, AB, LDAB, D, E, HOUS, LHOUS, WORK, LWORK, INFO)
DSYTRD_SB2ST reduces a real symmetric band matrix A to real symmetric tridiagonal form T
logical function lsamen(N, CA, CB)
LSAMEN
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
subroutine dstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEDC
subroutine dsterf(N, D, E, INFO)
DSTERF
subroutine dsptrd(UPLO, N, AP, D, E, TAU, INFO)
DSPTRD
subroutine dstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEIN
subroutine dopgtr(UPLO, N, AP, TAU, Q, LDQ, WORK, INFO)
DOPGTR
subroutine dsbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
DSBTRD
subroutine dormtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMTR
subroutine dorgtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
DORGTR
subroutine dopmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
DOPMTR
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 matrice...
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 matrice...
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 matrice...
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 matrice...
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 matrice...
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 matrice...
subroutine dsbevd_2stage(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSBEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
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 matrices
subroutine dsbev_2stage(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, LWORK, INFO)
DSBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER m...
subroutine dstev(JOBZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
subroutine dsbevx_2stage(JOBZ, RANGE, UPLO, N, KD, AB, LDAB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO)
DSBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER ...
subroutine dsbev(JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, WORK, INFO)
DSBEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices
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 matrice...
subroutine dpteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DPTEQR
subroutine dsytrd_sy2sb(UPLO, N, KD, A, LDA, AB, LDAB, TAU, WORK, LWORK, INFO)
DSYTRD_SY2SB
subroutine dsytrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
DSYTRD_2STAGE
subroutine dsytrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
DSYTRD
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
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
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
subroutine dsyevr_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSYEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...
subroutine dsyevx_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO)
DSYEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...
subroutine dsyev_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO)
DSYEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matr...
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
subroutine dsyevd_2stage(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO)
DSYEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...