453
454
455
456
457
458
459 INTEGER INFO, LDA, LDU, LIWORK, LWORK, NOUNIT, NSIZES,
460 $ NTYPES
461 DOUBLE PRECISION THRESH
462
463
464 LOGICAL DOTYPE( * )
465 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
466 DOUBLE PRECISION A( LDA, * ), D1( * ), D2( * ), D3( * ),
467 $ D4( * ), EVEIGS( * ), RESULT( * ), TAU( * ),
468 $ U( LDU, * ), V( LDU, * ), WA1( * ), WA2( * ),
469 $ WA3( * ), WORK( * ), Z( LDU, * )
470
471
472
473
474
475 DOUBLE PRECISION ZERO, ONE, TWO, TEN
476 parameter( zero = 0.0d0, one = 1.0d0, two = 2.0d0,
477 $ ten = 10.0d0 )
478 DOUBLE PRECISION HALF
479 parameter( half = 0.5d0 )
480 INTEGER MAXTYP
481 parameter( maxtyp = 18 )
482
483
484 LOGICAL BADNN
485 CHARACTER UPLO
486 INTEGER I, IDIAG, IHBW, IINFO, IL, IMODE, INDX, IROW,
487 $ ITEMP, ITYPE, IU, IUPLO, J, J1, J2, JCOL,
488 $ JSIZE, JTYPE, KD, LGN, LIWEDC, LWEDC, M, M2,
489 $ M3, MTYPES, N, NERRS, NMATS, NMAX, NTEST,
490 $ NTESTT
491 DOUBLE PRECISION ABSTOL, ANINV, ANORM, COND, OVFL, RTOVFL,
492 $ RTUNFL, TEMP1, TEMP2, TEMP3, ULP, ULPINV, UNFL,
493 $ VL, VU
494
495
496 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), ISEED2( 4 ),
497 $ ISEED3( 4 ), KMAGN( MAXTYP ), KMODE( MAXTYP ),
498 $ KTYPE( MAXTYP )
499
500
501 DOUBLE PRECISION DLAMCH, DLARND, DSXT1
503
504
513
514
515 CHARACTER*32 SRNAMT
516
517
518 COMMON / srnamc / srnamt
519
520
521 INTRINSIC abs, dble, int, log, max, min, sqrt
522
523
524 DATA ktype / 1, 2, 5*4, 5*5, 3*8, 3*9 /
525 DATA kmagn / 2*1, 1, 1, 1, 2, 3, 1, 1, 1, 2, 3, 1,
526 $ 2, 3, 1, 2, 3 /
527 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
528 $ 0, 0, 4, 4, 4 /
529
530
531
532
533
534 vl = zero
535 vu = zero
536
537
538
539 ntestt = 0
540 info = 0
541
542 badnn = .false.
543 nmax = 1
544 DO 10 j = 1, nsizes
545 nmax = max( nmax, nn( j ) )
546 IF( nn( j ).LT.0 )
547 $ badnn = .true.
548 10 CONTINUE
549
550
551
552 IF( nsizes.LT.0 ) THEN
553 info = -1
554 ELSE IF( badnn ) THEN
555 info = -2
556 ELSE IF( ntypes.LT.0 ) THEN
557 info = -3
558 ELSE IF( lda.LT.nmax ) THEN
559 info = -9
560 ELSE IF( ldu.LT.nmax ) THEN
561 info = -16
562 ELSE IF( 2*max( 2, nmax )**2.GT.lwork ) THEN
563 info = -21
564 END IF
565
566 IF( info.NE.0 ) THEN
567 CALL xerbla(
'DDRVST2STG', -info )
568 RETURN
569 END IF
570
571
572
573 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
574 $ RETURN
575
576
577
578 unfl =
dlamch(
'Safe minimum' )
579 ovfl =
dlamch(
'Overflow' )
581 ulpinv = one / ulp
582 rtunfl = sqrt( unfl )
583 rtovfl = sqrt( ovfl )
584
585
586
587 DO 20 i = 1, 4
588 iseed2( i ) = iseed( i )
589 iseed3( i ) = iseed( i )
590 20 CONTINUE
591
592 nerrs = 0
593 nmats = 0
594
595
596 DO 1740 jsize = 1, nsizes
597 n = nn( jsize )
598 IF( n.GT.0 ) THEN
599 lgn = int( log( dble( n ) ) / log( two ) )
600 IF( 2**lgn.LT.n )
601 $ lgn = lgn + 1
602 IF( 2**lgn.LT.n )
603 $ lgn = lgn + 1
604 lwedc = 1 + 4*n + 2*n*lgn + 4*n**2
605
606 liwedc = 3 + 5*n
607 ELSE
608 lwedc = 9
609
610 liwedc = 8
611 END IF
612 aninv = one / dble( max( 1, n ) )
613
614 IF( nsizes.NE.1 ) THEN
615 mtypes = min( maxtyp, ntypes )
616 ELSE
617 mtypes = min( maxtyp+1, ntypes )
618 END IF
619
620 DO 1730 jtype = 1, mtypes
621
622 IF( .NOT.dotype( jtype ) )
623 $ GO TO 1730
624 nmats = nmats + 1
625 ntest = 0
626
627 DO 30 j = 1, 4
628 ioldsd( j ) = iseed( j )
629 30 CONTINUE
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646 IF( mtypes.GT.maxtyp )
647 $ GO TO 110
648
649 itype = ktype( jtype )
650 imode = kmode( jtype )
651
652
653
654 GO TO ( 40, 50, 60 )kmagn( jtype )
655
656 40 CONTINUE
657 anorm = one
658 GO TO 70
659
660 50 CONTINUE
661 anorm = ( rtovfl*ulp )*aninv
662 GO TO 70
663
664 60 CONTINUE
665 anorm = rtunfl*n*ulpinv
666 GO TO 70
667
668 70 CONTINUE
669
670 CALL dlaset(
'Full', lda, n, zero, zero, a, lda )
671 iinfo = 0
672 cond = ulpinv
673
674
675
676
677
678 IF( itype.EQ.1 ) THEN
679 iinfo = 0
680
681 ELSE IF( itype.EQ.2 ) THEN
682
683
684
685 DO 80 jcol = 1, n
686 a( jcol, jcol ) = anorm
687 80 CONTINUE
688
689 ELSE IF( itype.EQ.4 ) THEN
690
691
692
693 CALL dlatms( n, n,
'S', iseed,
'S', work, imode, cond,
694 $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
695 $ iinfo )
696
697 ELSE IF( itype.EQ.5 ) THEN
698
699
700
701 CALL dlatms( n, n,
'S', iseed,
'S', work, imode, cond,
702 $ anorm, n, n, 'N', a, lda, work( n+1 ),
703 $ iinfo )
704
705 ELSE IF( itype.EQ.7 ) THEN
706
707
708
709 idumma( 1 ) = 1
710 CALL dlatmr( n, n,
'S', iseed,
'S', work, 6, one, one,
711 $ 'T', 'N', work( n+1 ), 1, one,
712 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
713 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
714
715 ELSE IF( itype.EQ.8 ) THEN
716
717
718
719 idumma( 1 ) = 1
720 CALL dlatmr( n, n,
'S', iseed,
'S', work, 6, one, one,
721 $ 'T', 'N', work( n+1 ), 1, one,
722 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
723 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
724
725 ELSE IF( itype.EQ.9 ) THEN
726
727
728
729 ihbw = int( ( n-1 )*
dlarnd( 1, iseed3 ) )
730 CALL dlatms( n, n,
'S', iseed,
'S', work, imode, cond,
731 $ anorm, ihbw, ihbw, 'Z', u, ldu, work( n+1 ),
732 $ iinfo )
733
734
735
736 CALL dlaset(
'Full', lda, n, zero, zero, a, lda )
737 DO 100 idiag = -ihbw, ihbw
738 irow = ihbw - idiag + 1
739 j1 = max( 1, idiag+1 )
740 j2 = min( n, n+idiag )
741 DO 90 j = j1, j2
742 i = j - idiag
743 a( i, j ) = u( irow, j )
744 90 CONTINUE
745 100 CONTINUE
746 ELSE
747 iinfo = 1
748 END IF
749
750 IF( iinfo.NE.0 ) THEN
751 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
752 $ ioldsd
753 info = abs( iinfo )
754 RETURN
755 END IF
756
757 110 CONTINUE
758
759 abstol = unfl + unfl
760 IF( n.LE.1 ) THEN
761 il = 1
762 iu = n
763 ELSE
764 il = 1 + int( ( n-1 )*
dlarnd( 1, iseed2 ) )
765 iu = 1 + int( ( n-1 )*
dlarnd( 1, iseed2 ) )
766 IF( il.GT.iu ) THEN
767 itemp = il
768 il = iu
769 iu = itemp
770 END IF
771 END IF
772
773
774
775 IF( jtype.LE.7 ) THEN
776 ntest = 1
777 DO 120 i = 1, n
778 d1( i ) = dble( a( i, i ) )
779 120 CONTINUE
780 DO 130 i = 1, n - 1
781 d2( i ) = dble( a( i+1, i ) )
782 130 CONTINUE
783 srnamt = 'DSTEV'
784 CALL dstev(
'V', n, d1, d2, z, ldu, work, iinfo )
785 IF( iinfo.NE.0 ) THEN
786 WRITE( nounit, fmt = 9999 )'DSTEV(V)', iinfo, n,
787 $ jtype, ioldsd
788 info = abs( iinfo )
789 IF( iinfo.LT.0 ) THEN
790 RETURN
791 ELSE
792 result( 1 ) = ulpinv
793 result( 2 ) = ulpinv
794 result( 3 ) = ulpinv
795 GO TO 180
796 END IF
797 END IF
798
799
800
801 DO 140 i = 1, n
802 d3( i ) = dble( a( i, i ) )
803 140 CONTINUE
804 DO 150 i = 1, n - 1
805 d4( i ) = dble( a( i+1, i ) )
806 150 CONTINUE
807 CALL dstt21( n, 0, d3, d4, d1, d2, z, ldu, work,
808 $ result( 1 ) )
809
810 ntest = 3
811 DO 160 i = 1, n - 1
812 d4( i ) = dble( a( i+1, i ) )
813 160 CONTINUE
814 srnamt = 'DSTEV'
815 CALL dstev(
'N', n, d3, d4, z, ldu, work, iinfo )
816 IF( iinfo.NE.0 ) THEN
817 WRITE( nounit, fmt = 9999 )'DSTEV(N)', iinfo, n,
818 $ jtype, ioldsd
819 info = abs( iinfo )
820 IF( iinfo.LT.0 ) THEN
821 RETURN
822 ELSE
823 result( 3 ) = ulpinv
824 GO TO 180
825 END IF
826 END IF
827
828
829
830 temp1 = zero
831 temp2 = zero
832 DO 170 j = 1, n
833 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
834 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
835 170 CONTINUE
836 result( 3 ) = temp2 / max( unfl,
837 $ ulp*max( temp1, temp2 ) )
838
839 180 CONTINUE
840
841 ntest = 4
842 DO 190 i = 1, n
843 eveigs( i ) = d3( i )
844 d1( i ) = dble( a( i, i ) )
845 190 CONTINUE
846 DO 200 i = 1, n - 1
847 d2( i ) = dble( a( i+1, i ) )
848 200 CONTINUE
849 srnamt = 'DSTEVX'
850 CALL dstevx(
'V',
'A', n, d1, d2, vl, vu, il, iu, abstol,
851 $ m, wa1, z, ldu, work, iwork, iwork( 5*n+1 ),
852 $ iinfo )
853 IF( iinfo.NE.0 ) THEN
854 WRITE( nounit, fmt = 9999 )'DSTEVX(V,A)', iinfo, n,
855 $ jtype, ioldsd
856 info = abs( iinfo )
857 IF( iinfo.LT.0 ) THEN
858 RETURN
859 ELSE
860 result( 4 ) = ulpinv
861 result( 5 ) = ulpinv
862 result( 6 ) = ulpinv
863 GO TO 250
864 END IF
865 END IF
866 IF( n.GT.0 ) THEN
867 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
868 ELSE
869 temp3 = zero
870 END IF
871
872
873
874 DO 210 i = 1, n
875 d3( i ) = dble( a( i, i ) )
876 210 CONTINUE
877 DO 220 i = 1, n - 1
878 d4( i ) = dble( a( i+1, i ) )
879 220 CONTINUE
880 CALL dstt21( n, 0, d3, d4, wa1, d2, z, ldu, work,
881 $ result( 4 ) )
882
883 ntest = 6
884 DO 230 i = 1, n - 1
885 d4( i ) = dble( a( i+1, i ) )
886 230 CONTINUE
887 srnamt = 'DSTEVX'
888 CALL dstevx(
'N',
'A', n, d3, d4, vl, vu, il, iu, abstol,
889 $ m2, wa2, z, ldu, work, iwork,
890 $ iwork( 5*n+1 ), iinfo )
891 IF( iinfo.NE.0 ) THEN
892 WRITE( nounit, fmt = 9999 )'DSTEVX(N,A)', iinfo, n,
893 $ jtype, ioldsd
894 info = abs( iinfo )
895 IF( iinfo.LT.0 ) THEN
896 RETURN
897 ELSE
898 result( 6 ) = ulpinv
899 GO TO 250
900 END IF
901 END IF
902
903
904
905 temp1 = zero
906 temp2 = zero
907 DO 240 j = 1, n
908 temp1 = max( temp1, abs( wa2( j ) ),
909 $ abs( eveigs( j ) ) )
910 temp2 = max( temp2, abs( wa2( j )-eveigs( j ) ) )
911 240 CONTINUE
912 result( 6 ) = temp2 / max( unfl,
913 $ ulp*max( temp1, temp2 ) )
914
915 250 CONTINUE
916
917 ntest = 7
918 DO 260 i = 1, n
919 d1( i ) = dble( a( i, i ) )
920 260 CONTINUE
921 DO 270 i = 1, n - 1
922 d2( i ) = dble( a( i+1, i ) )
923 270 CONTINUE
924 srnamt = 'DSTEVR'
925 CALL dstevr(
'V',
'A', n, d1, d2, vl, vu, il, iu, abstol,
926 $ m, wa1, z, ldu, iwork, work, lwork,
927 $ iwork(2*n+1), liwork-2*n, iinfo )
928 IF( iinfo.NE.0 ) THEN
929 WRITE( nounit, fmt = 9999 )'DSTEVR(V,A)', iinfo, n,
930 $ jtype, ioldsd
931 info = abs( iinfo )
932 IF( iinfo.LT.0 ) THEN
933 RETURN
934 ELSE
935 result( 7 ) = ulpinv
936 result( 8 ) = ulpinv
937 GO TO 320
938 END IF
939 END IF
940 IF( n.GT.0 ) THEN
941 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
942 ELSE
943 temp3 = zero
944 END IF
945
946
947
948 DO 280 i = 1, n
949 d3( i ) = dble( a( i, i ) )
950 280 CONTINUE
951 DO 290 i = 1, n - 1
952 d4( i ) = dble( a( i+1, i ) )
953 290 CONTINUE
954 CALL dstt21( n, 0, d3, d4, wa1, d2, z, ldu, work,
955 $ result( 7 ) )
956
957 ntest = 9
958 DO 300 i = 1, n - 1
959 d4( i ) = dble( a( i+1, i ) )
960 300 CONTINUE
961 srnamt = 'DSTEVR'
962 CALL dstevr(
'N',
'A', n, d3, d4, vl, vu, il, iu, abstol,
963 $ m2, wa2, z, ldu, iwork, work, lwork,
964 $ iwork(2*n+1), liwork-2*n, iinfo )
965 IF( iinfo.NE.0 ) THEN
966 WRITE( nounit, fmt = 9999 )'DSTEVR(N,A)', iinfo, n,
967 $ jtype, ioldsd
968 info = abs( iinfo )
969 IF( iinfo.LT.0 ) THEN
970 RETURN
971 ELSE
972 result( 9 ) = ulpinv
973 GO TO 320
974 END IF
975 END IF
976
977
978
979 temp1 = zero
980 temp2 = zero
981 DO 310 j = 1, n
982 temp1 = max( temp1, abs( wa2( j ) ),
983 $ abs( eveigs( j ) ) )
984 temp2 = max( temp2, abs( wa2( j )-eveigs( j ) ) )
985 310 CONTINUE
986 result( 9 ) = temp2 / max( unfl,
987 $ ulp*max( temp1, temp2 ) )
988
989 320 CONTINUE
990
991
992 ntest = 10
993 DO 330 i = 1, n
994 d1( i ) = dble( a( i, i ) )
995 330 CONTINUE
996 DO 340 i = 1, n - 1
997 d2( i ) = dble( a( i+1, i ) )
998 340 CONTINUE
999 srnamt = 'DSTEVX'
1000 CALL dstevx(
'V',
'I', n, d1, d2, vl, vu, il, iu, abstol,
1001 $ m2, wa2, z, ldu, work, iwork,
1002 $ iwork( 5*n+1 ), iinfo )
1003 IF( iinfo.NE.0 ) THEN
1004 WRITE( nounit, fmt = 9999 )'DSTEVX(V,I)', iinfo, n,
1005 $ jtype, ioldsd
1006 info = abs( iinfo )
1007 IF( iinfo.LT.0 ) THEN
1008 RETURN
1009 ELSE
1010 result( 10 ) = ulpinv
1011 result( 11 ) = ulpinv
1012 result( 12 ) = ulpinv
1013 GO TO 380
1014 END IF
1015 END IF
1016
1017
1018
1019 DO 350 i = 1, n
1020 d3( i ) = dble( a( i, i ) )
1021 350 CONTINUE
1022 DO 360 i = 1, n - 1
1023 d4( i ) = dble( a( i+1, i ) )
1024 360 CONTINUE
1025 CALL dstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1026 $ max( 1, m2 ), result( 10 ) )
1027
1028
1029 ntest = 12
1030 DO 370 i = 1, n - 1
1031 d4( i ) = dble( a( i+1, i ) )
1032 370 CONTINUE
1033 srnamt = 'DSTEVX'
1034 CALL dstevx(
'N',
'I', n, d3, d4, vl, vu, il, iu, abstol,
1035 $ m3, wa3, z, ldu, work, iwork,
1036 $ iwork( 5*n+1 ), iinfo )
1037 IF( iinfo.NE.0 ) THEN
1038 WRITE( nounit, fmt = 9999 )'DSTEVX(N,I)', iinfo, n,
1039 $ jtype, ioldsd
1040 info = abs( iinfo )
1041 IF( iinfo.LT.0 ) THEN
1042 RETURN
1043 ELSE
1044 result( 12 ) = ulpinv
1045 GO TO 380
1046 END IF
1047 END IF
1048
1049
1050
1051 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1052 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1053 result( 12 ) = ( temp1+temp2 ) / max( unfl, ulp*temp3 )
1054
1055 380 CONTINUE
1056
1057 ntest = 12
1058 IF( n.GT.0 ) THEN
1059 IF( il.NE.1 ) THEN
1060 vl = wa1( il ) - max( half*
1061 $ ( wa1( il )-wa1( il-1 ) ), ten*ulp*temp3,
1062 $ ten*rtunfl )
1063 ELSE
1064 vl = wa1( 1 ) - max( half*( wa1( n )-wa1( 1 ) ),
1065 $ ten*ulp*temp3, ten*rtunfl )
1066 END IF
1067 IF( iu.NE.n ) THEN
1068 vu = wa1( iu ) + max( half*
1069 $ ( wa1( iu+1 )-wa1( iu ) ), ten*ulp*temp3,
1070 $ ten*rtunfl )
1071 ELSE
1072 vu = wa1( n ) + max( half*( wa1( n )-wa1( 1 ) ),
1073 $ ten*ulp*temp3, ten*rtunfl )
1074 END IF
1075 ELSE
1076 vl = zero
1077 vu = one
1078 END IF
1079
1080 DO 390 i = 1, n
1081 d1( i ) = dble( a( i, i ) )
1082 390 CONTINUE
1083 DO 400 i = 1, n - 1
1084 d2( i ) = dble( a( i+1, i ) )
1085 400 CONTINUE
1086 srnamt = 'DSTEVX'
1087 CALL dstevx(
'V',
'V', n, d1, d2, vl, vu, il, iu, abstol,
1088 $ m2, wa2, z, ldu, work, iwork,
1089 $ iwork( 5*n+1 ), iinfo )
1090 IF( iinfo.NE.0 ) THEN
1091 WRITE( nounit, fmt = 9999 )'DSTEVX(V,V)', iinfo, n,
1092 $ jtype, ioldsd
1093 info = abs( iinfo )
1094 IF( iinfo.LT.0 ) THEN
1095 RETURN
1096 ELSE
1097 result( 13 ) = ulpinv
1098 result( 14 ) = ulpinv
1099 result( 15 ) = ulpinv
1100 GO TO 440
1101 END IF
1102 END IF
1103
1104 IF( m2.EQ.0 .AND. n.GT.0 ) THEN
1105 result( 13 ) = ulpinv
1106 result( 14 ) = ulpinv
1107 result( 15 ) = ulpinv
1108 GO TO 440
1109 END IF
1110
1111
1112
1113 DO 410 i = 1, n
1114 d3( i ) = dble( a( i, i ) )
1115 410 CONTINUE
1116 DO 420 i = 1, n - 1
1117 d4( i ) = dble( a( i+1, i ) )
1118 420 CONTINUE
1119 CALL dstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1120 $ max( 1, m2 ), result( 13 ) )
1121
1122 ntest = 15
1123 DO 430 i = 1, n - 1
1124 d4( i ) = dble( a( i+1, i ) )
1125 430 CONTINUE
1126 srnamt = 'DSTEVX'
1127 CALL dstevx(
'N',
'V', n, d3, d4, vl, vu, il, iu, abstol,
1128 $ m3, wa3, z, ldu, work, iwork,
1129 $ iwork( 5*n+1 ), iinfo )
1130 IF( iinfo.NE.0 ) THEN
1131 WRITE( nounit, fmt = 9999 )'DSTEVX(N,V)', iinfo, n,
1132 $ jtype, ioldsd
1133 info = abs( iinfo )
1134 IF( iinfo.LT.0 ) THEN
1135 RETURN
1136 ELSE
1137 result( 15 ) = ulpinv
1138 GO TO 440
1139 END IF
1140 END IF
1141
1142
1143
1144 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1145 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1146 result( 15 ) = ( temp1+temp2 ) / max( unfl, temp3*ulp )
1147
1148 440 CONTINUE
1149
1150 ntest = 16
1151 DO 450 i = 1, n
1152 d1( i ) = dble( a( i, i ) )
1153 450 CONTINUE
1154 DO 460 i = 1, n - 1
1155 d2( i ) = dble( a( i+1, i ) )
1156 460 CONTINUE
1157 srnamt = 'DSTEVD'
1158 CALL dstevd(
'V', n, d1, d2, z, ldu, work, lwedc, iwork,
1159 $ liwedc, iinfo )
1160 IF( iinfo.NE.0 ) THEN
1161 WRITE( nounit, fmt = 9999 )'DSTEVD(V)', iinfo, n,
1162 $ jtype, ioldsd
1163 info = abs( iinfo )
1164 IF( iinfo.LT.0 ) THEN
1165 RETURN
1166 ELSE
1167 result( 16 ) = ulpinv
1168 result( 17 ) = ulpinv
1169 result( 18 ) = ulpinv
1170 GO TO 510
1171 END IF
1172 END IF
1173
1174
1175
1176 DO 470 i = 1, n
1177 d3( i ) = dble( a( i, i ) )
1178 470 CONTINUE
1179 DO 480 i = 1, n - 1
1180 d4( i ) = dble( a( i+1, i ) )
1181 480 CONTINUE
1182 CALL dstt21( n, 0, d3, d4, d1, d2, z, ldu, work,
1183 $ result( 16 ) )
1184
1185 ntest = 18
1186 DO 490 i = 1, n - 1
1187 d4( i ) = dble( a( i+1, i ) )
1188 490 CONTINUE
1189 srnamt = 'DSTEVD'
1190 CALL dstevd(
'N', n, d3, d4, z, ldu, work, lwedc, iwork,
1191 $ liwedc, iinfo )
1192 IF( iinfo.NE.0 ) THEN
1193 WRITE( nounit, fmt = 9999 )'DSTEVD(N)', iinfo, n,
1194 $ jtype, ioldsd
1195 info = abs( iinfo )
1196 IF( iinfo.LT.0 ) THEN
1197 RETURN
1198 ELSE
1199 result( 18 ) = ulpinv
1200 GO TO 510
1201 END IF
1202 END IF
1203
1204
1205
1206 temp1 = zero
1207 temp2 = zero
1208 DO 500 j = 1, n
1209 temp1 = max( temp1, abs( eveigs( j ) ),
1210 $ abs( d3( j ) ) )
1211 temp2 = max( temp2, abs( eveigs( j )-d3( j ) ) )
1212 500 CONTINUE
1213 result( 18 ) = temp2 / max( unfl,
1214 $ ulp*max( temp1, temp2 ) )
1215
1216 510 CONTINUE
1217
1218 ntest = 19
1219 DO 520 i = 1, n
1220 d1( i ) = dble( a( i, i ) )
1221 520 CONTINUE
1222 DO 530 i = 1, n - 1
1223 d2( i ) = dble( a( i+1, i ) )
1224 530 CONTINUE
1225 srnamt = 'DSTEVR'
1226 CALL dstevr(
'V',
'I', n, d1, d2, vl, vu, il, iu, abstol,
1227 $ m2, wa2, z, ldu, iwork, work, lwork,
1228 $ iwork(2*n+1), liwork-2*n, iinfo )
1229 IF( iinfo.NE.0 ) THEN
1230 WRITE( nounit, fmt = 9999 )'DSTEVR(V,I)', iinfo, n,
1231 $ jtype, ioldsd
1232 info = abs( iinfo )
1233 IF( iinfo.LT.0 ) THEN
1234 RETURN
1235 ELSE
1236 result( 19 ) = ulpinv
1237 result( 20 ) = ulpinv
1238 result( 21 ) = ulpinv
1239 GO TO 570
1240 END IF
1241 END IF
1242
1243
1244
1245 DO 540 i = 1, n
1246 d3( i ) = dble( a( i, i ) )
1247 540 CONTINUE
1248 DO 550 i = 1, n - 1
1249 d4( i ) = dble( a( i+1, i ) )
1250 550 CONTINUE
1251 CALL dstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1252 $ max( 1, m2 ), result( 19 ) )
1253
1254
1255 ntest = 21
1256 DO 560 i = 1, n - 1
1257 d4( i ) = dble( a( i+1, i ) )
1258 560 CONTINUE
1259 srnamt = 'DSTEVR'
1260 CALL dstevr(
'N',
'I', n, d3, d4, vl, vu, il, iu, abstol,
1261 $ m3, wa3, z, ldu, iwork, work, lwork,
1262 $ iwork(2*n+1), liwork-2*n, iinfo )
1263 IF( iinfo.NE.0 ) THEN
1264 WRITE( nounit, fmt = 9999 )'DSTEVR(N,I)', iinfo, n,
1265 $ jtype, ioldsd
1266 info = abs( iinfo )
1267 IF( iinfo.LT.0 ) THEN
1268 RETURN
1269 ELSE
1270 result( 21 ) = ulpinv
1271 GO TO 570
1272 END IF
1273 END IF
1274
1275
1276
1277 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1278 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1279 result( 21 ) = ( temp1+temp2 ) / max( unfl, ulp*temp3 )
1280
1281 570 CONTINUE
1282
1283 ntest = 21
1284 IF( n.GT.0 ) THEN
1285 IF( il.NE.1 ) THEN
1286 vl = wa1( il ) - max( half*
1287 $ ( wa1( il )-wa1( il-1 ) ), ten*ulp*temp3,
1288 $ ten*rtunfl )
1289 ELSE
1290 vl = wa1( 1 ) - max( half*( wa1( n )-wa1( 1 ) ),
1291 $ ten*ulp*temp3, ten*rtunfl )
1292 END IF
1293 IF( iu.NE.n ) THEN
1294 vu = wa1( iu ) + max( half*
1295 $ ( wa1( iu+1 )-wa1( iu ) ), ten*ulp*temp3,
1296 $ ten*rtunfl )
1297 ELSE
1298 vu = wa1( n ) + max( half*( wa1( n )-wa1( 1 ) ),
1299 $ ten*ulp*temp3, ten*rtunfl )
1300 END IF
1301 ELSE
1302 vl = zero
1303 vu = one
1304 END IF
1305
1306 DO 580 i = 1, n
1307 d1( i ) = dble( a( i, i ) )
1308 580 CONTINUE
1309 DO 590 i = 1, n - 1
1310 d2( i ) = dble( a( i+1, i ) )
1311 590 CONTINUE
1312 srnamt = 'DSTEVR'
1313 CALL dstevr(
'V',
'V', n, d1, d2, vl, vu, il, iu, abstol,
1314 $ m2, wa2, z, ldu, iwork, work, lwork,
1315 $ iwork(2*n+1), liwork-2*n, iinfo )
1316 IF( iinfo.NE.0 ) THEN
1317 WRITE( nounit, fmt = 9999 )'DSTEVR(V,V)', iinfo, n,
1318 $ jtype, ioldsd
1319 info = abs( iinfo )
1320 IF( iinfo.LT.0 ) THEN
1321 RETURN
1322 ELSE
1323 result( 22 ) = ulpinv
1324 result( 23 ) = ulpinv
1325 result( 24 ) = ulpinv
1326 GO TO 630
1327 END IF
1328 END IF
1329
1330 IF( m2.EQ.0 .AND. n.GT.0 ) THEN
1331 result( 22 ) = ulpinv
1332 result( 23 ) = ulpinv
1333 result( 24 ) = ulpinv
1334 GO TO 630
1335 END IF
1336
1337
1338
1339 DO 600 i = 1, n
1340 d3( i ) = dble( a( i, i ) )
1341 600 CONTINUE
1342 DO 610 i = 1, n - 1
1343 d4( i ) = dble( a( i+1, i ) )
1344 610 CONTINUE
1345 CALL dstt22( n, m2, 0, d3, d4, wa2, d2, z, ldu, work,
1346 $ max( 1, m2 ), result( 22 ) )
1347
1348 ntest = 24
1349 DO 620 i = 1, n - 1
1350 d4( i ) = dble( a( i+1, i ) )
1351 620 CONTINUE
1352 srnamt = 'DSTEVR'
1353 CALL dstevr(
'N',
'V', n, d3, d4, vl, vu, il, iu, abstol,
1354 $ m3, wa3, z, ldu, iwork, work, lwork,
1355 $ iwork(2*n+1), liwork-2*n, iinfo )
1356 IF( iinfo.NE.0 ) THEN
1357 WRITE( nounit, fmt = 9999 )'DSTEVR(N,V)', iinfo, n,
1358 $ jtype, ioldsd
1359 info = abs( iinfo )
1360 IF( iinfo.LT.0 ) THEN
1361 RETURN
1362 ELSE
1363 result( 24 ) = ulpinv
1364 GO TO 630
1365 END IF
1366 END IF
1367
1368
1369
1370 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1371 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1372 result( 24 ) = ( temp1+temp2 ) / max( unfl, temp3*ulp )
1373
1374 630 CONTINUE
1375
1376
1377
1378 ELSE
1379
1380 DO 640 i = 1, 24
1381 result( i ) = zero
1382 640 CONTINUE
1383 ntest = 24
1384 END IF
1385
1386
1387
1388
1389 DO 1720 iuplo = 0, 1
1390 IF( iuplo.EQ.0 ) THEN
1391 uplo = 'L'
1392 ELSE
1393 uplo = 'U'
1394 END IF
1395
1396
1397
1398 CALL dlacpy(
' ', n, n, a, lda, v, ldu )
1399
1400 ntest = ntest + 1
1401 srnamt = 'DSYEV'
1402 CALL dsyev(
'V', uplo, n, a, ldu, d1, work, lwork,
1403 $ iinfo )
1404 IF( iinfo.NE.0 ) THEN
1405 WRITE( nounit, fmt = 9999 )'DSYEV(V,' // uplo // ')',
1406 $ iinfo, n, jtype, ioldsd
1407 info = abs( iinfo )
1408 IF( iinfo.LT.0 ) THEN
1409 RETURN
1410 ELSE
1411 result( ntest ) = ulpinv
1412 result( ntest+1 ) = ulpinv
1413 result( ntest+2 ) = ulpinv
1414 GO TO 660
1415 END IF
1416 END IF
1417
1418
1419
1420 CALL dsyt21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
1421 $ ldu, tau, work, result( ntest ) )
1422
1423 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1424
1425 ntest = ntest + 2
1426 srnamt = 'DSYEV_2STAGE'
1427 CALL dsyev_2stage(
'N', uplo, n, a, ldu, d3, work, lwork,
1428 $ iinfo )
1429 IF( iinfo.NE.0 ) THEN
1430 WRITE( nounit, fmt = 9999 )
1431 $ 'DSYEV_2STAGE(N,' // uplo // ')',
1432 $ iinfo, n, jtype, ioldsd
1433 info = abs( iinfo )
1434 IF( iinfo.LT.0 ) THEN
1435 RETURN
1436 ELSE
1437 result( ntest ) = ulpinv
1438 GO TO 660
1439 END IF
1440 END IF
1441
1442
1443
1444 temp1 = zero
1445 temp2 = zero
1446 DO 650 j = 1, n
1447 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1448 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1449 650 CONTINUE
1450 result( ntest ) = temp2 / max( unfl,
1451 $ ulp*max( temp1, temp2 ) )
1452
1453 660 CONTINUE
1454 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1455
1456 ntest = ntest + 1
1457
1458 IF( n.GT.0 ) THEN
1459 temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1460 IF( il.NE.1 ) THEN
1461 vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1462 $ ten*ulp*temp3, ten*rtunfl )
1463 ELSE IF( n.GT.0 ) THEN
1464 vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1465 $ ten*ulp*temp3, ten*rtunfl )
1466 END IF
1467 IF( iu.NE.n ) THEN
1468 vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1469 $ ten*ulp*temp3, ten*rtunfl )
1470 ELSE IF( n.GT.0 ) THEN
1471 vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1472 $ ten*ulp*temp3, ten*rtunfl )
1473 END IF
1474 ELSE
1475 temp3 = zero
1476 vl = zero
1477 vu = one
1478 END IF
1479
1480 srnamt = 'DSYEVX'
1481 CALL dsyevx(
'V',
'A', uplo, n, a, ldu, vl, vu, il, iu,
1482 $ abstol, m, wa1, z, ldu, work, lwork, iwork,
1483 $ iwork( 5*n+1 ), iinfo )
1484 IF( iinfo.NE.0 ) THEN
1485 WRITE( nounit, fmt = 9999 )'DSYEVX(V,A,' // uplo //
1486 $ ')', iinfo, n, jtype, ioldsd
1487 info = abs( iinfo )
1488 IF( iinfo.LT.0 ) THEN
1489 RETURN
1490 ELSE
1491 result( ntest ) = ulpinv
1492 result( ntest+1 ) = ulpinv
1493 result( ntest+2 ) = ulpinv
1494 GO TO 680
1495 END IF
1496 END IF
1497
1498
1499
1500 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1501
1502 CALL dsyt21( 1, uplo, n, 0, a, ldu, d1, d2, z, ldu, v,
1503 $ ldu, tau, work, result( ntest ) )
1504
1505 ntest = ntest + 2
1506 srnamt = 'DSYEVX_2STAGE'
1508 $ il, iu, abstol, m2, wa2, z, ldu, work,
1509 $ lwork, iwork, iwork( 5*n+1 ), iinfo )
1510 IF( iinfo.NE.0 ) THEN
1511 WRITE( nounit, fmt = 9999 )
1512 $ 'DSYEVX_2STAGE(N,A,' // uplo //
1513 $ ')', iinfo, n, jtype, ioldsd
1514 info = abs( iinfo )
1515 IF( iinfo.LT.0 ) THEN
1516 RETURN
1517 ELSE
1518 result( ntest ) = ulpinv
1519 GO TO 680
1520 END IF
1521 END IF
1522
1523
1524
1525 temp1 = zero
1526 temp2 = zero
1527 DO 670 j = 1, n
1528 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1529 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1530 670 CONTINUE
1531 result( ntest ) = temp2 / max( unfl,
1532 $ ulp*max( temp1, temp2 ) )
1533
1534 680 CONTINUE
1535
1536 ntest = ntest + 1
1537 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1538 srnamt = 'DSYEVX'
1539 CALL dsyevx(
'V',
'I', uplo, n, a, ldu, vl, vu, il, iu,
1540 $ abstol, m2, wa2, z, ldu, work, lwork, iwork,
1541 $ iwork( 5*n+1 ), iinfo )
1542 IF( iinfo.NE.0 ) THEN
1543 WRITE( nounit, fmt = 9999 )'DSYEVX(V,I,' // uplo //
1544 $ ')', iinfo, n, jtype, ioldsd
1545 info = abs( iinfo )
1546 IF( iinfo.LT.0 ) THEN
1547 RETURN
1548 ELSE
1549 result( ntest ) = ulpinv
1550 result( ntest+1 ) = ulpinv
1551 result( ntest+2 ) = ulpinv
1552 GO TO 690
1553 END IF
1554 END IF
1555
1556
1557
1558 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1559
1560 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1561 $ v, ldu, tau, work, result( ntest ) )
1562
1563 ntest = ntest + 2
1564 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1565 srnamt = 'DSYEVX_2STAGE'
1567 $ il, iu, abstol, m3, wa3, z, ldu, work,
1568 $ lwork, iwork, iwork( 5*n+1 ), iinfo )
1569 IF( iinfo.NE.0 ) THEN
1570 WRITE( nounit, fmt = 9999 )
1571 $ 'DSYEVX_2STAGE(N,I,' // uplo //
1572 $ ')', iinfo, n, jtype, ioldsd
1573 info = abs( iinfo )
1574 IF( iinfo.LT.0 ) THEN
1575 RETURN
1576 ELSE
1577 result( ntest ) = ulpinv
1578 GO TO 690
1579 END IF
1580 END IF
1581
1582
1583
1584 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1585 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1586 result( ntest ) = ( temp1+temp2 ) /
1587 $ max( unfl, ulp*temp3 )
1588 690 CONTINUE
1589
1590 ntest = ntest + 1
1591 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1592 srnamt = 'DSYEVX'
1593 CALL dsyevx(
'V',
'V', uplo, n, a, ldu, vl, vu, il, iu,
1594 $ abstol, m2, wa2, z, ldu, work, lwork, iwork,
1595 $ iwork( 5*n+1 ), iinfo )
1596 IF( iinfo.NE.0 ) THEN
1597 WRITE( nounit, fmt = 9999 )'DSYEVX(V,V,' // uplo //
1598 $ ')', iinfo, n, jtype, ioldsd
1599 info = abs( iinfo )
1600 IF( iinfo.LT.0 ) THEN
1601 RETURN
1602 ELSE
1603 result( ntest ) = ulpinv
1604 result( ntest+1 ) = ulpinv
1605 result( ntest+2 ) = ulpinv
1606 GO TO 700
1607 END IF
1608 END IF
1609
1610
1611
1612 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1613
1614 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1615 $ v, ldu, tau, work, result( ntest ) )
1616
1617 ntest = ntest + 2
1618 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1619 srnamt = 'DSYEVX_2STAGE'
1621 $ il, iu, abstol, m3, wa3, z, ldu, work,
1622 $ lwork, iwork, iwork( 5*n+1 ), iinfo )
1623 IF( iinfo.NE.0 ) THEN
1624 WRITE( nounit, fmt = 9999 )
1625 $ 'DSYEVX_2STAGE(N,V,' // uplo //
1626 $ ')', iinfo, n, jtype, ioldsd
1627 info = abs( iinfo )
1628 IF( iinfo.LT.0 ) THEN
1629 RETURN
1630 ELSE
1631 result( ntest ) = ulpinv
1632 GO TO 700
1633 END IF
1634 END IF
1635
1636 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1637 result( ntest ) = ulpinv
1638 GO TO 700
1639 END IF
1640
1641
1642
1643 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1644 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1645 IF( n.GT.0 ) THEN
1646 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1647 ELSE
1648 temp3 = zero
1649 END IF
1650 result( ntest ) = ( temp1+temp2 ) /
1651 $ max( unfl, temp3*ulp )
1652
1653 700 CONTINUE
1654
1655
1656
1657 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
1658
1659
1660
1661
1662 IF( iuplo.EQ.1 ) THEN
1663 indx = 1
1664 DO 720 j = 1, n
1665 DO 710 i = 1, j
1666 work( indx ) = a( i, j )
1667 indx = indx + 1
1668 710 CONTINUE
1669 720 CONTINUE
1670 ELSE
1671 indx = 1
1672 DO 740 j = 1, n
1673 DO 730 i = j, n
1674 work( indx ) = a( i, j )
1675 indx = indx + 1
1676 730 CONTINUE
1677 740 CONTINUE
1678 END IF
1679
1680 ntest = ntest + 1
1681 srnamt = 'DSPEV'
1682 CALL dspev(
'V', uplo, n, work, d1, z, ldu, v, iinfo )
1683 IF( iinfo.NE.0 ) THEN
1684 WRITE( nounit, fmt = 9999 )'DSPEV(V,' // uplo // ')',
1685 $ iinfo, n, jtype, ioldsd
1686 info = abs( iinfo )
1687 IF( iinfo.LT.0 ) THEN
1688 RETURN
1689 ELSE
1690 result( ntest ) = ulpinv
1691 result( ntest+1 ) = ulpinv
1692 result( ntest+2 ) = ulpinv
1693 GO TO 800
1694 END IF
1695 END IF
1696
1697
1698
1699 CALL dsyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
1700 $ ldu, tau, work, result( ntest ) )
1701
1702 IF( iuplo.EQ.1 ) THEN
1703 indx = 1
1704 DO 760 j = 1, n
1705 DO 750 i = 1, j
1706 work( indx ) = a( i, j )
1707 indx = indx + 1
1708 750 CONTINUE
1709 760 CONTINUE
1710 ELSE
1711 indx = 1
1712 DO 780 j = 1, n
1713 DO 770 i = j, n
1714 work( indx ) = a( i, j )
1715 indx = indx + 1
1716 770 CONTINUE
1717 780 CONTINUE
1718 END IF
1719
1720 ntest = ntest + 2
1721 srnamt = 'DSPEV'
1722 CALL dspev(
'N', uplo, n, work, d3, z, ldu, v, iinfo )
1723 IF( iinfo.NE.0 ) THEN
1724 WRITE( nounit, fmt = 9999 )'DSPEV(N,' // uplo // ')',
1725 $ iinfo, n, jtype, ioldsd
1726 info = abs( iinfo )
1727 IF( iinfo.LT.0 ) THEN
1728 RETURN
1729 ELSE
1730 result( ntest ) = ulpinv
1731 GO TO 800
1732 END IF
1733 END IF
1734
1735
1736
1737 temp1 = zero
1738 temp2 = zero
1739 DO 790 j = 1, n
1740 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
1741 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
1742 790 CONTINUE
1743 result( ntest ) = temp2 / max( unfl,
1744 $ ulp*max( temp1, temp2 ) )
1745
1746
1747
1748
1749 800 CONTINUE
1750 IF( iuplo.EQ.1 ) THEN
1751 indx = 1
1752 DO 820 j = 1, n
1753 DO 810 i = 1, j
1754 work( indx ) = a( i, j )
1755 indx = indx + 1
1756 810 CONTINUE
1757 820 CONTINUE
1758 ELSE
1759 indx = 1
1760 DO 840 j = 1, n
1761 DO 830 i = j, n
1762 work( indx ) = a( i, j )
1763 indx = indx + 1
1764 830 CONTINUE
1765 840 CONTINUE
1766 END IF
1767
1768 ntest = ntest + 1
1769
1770 IF( n.GT.0 ) THEN
1771 temp3 = max( abs( d1( 1 ) ), abs( d1( n ) ) )
1772 IF( il.NE.1 ) THEN
1773 vl = d1( il ) - max( half*( d1( il )-d1( il-1 ) ),
1774 $ ten*ulp*temp3, ten*rtunfl )
1775 ELSE IF( n.GT.0 ) THEN
1776 vl = d1( 1 ) - max( half*( d1( n )-d1( 1 ) ),
1777 $ ten*ulp*temp3, ten*rtunfl )
1778 END IF
1779 IF( iu.NE.n ) THEN
1780 vu = d1( iu ) + max( half*( d1( iu+1 )-d1( iu ) ),
1781 $ ten*ulp*temp3, ten*rtunfl )
1782 ELSE IF( n.GT.0 ) THEN
1783 vu = d1( n ) + max( half*( d1( n )-d1( 1 ) ),
1784 $ ten*ulp*temp3, ten*rtunfl )
1785 END IF
1786 ELSE
1787 temp3 = zero
1788 vl = zero
1789 vu = one
1790 END IF
1791
1792 srnamt = 'DSPEVX'
1793 CALL dspevx(
'V',
'A', uplo, n, work, vl, vu, il, iu,
1794 $ abstol, m, wa1, z, ldu, v, iwork,
1795 $ iwork( 5*n+1 ), iinfo )
1796 IF( iinfo.NE.0 ) THEN
1797 WRITE( nounit, fmt = 9999 )'DSPEVX(V,A,' // uplo //
1798 $ ')', iinfo, n, jtype, ioldsd
1799 info = abs( iinfo )
1800 IF( iinfo.LT.0 ) THEN
1801 RETURN
1802 ELSE
1803 result( ntest ) = ulpinv
1804 result( ntest+1 ) = ulpinv
1805 result( ntest+2 ) = ulpinv
1806 GO TO 900
1807 END IF
1808 END IF
1809
1810
1811
1812 CALL dsyt21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
1813 $ ldu, tau, work, result( ntest ) )
1814
1815 ntest = ntest + 2
1816
1817 IF( iuplo.EQ.1 ) THEN
1818 indx = 1
1819 DO 860 j = 1, n
1820 DO 850 i = 1, j
1821 work( indx ) = a( i, j )
1822 indx = indx + 1
1823 850 CONTINUE
1824 860 CONTINUE
1825 ELSE
1826 indx = 1
1827 DO 880 j = 1, n
1828 DO 870 i = j, n
1829 work( indx ) = a( i, j )
1830 indx = indx + 1
1831 870 CONTINUE
1832 880 CONTINUE
1833 END IF
1834
1835 srnamt = 'DSPEVX'
1836 CALL dspevx(
'N',
'A', uplo, n, work, vl, vu, il, iu,
1837 $ abstol, m2, wa2, z, ldu, v, iwork,
1838 $ iwork( 5*n+1 ), iinfo )
1839 IF( iinfo.NE.0 ) THEN
1840 WRITE( nounit, fmt = 9999 )'DSPEVX(N,A,' // uplo //
1841 $ ')', iinfo, n, jtype, ioldsd
1842 info = abs( iinfo )
1843 IF( iinfo.LT.0 ) THEN
1844 RETURN
1845 ELSE
1846 result( ntest ) = ulpinv
1847 GO TO 900
1848 END IF
1849 END IF
1850
1851
1852
1853 temp1 = zero
1854 temp2 = zero
1855 DO 890 j = 1, n
1856 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
1857 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
1858 890 CONTINUE
1859 result( ntest ) = temp2 / max( unfl,
1860 $ ulp*max( temp1, temp2 ) )
1861
1862 900 CONTINUE
1863 IF( iuplo.EQ.1 ) THEN
1864 indx = 1
1865 DO 920 j = 1, n
1866 DO 910 i = 1, j
1867 work( indx ) = a( i, j )
1868 indx = indx + 1
1869 910 CONTINUE
1870 920 CONTINUE
1871 ELSE
1872 indx = 1
1873 DO 940 j = 1, n
1874 DO 930 i = j, n
1875 work( indx ) = a( i, j )
1876 indx = indx + 1
1877 930 CONTINUE
1878 940 CONTINUE
1879 END IF
1880
1881 ntest = ntest + 1
1882
1883 srnamt = 'DSPEVX'
1884 CALL dspevx(
'V',
'I', uplo, n, work, vl, vu, il, iu,
1885 $ abstol, m2, wa2, z, ldu, v, iwork,
1886 $ iwork( 5*n+1 ), iinfo )
1887 IF( iinfo.NE.0 ) THEN
1888 WRITE( nounit, fmt = 9999 )'DSPEVX(V,I,' // uplo //
1889 $ ')', iinfo, n, jtype, ioldsd
1890 info = abs( iinfo )
1891 IF( iinfo.LT.0 ) THEN
1892 RETURN
1893 ELSE
1894 result( ntest ) = ulpinv
1895 result( ntest+1 ) = ulpinv
1896 result( ntest+2 ) = ulpinv
1897 GO TO 990
1898 END IF
1899 END IF
1900
1901
1902
1903 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
1904 $ v, ldu, tau, work, result( ntest ) )
1905
1906 ntest = ntest + 2
1907
1908 IF( iuplo.EQ.1 ) THEN
1909 indx = 1
1910 DO 960 j = 1, n
1911 DO 950 i = 1, j
1912 work( indx ) = a( i, j )
1913 indx = indx + 1
1914 950 CONTINUE
1915 960 CONTINUE
1916 ELSE
1917 indx = 1
1918 DO 980 j = 1, n
1919 DO 970 i = j, n
1920 work( indx ) = a( i, j )
1921 indx = indx + 1
1922 970 CONTINUE
1923 980 CONTINUE
1924 END IF
1925
1926 srnamt = 'DSPEVX'
1927 CALL dspevx(
'N',
'I', uplo, n, work, vl, vu, il, iu,
1928 $ abstol, m3, wa3, z, ldu, v, iwork,
1929 $ iwork( 5*n+1 ), iinfo )
1930 IF( iinfo.NE.0 ) THEN
1931 WRITE( nounit, fmt = 9999 )'DSPEVX(N,I,' // uplo //
1932 $ ')', iinfo, n, jtype, ioldsd
1933 info = abs( iinfo )
1934 IF( iinfo.LT.0 ) THEN
1935 RETURN
1936 ELSE
1937 result( ntest ) = ulpinv
1938 GO TO 990
1939 END IF
1940 END IF
1941
1942 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
1943 result( ntest ) = ulpinv
1944 GO TO 990
1945 END IF
1946
1947
1948
1949 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
1950 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
1951 IF( n.GT.0 ) THEN
1952 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
1953 ELSE
1954 temp3 = zero
1955 END IF
1956 result( ntest ) = ( temp1+temp2 ) /
1957 $ max( unfl, temp3*ulp )
1958
1959 990 CONTINUE
1960 IF( iuplo.EQ.1 ) THEN
1961 indx = 1
1962 DO 1010 j = 1, n
1963 DO 1000 i = 1, j
1964 work( indx ) = a( i, j )
1965 indx = indx + 1
1966 1000 CONTINUE
1967 1010 CONTINUE
1968 ELSE
1969 indx = 1
1970 DO 1030 j = 1, n
1971 DO 1020 i = j, n
1972 work( indx ) = a( i, j )
1973 indx = indx + 1
1974 1020 CONTINUE
1975 1030 CONTINUE
1976 END IF
1977
1978 ntest = ntest + 1
1979
1980 srnamt = 'DSPEVX'
1981 CALL dspevx(
'V',
'V', uplo, n, work, vl, vu, il, iu,
1982 $ abstol, m2, wa2, z, ldu, v, iwork,
1983 $ iwork( 5*n+1 ), iinfo )
1984 IF( iinfo.NE.0 ) THEN
1985 WRITE( nounit, fmt = 9999 )'DSPEVX(V,V,' // uplo //
1986 $ ')', iinfo, n, jtype, ioldsd
1987 info = abs( iinfo )
1988 IF( iinfo.LT.0 ) THEN
1989 RETURN
1990 ELSE
1991 result( ntest ) = ulpinv
1992 result( ntest+1 ) = ulpinv
1993 result( ntest+2 ) = ulpinv
1994 GO TO 1080
1995 END IF
1996 END IF
1997
1998
1999
2000 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2001 $ v, ldu, tau, work, result( ntest ) )
2002
2003 ntest = ntest + 2
2004
2005 IF( iuplo.EQ.1 ) THEN
2006 indx = 1
2007 DO 1050 j = 1, n
2008 DO 1040 i = 1, j
2009 work( indx ) = a( i, j )
2010 indx = indx + 1
2011 1040 CONTINUE
2012 1050 CONTINUE
2013 ELSE
2014 indx = 1
2015 DO 1070 j = 1, n
2016 DO 1060 i = j, n
2017 work( indx ) = a( i, j )
2018 indx = indx + 1
2019 1060 CONTINUE
2020 1070 CONTINUE
2021 END IF
2022
2023 srnamt = 'DSPEVX'
2024 CALL dspevx(
'N',
'V', uplo, n, work, vl, vu, il, iu,
2025 $ abstol, m3, wa3, z, ldu, v, iwork,
2026 $ iwork( 5*n+1 ), iinfo )
2027 IF( iinfo.NE.0 ) THEN
2028 WRITE( nounit, fmt = 9999 )'DSPEVX(N,V,' // uplo //
2029 $ ')', iinfo, n, jtype, ioldsd
2030 info = abs( iinfo )
2031 IF( iinfo.LT.0 ) THEN
2032 RETURN
2033 ELSE
2034 result( ntest ) = ulpinv
2035 GO TO 1080
2036 END IF
2037 END IF
2038
2039 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2040 result( ntest ) = ulpinv
2041 GO TO 1080
2042 END IF
2043
2044
2045
2046 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2047 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2048 IF( n.GT.0 ) THEN
2049 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2050 ELSE
2051 temp3 = zero
2052 END IF
2053 result( ntest ) = ( temp1+temp2 ) /
2054 $ max( unfl, temp3*ulp )
2055
2056 1080 CONTINUE
2057
2058
2059
2060 IF( jtype.LE.7 ) THEN
2061 kd = 1
2062 ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
2063 kd = max( n-1, 0 )
2064 ELSE
2065 kd = ihbw
2066 END IF
2067
2068
2069
2070
2071 IF( iuplo.EQ.1 ) THEN
2072 DO 1100 j = 1, n
2073 DO 1090 i = max( 1, j-kd ), j
2074 v( kd+1+i-j, j ) = a( i, j )
2075 1090 CONTINUE
2076 1100 CONTINUE
2077 ELSE
2078 DO 1120 j = 1, n
2079 DO 1110 i = j, min( n, j+kd )
2080 v( 1+i-j, j ) = a( i, j )
2081 1110 CONTINUE
2082 1120 CONTINUE
2083 END IF
2084
2085 ntest = ntest + 1
2086 srnamt = 'DSBEV'
2087 CALL dsbev(
'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
2088 $ iinfo )
2089 IF( iinfo.NE.0 ) THEN
2090 WRITE( nounit, fmt = 9999 )'DSBEV(V,' // uplo // ')',
2091 $ iinfo, n, jtype, ioldsd
2092 info = abs( iinfo )
2093 IF( iinfo.LT.0 ) THEN
2094 RETURN
2095 ELSE
2096 result( ntest ) = ulpinv
2097 result( ntest+1 ) = ulpinv
2098 result( ntest+2 ) = ulpinv
2099 GO TO 1180
2100 END IF
2101 END IF
2102
2103
2104
2105 CALL dsyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
2106 $ ldu, tau, work, result( ntest ) )
2107
2108 IF( iuplo.EQ.1 ) THEN
2109 DO 1140 j = 1, n
2110 DO 1130 i = max( 1, j-kd ), j
2111 v( kd+1+i-j, j ) = a( i, j )
2112 1130 CONTINUE
2113 1140 CONTINUE
2114 ELSE
2115 DO 1160 j = 1, n
2116 DO 1150 i = j, min( n, j+kd )
2117 v( 1+i-j, j ) = a( i, j )
2118 1150 CONTINUE
2119 1160 CONTINUE
2120 END IF
2121
2122 ntest = ntest + 2
2123 srnamt = 'DSBEV_2STAGE'
2124 CALL dsbev_2stage(
'N', uplo, n, kd, v, ldu, d3, z, ldu,
2125 $ work, lwork, iinfo )
2126 IF( iinfo.NE.0 ) THEN
2127 WRITE( nounit, fmt = 9999 )
2128 $ 'DSBEV_2STAGE(N,' // uplo // ')',
2129 $ iinfo, n, jtype, ioldsd
2130 info = abs( iinfo )
2131 IF( iinfo.LT.0 ) THEN
2132 RETURN
2133 ELSE
2134 result( ntest ) = ulpinv
2135 GO TO 1180
2136 END IF
2137 END IF
2138
2139
2140
2141 temp1 = zero
2142 temp2 = zero
2143 DO 1170 j = 1, n
2144 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2145 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2146 1170 CONTINUE
2147 result( ntest ) = temp2 / max( unfl,
2148 $ ulp*max( temp1, temp2 ) )
2149
2150
2151
2152
2153 1180 CONTINUE
2154 IF( iuplo.EQ.1 ) THEN
2155 DO 1200 j = 1, n
2156 DO 1190 i = max( 1, j-kd ), j
2157 v( kd+1+i-j, j ) = a( i, j )
2158 1190 CONTINUE
2159 1200 CONTINUE
2160 ELSE
2161 DO 1220 j = 1, n
2162 DO 1210 i = j, min( n, j+kd )
2163 v( 1+i-j, j ) = a( i, j )
2164 1210 CONTINUE
2165 1220 CONTINUE
2166 END IF
2167
2168 ntest = ntest + 1
2169 srnamt = 'DSBEVX'
2170 CALL dsbevx(
'V',
'A', uplo, n, kd, v, ldu, u, ldu, vl,
2171 $ vu, il, iu, abstol, m, wa2, z, ldu, work,
2172 $ iwork, iwork( 5*n+1 ), iinfo )
2173 IF( iinfo.NE.0 ) THEN
2174 WRITE( nounit, fmt = 9999 )'DSBEVX(V,A,' // uplo //
2175 $ ')', iinfo, n, jtype, ioldsd
2176 info = abs( iinfo )
2177 IF( iinfo.LT.0 ) THEN
2178 RETURN
2179 ELSE
2180 result( ntest ) = ulpinv
2181 result( ntest+1 ) = ulpinv
2182 result( ntest+2 ) = ulpinv
2183 GO TO 1280
2184 END IF
2185 END IF
2186
2187
2188
2189 CALL dsyt21( 1, uplo, n, 0, a, ldu, wa2, d2, z, ldu, v,
2190 $ ldu, tau, work, result( ntest ) )
2191
2192 ntest = ntest + 2
2193
2194 IF( iuplo.EQ.1 ) THEN
2195 DO 1240 j = 1, n
2196 DO 1230 i = max( 1, j-kd ), j
2197 v( kd+1+i-j, j ) = a( i, j )
2198 1230 CONTINUE
2199 1240 CONTINUE
2200 ELSE
2201 DO 1260 j = 1, n
2202 DO 1250 i = j, min( n, j+kd )
2203 v( 1+i-j, j ) = a( i, j )
2204 1250 CONTINUE
2205 1260 CONTINUE
2206 END IF
2207
2208 srnamt = 'DSBEVX_2STAGE'
2210 $ u, ldu, vl, vu, il, iu, abstol, m3, wa3,
2211 $ z, ldu, work, lwork, iwork, iwork( 5*n+1 ),
2212 $ iinfo )
2213 IF( iinfo.NE.0 ) THEN
2214 WRITE( nounit, fmt = 9999 )
2215 $ 'DSBEVX_2STAGE(N,A,' // uplo //
2216 $ ')', iinfo, n, jtype, ioldsd
2217 info = abs( iinfo )
2218 IF( iinfo.LT.0 ) THEN
2219 RETURN
2220 ELSE
2221 result( ntest ) = ulpinv
2222 GO TO 1280
2223 END IF
2224 END IF
2225
2226
2227
2228 temp1 = zero
2229 temp2 = zero
2230 DO 1270 j = 1, n
2231 temp1 = max( temp1, abs( wa2( j ) ), abs( wa3( j ) ) )
2232 temp2 = max( temp2, abs( wa2( j )-wa3( j ) ) )
2233 1270 CONTINUE
2234 result( ntest ) = temp2 / max( unfl,
2235 $ ulp*max( temp1, temp2 ) )
2236
2237 1280 CONTINUE
2238 ntest = ntest + 1
2239 IF( iuplo.EQ.1 ) THEN
2240 DO 1300 j = 1, n
2241 DO 1290 i = max( 1, j-kd ), j
2242 v( kd+1+i-j, j ) = a( i, j )
2243 1290 CONTINUE
2244 1300 CONTINUE
2245 ELSE
2246 DO 1320 j = 1, n
2247 DO 1310 i = j, min( n, j+kd )
2248 v( 1+i-j, j ) = a( i, j )
2249 1310 CONTINUE
2250 1320 CONTINUE
2251 END IF
2252
2253 srnamt = 'DSBEVX'
2254 CALL dsbevx(
'V',
'I', uplo, n, kd, v, ldu, u, ldu, vl,
2255 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
2256 $ iwork, iwork( 5*n+1 ), iinfo )
2257 IF( iinfo.NE.0 ) THEN
2258 WRITE( nounit, fmt = 9999 )'DSBEVX(V,I,' // uplo //
2259 $ ')', iinfo, n, jtype, ioldsd
2260 info = abs( iinfo )
2261 IF( iinfo.LT.0 ) THEN
2262 RETURN
2263 ELSE
2264 result( ntest ) = ulpinv
2265 result( ntest+1 ) = ulpinv
2266 result( ntest+2 ) = ulpinv
2267 GO TO 1370
2268 END IF
2269 END IF
2270
2271
2272
2273 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2274 $ v, ldu, tau, work, result( ntest ) )
2275
2276 ntest = ntest + 2
2277
2278 IF( iuplo.EQ.1 ) THEN
2279 DO 1340 j = 1, n
2280 DO 1330 i = max( 1, j-kd ), j
2281 v( kd+1+i-j, j ) = a( i, j )
2282 1330 CONTINUE
2283 1340 CONTINUE
2284 ELSE
2285 DO 1360 j = 1, n
2286 DO 1350 i = j, min( n, j+kd )
2287 v( 1+i-j, j ) = a( i, j )
2288 1350 CONTINUE
2289 1360 CONTINUE
2290 END IF
2291
2292 srnamt = 'DSBEVX_2STAGE'
2294 $ u, ldu, vl, vu, il, iu, abstol, m3, wa3,
2295 $ z, ldu, work, lwork, iwork, iwork( 5*n+1 ),
2296 $ iinfo )
2297 IF( iinfo.NE.0 ) THEN
2298 WRITE( nounit, fmt = 9999 )
2299 $ 'DSBEVX_2STAGE(N,I,' // uplo //
2300 $ ')', iinfo, n, jtype, ioldsd
2301 info = abs( iinfo )
2302 IF( iinfo.LT.0 ) THEN
2303 RETURN
2304 ELSE
2305 result( ntest ) = ulpinv
2306 GO TO 1370
2307 END IF
2308 END IF
2309
2310
2311
2312 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2313 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2314 IF( n.GT.0 ) THEN
2315 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2316 ELSE
2317 temp3 = zero
2318 END IF
2319 result( ntest ) = ( temp1+temp2 ) /
2320 $ max( unfl, temp3*ulp )
2321
2322 1370 CONTINUE
2323 ntest = ntest + 1
2324 IF( iuplo.EQ.1 ) THEN
2325 DO 1390 j = 1, n
2326 DO 1380 i = max( 1, j-kd ), j
2327 v( kd+1+i-j, j ) = a( i, j )
2328 1380 CONTINUE
2329 1390 CONTINUE
2330 ELSE
2331 DO 1410 j = 1, n
2332 DO 1400 i = j, min( n, j+kd )
2333 v( 1+i-j, j ) = a( i, j )
2334 1400 CONTINUE
2335 1410 CONTINUE
2336 END IF
2337
2338 srnamt = 'DSBEVX'
2339 CALL dsbevx(
'V',
'V', uplo, n, kd, v, ldu, u, ldu, vl,
2340 $ vu, il, iu, abstol, m2, wa2, z, ldu, work,
2341 $ iwork, iwork( 5*n+1 ), iinfo )
2342 IF( iinfo.NE.0 ) THEN
2343 WRITE( nounit, fmt = 9999 )'DSBEVX(V,V,' // uplo //
2344 $ ')', iinfo, n, jtype, ioldsd
2345 info = abs( iinfo )
2346 IF( iinfo.LT.0 ) THEN
2347 RETURN
2348 ELSE
2349 result( ntest ) = ulpinv
2350 result( ntest+1 ) = ulpinv
2351 result( ntest+2 ) = ulpinv
2352 GO TO 1460
2353 END IF
2354 END IF
2355
2356
2357
2358 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2359 $ v, ldu, tau, work, result( ntest ) )
2360
2361 ntest = ntest + 2
2362
2363 IF( iuplo.EQ.1 ) THEN
2364 DO 1430 j = 1, n
2365 DO 1420 i = max( 1, j-kd ), j
2366 v( kd+1+i-j, j ) = a( i, j )
2367 1420 CONTINUE
2368 1430 CONTINUE
2369 ELSE
2370 DO 1450 j = 1, n
2371 DO 1440 i = j, min( n, j+kd )
2372 v( 1+i-j, j ) = a( i, j )
2373 1440 CONTINUE
2374 1450 CONTINUE
2375 END IF
2376
2377 srnamt = 'DSBEVX_2STAGE'
2379 $ u, ldu, vl, vu, il, iu, abstol, m3, wa3,
2380 $ z, ldu, work, lwork, iwork, iwork( 5*n+1 ),
2381 $ iinfo )
2382 IF( iinfo.NE.0 ) THEN
2383 WRITE( nounit, fmt = 9999 )
2384 $ 'DSBEVX_2STAGE(N,V,' // uplo //
2385 $ ')', iinfo, n, jtype, ioldsd
2386 info = abs( iinfo )
2387 IF( iinfo.LT.0 ) THEN
2388 RETURN
2389 ELSE
2390 result( ntest ) = ulpinv
2391 GO TO 1460
2392 END IF
2393 END IF
2394
2395 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2396 result( ntest ) = ulpinv
2397 GO TO 1460
2398 END IF
2399
2400
2401
2402 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2403 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2404 IF( n.GT.0 ) THEN
2405 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2406 ELSE
2407 temp3 = zero
2408 END IF
2409 result( ntest ) = ( temp1+temp2 ) /
2410 $ max( unfl, temp3*ulp )
2411
2412 1460 CONTINUE
2413
2414
2415
2416 CALL dlacpy(
' ', n, n, a, lda, v, ldu )
2417
2418 ntest = ntest + 1
2419 srnamt = 'DSYEVD'
2420 CALL dsyevd(
'V', uplo, n, a, ldu, d1, work, lwedc,
2421 $ iwork, liwedc, iinfo )
2422 IF( iinfo.NE.0 ) THEN
2423 WRITE( nounit, fmt = 9999 )'DSYEVD(V,' // uplo //
2424 $ ')', iinfo, n, jtype, ioldsd
2425 info = abs( iinfo )
2426 IF( iinfo.LT.0 ) THEN
2427 RETURN
2428 ELSE
2429 result( ntest ) = ulpinv
2430 result( ntest+1 ) = ulpinv
2431 result( ntest+2 ) = ulpinv
2432 GO TO 1480
2433 END IF
2434 END IF
2435
2436
2437
2438 CALL dsyt21( 1, uplo, n, 0, v, ldu, d1, d2, a, ldu, z,
2439 $ ldu, tau, work, result( ntest ) )
2440
2441 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2442
2443 ntest = ntest + 2
2444 srnamt = 'DSYEVD_2STAGE'
2446 $ lwork, iwork, liwedc, iinfo )
2447 IF( iinfo.NE.0 ) THEN
2448 WRITE( nounit, fmt = 9999 )
2449 $ 'DSYEVD_2STAGE(N,' // uplo //
2450 $ ')', iinfo, n, jtype, ioldsd
2451 info = abs( iinfo )
2452 IF( iinfo.LT.0 ) THEN
2453 RETURN
2454 ELSE
2455 result( ntest ) = ulpinv
2456 GO TO 1480
2457 END IF
2458 END IF
2459
2460
2461
2462 temp1 = zero
2463 temp2 = zero
2464 DO 1470 j = 1, n
2465 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2466 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2467 1470 CONTINUE
2468 result( ntest ) = temp2 / max( unfl,
2469 $ ulp*max( temp1, temp2 ) )
2470
2471 1480 CONTINUE
2472
2473
2474
2475 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2476
2477
2478
2479
2480 IF( iuplo.EQ.1 ) THEN
2481 indx = 1
2482 DO 1500 j = 1, n
2483 DO 1490 i = 1, j
2484 work( indx ) = a( i, j )
2485 indx = indx + 1
2486 1490 CONTINUE
2487 1500 CONTINUE
2488 ELSE
2489 indx = 1
2490 DO 1520 j = 1, n
2491 DO 1510 i = j, n
2492 work( indx ) = a( i, j )
2493 indx = indx + 1
2494 1510 CONTINUE
2495 1520 CONTINUE
2496 END IF
2497
2498 ntest = ntest + 1
2499 srnamt = 'DSPEVD'
2500 CALL dspevd(
'V', uplo, n, work, d1, z, ldu,
2501 $ work( indx ), lwedc-indx+1, iwork, liwedc,
2502 $ iinfo )
2503 IF( iinfo.NE.0 ) THEN
2504 WRITE( nounit, fmt = 9999 )'DSPEVD(V,' // uplo //
2505 $ ')', iinfo, n, jtype, ioldsd
2506 info = abs( iinfo )
2507 IF( iinfo.LT.0 ) THEN
2508 RETURN
2509 ELSE
2510 result( ntest ) = ulpinv
2511 result( ntest+1 ) = ulpinv
2512 result( ntest+2 ) = ulpinv
2513 GO TO 1580
2514 END IF
2515 END IF
2516
2517
2518
2519 CALL dsyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
2520 $ ldu, tau, work, result( ntest ) )
2521
2522 IF( iuplo.EQ.1 ) THEN
2523 indx = 1
2524 DO 1540 j = 1, n
2525 DO 1530 i = 1, j
2526
2527 work( indx ) = a( i, j )
2528 indx = indx + 1
2529 1530 CONTINUE
2530 1540 CONTINUE
2531 ELSE
2532 indx = 1
2533 DO 1560 j = 1, n
2534 DO 1550 i = j, n
2535 work( indx ) = a( i, j )
2536 indx = indx + 1
2537 1550 CONTINUE
2538 1560 CONTINUE
2539 END IF
2540
2541 ntest = ntest + 2
2542 srnamt = 'DSPEVD'
2543 CALL dspevd(
'N', uplo, n, work, d3, z, ldu,
2544 $ work( indx ), lwedc-indx+1, iwork, liwedc,
2545 $ iinfo )
2546 IF( iinfo.NE.0 ) THEN
2547 WRITE( nounit, fmt = 9999 )'DSPEVD(N,' // uplo //
2548 $ ')', iinfo, n, jtype, ioldsd
2549 info = abs( iinfo )
2550 IF( iinfo.LT.0 ) THEN
2551 RETURN
2552 ELSE
2553 result( ntest ) = ulpinv
2554 GO TO 1580
2555 END IF
2556 END IF
2557
2558
2559
2560 temp1 = zero
2561 temp2 = zero
2562 DO 1570 j = 1, n
2563 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2564 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2565 1570 CONTINUE
2566 result( ntest ) = temp2 / max( unfl,
2567 $ ulp*max( temp1, temp2 ) )
2568 1580 CONTINUE
2569
2570
2571
2572 IF( jtype.LE.7 ) THEN
2573 kd = 1
2574 ELSE IF( jtype.GE.8 .AND. jtype.LE.15 ) THEN
2575 kd = max( n-1, 0 )
2576 ELSE
2577 kd = ihbw
2578 END IF
2579
2580
2581
2582
2583 IF( iuplo.EQ.1 ) THEN
2584 DO 1600 j = 1, n
2585 DO 1590 i = max( 1, j-kd ), j
2586 v( kd+1+i-j, j ) = a( i, j )
2587 1590 CONTINUE
2588 1600 CONTINUE
2589 ELSE
2590 DO 1620 j = 1, n
2591 DO 1610 i = j, min( n, j+kd )
2592 v( 1+i-j, j ) = a( i, j )
2593 1610 CONTINUE
2594 1620 CONTINUE
2595 END IF
2596
2597 ntest = ntest + 1
2598 srnamt = 'DSBEVD'
2599 CALL dsbevd(
'V', uplo, n, kd, v, ldu, d1, z, ldu, work,
2600 $ lwedc, iwork, liwedc, iinfo )
2601 IF( iinfo.NE.0 ) THEN
2602 WRITE( nounit, fmt = 9999 )'DSBEVD(V,' // uplo //
2603 $ ')', iinfo, n, jtype, ioldsd
2604 info = abs( iinfo )
2605 IF( iinfo.LT.0 ) THEN
2606 RETURN
2607 ELSE
2608 result( ntest ) = ulpinv
2609 result( ntest+1 ) = ulpinv
2610 result( ntest+2 ) = ulpinv
2611 GO TO 1680
2612 END IF
2613 END IF
2614
2615
2616
2617 CALL dsyt21( 1, uplo, n, 0, a, lda, d1, d2, z, ldu, v,
2618 $ ldu, tau, work, result( ntest ) )
2619
2620 IF( iuplo.EQ.1 ) THEN
2621 DO 1640 j = 1, n
2622 DO 1630 i = max( 1, j-kd ), j
2623 v( kd+1+i-j, j ) = a( i, j )
2624 1630 CONTINUE
2625 1640 CONTINUE
2626 ELSE
2627 DO 1660 j = 1, n
2628 DO 1650 i = j, min( n, j+kd )
2629 v( 1+i-j, j ) = a( i, j )
2630 1650 CONTINUE
2631 1660 CONTINUE
2632 END IF
2633
2634 ntest = ntest + 2
2635 srnamt = 'DSBEVD_2STAGE'
2637 $ work, lwork, iwork, liwedc, iinfo )
2638 IF( iinfo.NE.0 ) THEN
2639 WRITE( nounit, fmt = 9999 )
2640 $ 'DSBEVD_2STAGE(N,' // uplo //
2641 $ ')', iinfo, n, jtype, ioldsd
2642 info = abs( iinfo )
2643 IF( iinfo.LT.0 ) THEN
2644 RETURN
2645 ELSE
2646 result( ntest ) = ulpinv
2647 GO TO 1680
2648 END IF
2649 END IF
2650
2651
2652
2653 temp1 = zero
2654 temp2 = zero
2655 DO 1670 j = 1, n
2656 temp1 = max( temp1, abs( d1( j ) ), abs( d3( j ) ) )
2657 temp2 = max( temp2, abs( d1( j )-d3( j ) ) )
2658 1670 CONTINUE
2659 result( ntest ) = temp2 / max( unfl,
2660 $ ulp*max( temp1, temp2 ) )
2661
2662 1680 CONTINUE
2663
2664
2665 CALL dlacpy(
' ', n, n, a, lda, v, ldu )
2666 ntest = ntest + 1
2667 srnamt = 'DSYEVR'
2668 CALL dsyevr(
'V',
'A', uplo, n, a, ldu, vl, vu, il, iu,
2669 $ abstol, m, wa1, z, ldu, iwork, work, lwork,
2670 $ iwork(2*n+1), liwork-2*n, iinfo )
2671 IF( iinfo.NE.0 ) THEN
2672 WRITE( nounit, fmt = 9999 )'DSYEVR(V,A,' // uplo //
2673 $ ')', iinfo, n, jtype, ioldsd
2674 info = abs( iinfo )
2675 IF( iinfo.LT.0 ) THEN
2676 RETURN
2677 ELSE
2678 result( ntest ) = ulpinv
2679 result( ntest+1 ) = ulpinv
2680 result( ntest+2 ) = ulpinv
2681 GO TO 1700
2682 END IF
2683 END IF
2684
2685
2686
2687 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2688
2689 CALL dsyt21( 1, uplo, n, 0, a, ldu, wa1, d2, z, ldu, v,
2690 $ ldu, tau, work, result( ntest ) )
2691
2692 ntest = ntest + 2
2693 srnamt = 'DSYEVR_2STAGE'
2695 $ il, iu, abstol, m2, wa2, z, ldu, iwork,
2696 $ work, lwork, iwork(2*n+1), liwork-2*n,
2697 $ iinfo )
2698 IF( iinfo.NE.0 ) THEN
2699 WRITE( nounit, fmt = 9999 )
2700 $ 'DSYEVR_2STAGE(N,A,' // uplo //
2701 $ ')', iinfo, n, jtype, ioldsd
2702 info = abs( iinfo )
2703 IF( iinfo.LT.0 ) THEN
2704 RETURN
2705 ELSE
2706 result( ntest ) = ulpinv
2707 GO TO 1700
2708 END IF
2709 END IF
2710
2711
2712
2713 temp1 = zero
2714 temp2 = zero
2715 DO 1690 j = 1, n
2716 temp1 = max( temp1, abs( wa1( j ) ), abs( wa2( j ) ) )
2717 temp2 = max( temp2, abs( wa1( j )-wa2( j ) ) )
2718 1690 CONTINUE
2719 result( ntest ) = temp2 / max( unfl,
2720 $ ulp*max( temp1, temp2 ) )
2721
2722 1700 CONTINUE
2723
2724 ntest = ntest + 1
2725 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2726 srnamt = 'DSYEVR'
2727 CALL dsyevr(
'V',
'I', uplo, n, a, ldu, vl, vu, il, iu,
2728 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
2729 $ iwork(2*n+1), liwork-2*n, iinfo )
2730 IF( iinfo.NE.0 ) THEN
2731 WRITE( nounit, fmt = 9999 )'DSYEVR(V,I,' // uplo //
2732 $ ')', iinfo, n, jtype, ioldsd
2733 info = abs( iinfo )
2734 IF( iinfo.LT.0 ) THEN
2735 RETURN
2736 ELSE
2737 result( ntest ) = ulpinv
2738 result( ntest+1 ) = ulpinv
2739 result( ntest+2 ) = ulpinv
2740 GO TO 1710
2741 END IF
2742 END IF
2743
2744
2745
2746 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2747
2748 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2749 $ v, ldu, tau, work, result( ntest ) )
2750
2751 ntest = ntest + 2
2752 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2753 srnamt = 'DSYEVR_2STAGE'
2755 $ il, iu, abstol, m3, wa3, z, ldu, iwork,
2756 $ work, lwork, iwork(2*n+1), liwork-2*n,
2757 $ iinfo )
2758 IF( iinfo.NE.0 ) THEN
2759 WRITE( nounit, fmt = 9999 )
2760 $ 'DSYEVR_2STAGE(N,I,' // uplo //
2761 $ ')', iinfo, n, jtype, ioldsd
2762 info = abs( iinfo )
2763 IF( iinfo.LT.0 ) THEN
2764 RETURN
2765 ELSE
2766 result( ntest ) = ulpinv
2767 GO TO 1710
2768 END IF
2769 END IF
2770
2771
2772
2773 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2774 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2775 result( ntest ) = ( temp1+temp2 ) /
2776 $ max( unfl, ulp*temp3 )
2777 1710 CONTINUE
2778
2779 ntest = ntest + 1
2780 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2781 srnamt = 'DSYEVR'
2782 CALL dsyevr(
'V',
'V', uplo, n, a, ldu, vl, vu, il, iu,
2783 $ abstol, m2, wa2, z, ldu, iwork, work, lwork,
2784 $ iwork(2*n+1), liwork-2*n, iinfo )
2785 IF( iinfo.NE.0 ) THEN
2786 WRITE( nounit, fmt = 9999 )'DSYEVR(V,V,' // uplo //
2787 $ ')', iinfo, n, jtype, ioldsd
2788 info = abs( iinfo )
2789 IF( iinfo.LT.0 ) THEN
2790 RETURN
2791 ELSE
2792 result( ntest ) = ulpinv
2793 result( ntest+1 ) = ulpinv
2794 result( ntest+2 ) = ulpinv
2795 GO TO 1720
2796 END IF
2797 END IF
2798
2799
2800
2801 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2802
2803 CALL dsyt22( 1, uplo, n, m2, 0, a, ldu, wa2, d2, z, ldu,
2804 $ v, ldu, tau, work, result( ntest ) )
2805
2806 ntest = ntest + 2
2807 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2808 srnamt = 'DSYEVR_2STAGE'
2810 $ il, iu, abstol, m3, wa3, z, ldu, iwork,
2811 $ work, lwork, iwork(2*n+1), liwork-2*n,
2812 $ iinfo )
2813 IF( iinfo.NE.0 ) THEN
2814 WRITE( nounit, fmt = 9999 )
2815 $ 'DSYEVR_2STAGE(N,V,' // uplo //
2816 $ ')', iinfo, n, jtype, ioldsd
2817 info = abs( iinfo )
2818 IF( iinfo.LT.0 ) THEN
2819 RETURN
2820 ELSE
2821 result( ntest ) = ulpinv
2822 GO TO 1720
2823 END IF
2824 END IF
2825
2826 IF( m3.EQ.0 .AND. n.GT.0 ) THEN
2827 result( ntest ) = ulpinv
2828 GO TO 1720
2829 END IF
2830
2831
2832
2833 temp1 =
dsxt1( 1, wa2, m2, wa3, m3, abstol, ulp, unfl )
2834 temp2 =
dsxt1( 1, wa3, m3, wa2, m2, abstol, ulp, unfl )
2835 IF( n.GT.0 ) THEN
2836 temp3 = max( abs( wa1( 1 ) ), abs( wa1( n ) ) )
2837 ELSE
2838 temp3 = zero
2839 END IF
2840 result( ntest ) = ( temp1+temp2 ) /
2841 $ max( unfl, temp3*ulp )
2842
2843 CALL dlacpy(
' ', n, n, v, ldu, a, lda )
2844
2845 1720 CONTINUE
2846
2847
2848
2849 ntestt = ntestt + ntest
2850
2851 CALL dlafts(
'DST', n, n, jtype, ntest, result, ioldsd,
2852 $ thresh, nounit, nerrs )
2853
2854 1730 CONTINUE
2855 1740 CONTINUE
2856
2857
2858
2859 CALL alasvm(
'DST', nounit, nerrs, ntestt, 0 )
2860
2861 9999 FORMAT( ' DDRVST2STG: ', a, ' returned INFO=', i6, '.', / 9x,
2862 $ 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
2863
2864 RETURN
2865
2866
2867
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
subroutine xerbla(srname, info)
subroutine dlafts(type, m, n, imat, ntests, result, iseed, thresh, iounit, ie)
DLAFTS
double precision function dlarnd(idist, iseed)
DLARND
subroutine dlatmr(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)
DLATMR
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
subroutine dstt21(n, kband, ad, ae, sd, se, u, ldu, work, result)
DSTT21
subroutine dstt22(n, m, kband, ad, ae, sd, se, u, ldu, work, ldwork, result)
DSTT22
double precision function dsxt1(ijob, d1, n1, d2, n2, abstol, ulp, unfl)
DSXT1
subroutine dsyt21(itype, uplo, n, kband, a, lda, d, e, u, ldu, v, ldv, tau, work, result)
DSYT21
subroutine dsyt22(itype, uplo, n, m, kband, a, lda, d, e, u, ldu, v, ldv, tau, work, result)
DSYT22
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 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_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 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 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 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 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 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 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...
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 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 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 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 dsytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
DSYTRD_2STAGE
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
subroutine dsytrd_sy2sb(uplo, n, kd, a, lda, ab, ldab, tau, work, lwork, info)
DSYTRD_SY2SB
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 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 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 dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
double precision function dlamch(cmach)
DLAMCH
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
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 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 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...