LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ schkbd()

subroutine schkbd ( integer  NSIZES,
integer, dimension( * )  MVAL,
integer, dimension( * )  NVAL,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer  NRHS,
integer, dimension( 4 )  ISEED,
real  THRESH,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  BD,
real, dimension( * )  BE,
real, dimension( * )  S1,
real, dimension( * )  S2,
real, dimension( ldx, * )  X,
integer  LDX,
real, dimension( ldx, * )  Y,
real, dimension( ldx, * )  Z,
real, dimension( ldq, * )  Q,
integer  LDQ,
real, dimension( ldpt, * )  PT,
integer  LDPT,
real, dimension( ldpt, * )  U,
real, dimension( ldpt, * )  VT,
real, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  NOUT,
integer  INFO 
)

SCHKBD

Purpose:
 SCHKBD checks the singular value decomposition (SVD) routines.

 SGEBRD reduces a real general m by n matrix A to upper or lower
 bidiagonal form B by an orthogonal transformation:  Q' * A * P = B
 (or A = Q * B * P').  The matrix B is upper bidiagonal if m >= n
 and lower bidiagonal if m < n.

 SORGBR generates the orthogonal matrices Q and P' from SGEBRD.
 Note that Q and P are not necessarily square.

 SBDSQR computes the singular value decomposition of the bidiagonal
 matrix B as B = U S V'.  It is called three times to compute
    1)  B = U S1 V', where S1 is the diagonal matrix of singular
        values and the columns of the matrices U and V are the left
        and right singular vectors, respectively, of B.
    2)  Same as 1), but the singular values are stored in S2 and the
        singular vectors are not computed.
    3)  A = (UQ) S (P'V'), the SVD of the original matrix A.
 In addition, SBDSQR has an option to apply the left orthogonal matrix
 U to a matrix X, useful in least squares applications.

 SBDSDC computes the singular value decomposition of the bidiagonal
 matrix B as B = U S V' using divide-and-conquer. It is called twice
 to compute
    1) B = U S1 V', where S1 is the diagonal matrix of singular
        values and the columns of the matrices U and V are the left
        and right singular vectors, respectively, of B.
    2) Same as 1), but the singular values are stored in S2 and the
        singular vectors are not computed.

  SBDSVDX computes the singular value decomposition of the bidiagonal
  matrix B as B = U S V' using bisection and inverse iteration. It is
  called six times to compute
     1) B = U S1 V', RANGE='A', where S1 is the diagonal matrix of singular
         values and the columns of the matrices U and V are the left
         and right singular vectors, respectively, of B.
     2) Same as 1), but the singular values are stored in S2 and the
         singular vectors are not computed.
     3) B = U S1 V', RANGE='I', with where S1 is the diagonal matrix of singular
         values and the columns of the matrices U and V are the left
         and right singular vectors, respectively, of B
     4) Same as 3), but the singular values are stored in S2 and the
         singular vectors are not computed.
     5) B = U S1 V', RANGE='V', with where S1 is the diagonal matrix of singular
         values and the columns of the matrices U and V are the left
         and right singular vectors, respectively, of B
     6) Same as 5), but the singular values are stored in S2 and the
         singular vectors are not computed.

 For each pair of matrix dimensions (M,N) and each selected matrix
 type, an M by N matrix A and an M by NRHS matrix X are generated.
 The problem dimensions are as follows
    A:          M x N
    Q:          M x min(M,N) (but M x M if NRHS > 0)
    P:          min(M,N) x N
    B:          min(M,N) x min(M,N)
    U, V:       min(M,N) x min(M,N)
    S1, S2      diagonal, order min(M,N)
    X:          M x NRHS

 For each generated matrix, 14 tests are performed:

 Test SGEBRD and SORGBR

 (1)   | A - Q B PT | / ( |A| max(M,N) ulp ), PT = P'

 (2)   | I - Q' Q | / ( M ulp )

 (3)   | I - PT PT' | / ( N ulp )

 Test SBDSQR on bidiagonal matrix B

 (4)   | B - U S1 VT | / ( |B| min(M,N) ulp ), VT = V'

 (5)   | Y - U Z | / ( |Y| max(min(M,N),k) ulp ), where Y = Q' X
                                                  and   Z = U' Y.
 (6)   | I - U' U | / ( min(M,N) ulp )

 (7)   | I - VT VT' | / ( min(M,N) ulp )

 (8)   S1 contains min(M,N) nonnegative values in decreasing order.
       (Return 0 if true, 1/ULP if false.)

 (9)   | S1 - S2 | / ( |S1| ulp ), where S2 is computed without
                                   computing U and V.

 (10)  0 if the true singular values of B are within THRESH of
       those in S1.  2*THRESH if they are not.  (Tested using
       SSVDCH)

 Test SBDSQR on matrix A

 (11)  | A - (QU) S (VT PT) | / ( |A| max(M,N) ulp )

 (12)  | X - (QU) Z | / ( |X| max(M,k) ulp )

 (13)  | I - (QU)'(QU) | / ( M ulp )

 (14)  | I - (VT PT) (PT'VT') | / ( N ulp )

 Test SBDSDC on bidiagonal matrix B

 (15)  | B - U S1 VT | / ( |B| min(M,N) ulp ), VT = V'

 (16)  | I - U' U | / ( min(M,N) ulp )

 (17)  | I - VT VT' | / ( min(M,N) ulp )

 (18)  S1 contains min(M,N) nonnegative values in decreasing order.
       (Return 0 if true, 1/ULP if false.)

 (19)  | S1 - S2 | / ( |S1| ulp ), where S2 is computed without
                                   computing U and V.
  Test SBDSVDX on bidiagonal matrix B

  (20)  | B - U S1 VT | / ( |B| min(M,N) ulp ), VT = V'

  (21)  | I - U' U | / ( min(M,N) ulp )

  (22)  | I - VT VT' | / ( min(M,N) ulp )

  (23)  S1 contains min(M,N) nonnegative values in decreasing order.
        (Return 0 if true, 1/ULP if false.)

  (24)  | S1 - S2 | / ( |S1| ulp ), where S2 is computed without
                                    computing U and V.

  (25)  | S1 - U' B VT' | / ( |S| n ulp )    SBDSVDX('V', 'I')

  (26)  | I - U' U | / ( min(M,N) ulp )

  (27)  | I - VT VT' | / ( min(M,N) ulp )

  (28)  S1 contains min(M,N) nonnegative values in decreasing order.
        (Return 0 if true, 1/ULP if false.)

  (29)  | S1 - S2 | / ( |S1| ulp ), where S2 is computed without
                                    computing U and V.

  (30)  | S1 - U' B VT' | / ( |S1| n ulp )   SBDSVDX('V', 'V')

  (31)  | I - U' U | / ( min(M,N) ulp )

  (32)  | I - VT VT' | / ( min(M,N) ulp )

  (33)  S1 contains min(M,N) nonnegative values in decreasing order.
        (Return 0 if true, 1/ULP if false.)

  (34)  | S1 - S2 | / ( |S1| ulp ), where S2 is computed without
                                    computing U and V.

 The possible matrix types are

 (1)  The zero matrix.
 (2)  The identity matrix.

 (3)  A diagonal matrix with evenly spaced entries
      1, ..., ULP  and random signs.
      (ULP = (first number larger than 1) - 1 )
 (4)  A diagonal matrix with geometrically spaced entries
      1, ..., ULP  and random signs.
 (5)  A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
      and random signs.

 (6)  Same as (3), but multiplied by SQRT( overflow threshold )
 (7)  Same as (3), but multiplied by SQRT( underflow threshold )

 (8)  A matrix of the form  U D V, where U and V are orthogonal and
      D has evenly spaced entries 1, ..., ULP with random signs
      on the diagonal.

 (9)  A matrix of the form  U D V, where U and V are orthogonal and
      D has geometrically spaced entries 1, ..., ULP with random
      signs on the diagonal.

 (10) A matrix of the form  U D V, where U and V are orthogonal and
      D has "clustered" entries 1, ULP,..., ULP with random
      signs on the diagonal.

 (11) Same as (8), but multiplied by SQRT( overflow threshold )
 (12) Same as (8), but multiplied by SQRT( underflow threshold )

 (13) Rectangular matrix with random entries chosen from (-1,1).
 (14) Same as (13), but multiplied by SQRT( overflow threshold )
 (15) Same as (13), but multiplied by SQRT( underflow threshold )

 Special case:
 (16) A bidiagonal matrix with random entries chosen from a
      logarithmic distribution on [ulp^2,ulp^(-2)]  (I.e., each
      entry is  e^x, where x is chosen uniformly on
      [ 2 log(ulp), -2 log(ulp) ] .)  For *this* type:
      (a) SGEBRD is not called to reduce it to bidiagonal form.
      (b) the bidiagonal is  min(M,N) x min(M,N); if M<N, the
          matrix will be lower bidiagonal, otherwise upper.
      (c) only tests 5--8 and 14 are performed.

 A subset of the full set of matrix types may be selected through
 the logical array DOTYPE.
Parameters
[in]NSIZES
          NSIZES is INTEGER
          The number of values of M and N contained in the vectors
          MVAL and NVAL.  The matrix sizes are used in pairs (M,N).
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[in]NVAL
          NVAL is INTEGER array, dimension (NM)
          The values of the matrix column dimension N.
[in]NTYPES
          NTYPES is INTEGER
          The number of elements in DOTYPE.   If it is zero, SCHKBD
          does nothing.  It must be at least zero.  If it is MAXTYP+1
          and NSIZES is 1, then an additional type, MAXTYP+1 is
          defined, which is to use whatever matrices are in A and B.
          This is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
          DOTYPE(MAXTYP+1) is .TRUE. .
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          If DOTYPE(j) is .TRUE., then for each size (m,n), a matrix
          of type j will be generated.  If NTYPES is smaller than the
          maximum number of types defined (PARAMETER MAXTYP), then
          types NTYPES+1 through MAXTYP will not be generated.  If
          NTYPES is larger than MAXTYP, DOTYPE(MAXTYP+1) through
          DOTYPE(NTYPES) will be ignored.
[in]NRHS
          NRHS is INTEGER
          The number of columns in the "right-hand side" matrices X, Y,
          and Z, used in testing SBDSQR.  If NRHS = 0, then the
          operations on the right-hand side will not be tested.
          NRHS must be at least 0.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          On entry ISEED specifies the seed of the random number
          generator. The array elements should be between 0 and 4095;
          if not they will be reduced mod 4096.  Also, ISEED(4) must
          be odd.  The values of ISEED are changed on exit, and can be
          used in the next call to SCHKBD to continue the same random
          number sequence.
[in]THRESH
          THRESH is REAL
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.  Note that the
          expected value of the test ratios is O(1), so THRESH should
          be a reasonably small multiple of 1, e.g., 10 or 100.
[out]A
          A is REAL array, dimension (LDA,NMAX)
          where NMAX is the maximum value of N in NVAL.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,MMAX),
          where MMAX is the maximum value of M in MVAL.
[out]BD
          BD is REAL array, dimension
                      (max(min(MVAL(j),NVAL(j))))
[out]BE
          BE is REAL array, dimension
                      (max(min(MVAL(j),NVAL(j))))
[out]S1
          S1 is REAL array, dimension
                      (max(min(MVAL(j),NVAL(j))))
[out]S2
          S2 is REAL array, dimension
                      (max(min(MVAL(j),NVAL(j))))
[out]X
          X is REAL array, dimension (LDX,NRHS)
[in]LDX
          LDX is INTEGER
          The leading dimension of the arrays X, Y, and Z.
          LDX >= max(1,MMAX)
[out]Y
          Y is REAL array, dimension (LDX,NRHS)
[out]Z
          Z is REAL array, dimension (LDX,NRHS)
[out]Q
          Q is REAL array, dimension (LDQ,MMAX)
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  LDQ >= max(1,MMAX).
[out]PT
          PT is REAL array, dimension (LDPT,NMAX)
[in]LDPT
          LDPT is INTEGER
          The leading dimension of the arrays PT, U, and V.
          LDPT >= max(1, max(min(MVAL(j),NVAL(j)))).
[out]U
          U is REAL array, dimension
                      (LDPT,max(min(MVAL(j),NVAL(j))))
[out]VT
          VT is REAL array, dimension
                      (LDPT,max(min(MVAL(j),NVAL(j))))
[out]WORK
          WORK is REAL array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The number of entries in WORK.  This must be at least
          3(M+N) and  M(M + max(M,N,k) + 1) + N*min(M,N)  for all
          pairs  (M,N)=(MM(j),NN(j))
[out]IWORK
          IWORK is INTEGER array, dimension at least 8*min(M,N)
[in]NOUT
          NOUT is INTEGER
          The FORTRAN unit number for printing out error messages
          (e.g., if a routine returns IINFO not equal to 0.)
[out]INFO
          INFO is INTEGER
          If 0, then everything ran OK.
           -1: NSIZES < 0
           -2: Some MM(j) < 0
           -3: Some NN(j) < 0
           -4: NTYPES < 0
           -6: NRHS  < 0
           -8: THRESH < 0
          -11: LDA < 1 or LDA < MMAX, where MMAX is max( MM(j) ).
          -17: LDB < 1 or LDB < MMAX.
          -21: LDQ < 1 or LDQ < MMAX.
          -23: LDPT< 1 or LDPT< MNMAX.
          -27: LWORK too small.
          If  SLATMR, SLATMS, SGEBRD, SORGBR, or SBDSQR,
              returns an error code, the
              absolute value of it is returned.

-----------------------------------------------------------------------

     Some Local Variables and Parameters:
     ---- ----- --------- --- ----------

     ZERO, ONE       Real 0 and 1.
     MAXTYP          The number of types defined.
     NTEST           The number of tests performed, or which can
                     be performed so far, for the current matrix.
     MMAX            Largest value in NN.
     NMAX            Largest value in NN.
     MNMIN           min(MM(j), NN(j)) (the dimension of the bidiagonal
                     matrix.)
     MNMAX           The maximum value of MNMIN for j=1,...,NSIZES.
     NFAIL           The number of tests which have exceeded THRESH
     COND, IMODE     Values to be passed to the matrix generators.
     ANORM           Norm of A; passed to matrix generators.

     OVFL, UNFL      Overflow and underflow thresholds.
     RTOVFL, RTUNFL  Square roots of the previous 2 values.
     ULP, ULPINV     Finest relative precision and its inverse.

             The following four arrays decode JTYPE:
     KTYPE(j)        The general type (1-10) for type "j".
     KMODE(j)        The MODE value to be passed to the matrix
                     generator for type "j".
     KMAGN(j)        The order of magnitude ( O(1),
                     O(overflow^(1/2) ), O(underflow^(1/2) )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 489 of file schkbd.f.

493 *
494 * -- LAPACK test routine --
495 * -- LAPACK is a software package provided by Univ. of Tennessee, --
496 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
497 *
498 * .. Scalar Arguments ..
499  INTEGER INFO, LDA, LDPT, LDQ, LDX, LWORK, NOUT, NRHS,
500  $ NSIZES, NTYPES
501  REAL THRESH
502 * ..
503 * .. Array Arguments ..
504  LOGICAL DOTYPE( * )
505  INTEGER ISEED( 4 ), IWORK( * ), MVAL( * ), NVAL( * )
506  REAL A( LDA, * ), BD( * ), BE( * ), PT( LDPT, * ),
507  $ Q( LDQ, * ), S1( * ), S2( * ), U( LDPT, * ),
508  $ VT( LDPT, * ), WORK( * ), X( LDX, * ),
509  $ Y( LDX, * ), Z( LDX, * )
510 * ..
511 *
512 * ======================================================================
513 *
514 * .. Parameters ..
515  REAL ZERO, ONE, TWO, HALF
516  parameter( zero = 0.0e0, one = 1.0e0, two = 2.0e0,
517  $ half = 0.5e0 )
518  INTEGER MAXTYP
519  parameter( maxtyp = 16 )
520 * ..
521 * .. Local Scalars ..
522  LOGICAL BADMM, BADNN, BIDIAG
523  CHARACTER UPLO
524  CHARACTER*3 PATH
525  INTEGER I, IINFO, IL, IMODE, ITEMP, ITYPE, IU, IWBD,
526  $ IWBE, IWBS, IWBZ, IWWORK, J, JCOL, JSIZE,
527  $ JTYPE, LOG2UI, M, MINWRK, MMAX, MNMAX, MNMIN,
528  $ MNMIN2, MQ, MTYPES, N, NFAIL, NMAX,
529  $ NS1, NS2, NTEST
530  REAL ABSTOL, AMNINV, ANORM, COND, OVFL, RTOVFL,
531  $ RTUNFL, TEMP1, TEMP2, ULP, ULPINV, UNFL,
532  $ VL, VU
533 * ..
534 * .. Local Arrays ..
535  INTEGER IDUM( 1 ), IOLDSD( 4 ), ISEED2( 4 ),
536  $ KMAGN( MAXTYP ), KMODE( MAXTYP ),
537  $ KTYPE( MAXTYP )
538  REAL DUM( 1 ), DUMMA( 1 ), RESULT( 40 )
539 * ..
540 * .. External Functions ..
541  REAL SLAMCH, SLARND, SSXT1
542  EXTERNAL slamch, slarnd, ssxt1
543 * ..
544 * .. External Subroutines ..
545  EXTERNAL alasum, sbdsdc, sbdsqr, sbdsvdx, sbdt01,
549 * ..
550 * .. Intrinsic Functions ..
551  INTRINSIC abs, exp, int, log, max, min, sqrt
552 * ..
553 * .. Scalars in Common ..
554  LOGICAL LERR, OK
555  CHARACTER*32 SRNAMT
556  INTEGER INFOT, NUNIT
557 * ..
558 * .. Common blocks ..
559  COMMON / infoc / infot, nunit, ok, lerr
560  COMMON / srnamc / srnamt
561 * ..
562 * .. Data statements ..
563  DATA ktype / 1, 2, 5*4, 5*6, 3*9, 10 /
564  DATA kmagn / 2*1, 3*1, 2, 3, 3*1, 2, 3, 1, 2, 3, 0 /
565  DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
566  $ 0, 0, 0 /
567 * ..
568 * .. Executable Statements ..
569 *
570 * Check for errors
571 *
572  info = 0
573 *
574  badmm = .false.
575  badnn = .false.
576  mmax = 1
577  nmax = 1
578  mnmax = 1
579  minwrk = 1
580  DO 10 j = 1, nsizes
581  mmax = max( mmax, mval( j ) )
582  IF( mval( j ).LT.0 )
583  $ badmm = .true.
584  nmax = max( nmax, nval( j ) )
585  IF( nval( j ).LT.0 )
586  $ badnn = .true.
587  mnmax = max( mnmax, min( mval( j ), nval( j ) ) )
588  minwrk = max( minwrk, 3*( mval( j )+nval( j ) ),
589  $ mval( j )*( mval( j )+max( mval( j ), nval( j ),
590  $ nrhs )+1 )+nval( j )*min( nval( j ), mval( j ) ) )
591  10 CONTINUE
592 *
593 * Check for errors
594 *
595  IF( nsizes.LT.0 ) THEN
596  info = -1
597  ELSE IF( badmm ) THEN
598  info = -2
599  ELSE IF( badnn ) THEN
600  info = -3
601  ELSE IF( ntypes.LT.0 ) THEN
602  info = -4
603  ELSE IF( nrhs.LT.0 ) THEN
604  info = -6
605  ELSE IF( lda.LT.mmax ) THEN
606  info = -11
607  ELSE IF( ldx.LT.mmax ) THEN
608  info = -17
609  ELSE IF( ldq.LT.mmax ) THEN
610  info = -21
611  ELSE IF( ldpt.LT.mnmax ) THEN
612  info = -23
613  ELSE IF( minwrk.GT.lwork ) THEN
614  info = -27
615  END IF
616 *
617  IF( info.NE.0 ) THEN
618  CALL xerbla( 'SCHKBD', -info )
619  RETURN
620  END IF
621 *
622 * Initialize constants
623 *
624  path( 1: 1 ) = 'Single precision'
625  path( 2: 3 ) = 'BD'
626  nfail = 0
627  ntest = 0
628  unfl = slamch( 'Safe minimum' )
629  ovfl = slamch( 'Overflow' )
630  CALL slabad( unfl, ovfl )
631  ulp = slamch( 'Precision' )
632  ulpinv = one / ulp
633  log2ui = int( log( ulpinv ) / log( two ) )
634  rtunfl = sqrt( unfl )
635  rtovfl = sqrt( ovfl )
636  infot = 0
637  abstol = 2*unfl
638 *
639 * Loop over sizes, types
640 *
641  DO 300 jsize = 1, nsizes
642  m = mval( jsize )
643  n = nval( jsize )
644  mnmin = min( m, n )
645  amninv = one / max( m, n, 1 )
646 *
647  IF( nsizes.NE.1 ) THEN
648  mtypes = min( maxtyp, ntypes )
649  ELSE
650  mtypes = min( maxtyp+1, ntypes )
651  END IF
652 *
653  DO 290 jtype = 1, mtypes
654  IF( .NOT.dotype( jtype ) )
655  $ GO TO 290
656 *
657  DO 20 j = 1, 4
658  ioldsd( j ) = iseed( j )
659  20 CONTINUE
660 *
661  DO 30 j = 1, 34
662  result( j ) = -one
663  30 CONTINUE
664 *
665  uplo = ' '
666 *
667 * Compute "A"
668 *
669 * Control parameters:
670 *
671 * KMAGN KMODE KTYPE
672 * =1 O(1) clustered 1 zero
673 * =2 large clustered 2 identity
674 * =3 small exponential (none)
675 * =4 arithmetic diagonal, (w/ eigenvalues)
676 * =5 random symmetric, w/ eigenvalues
677 * =6 nonsymmetric, w/ singular values
678 * =7 random diagonal
679 * =8 random symmetric
680 * =9 random nonsymmetric
681 * =10 random bidiagonal (log. distrib.)
682 *
683  IF( mtypes.GT.maxtyp )
684  $ GO TO 100
685 *
686  itype = ktype( jtype )
687  imode = kmode( jtype )
688 *
689 * Compute norm
690 *
691  GO TO ( 40, 50, 60 )kmagn( jtype )
692 *
693  40 CONTINUE
694  anorm = one
695  GO TO 70
696 *
697  50 CONTINUE
698  anorm = ( rtovfl*ulp )*amninv
699  GO TO 70
700 *
701  60 CONTINUE
702  anorm = rtunfl*max( m, n )*ulpinv
703  GO TO 70
704 *
705  70 CONTINUE
706 *
707  CALL slaset( 'Full', lda, n, zero, zero, a, lda )
708  iinfo = 0
709  cond = ulpinv
710 *
711  bidiag = .false.
712  IF( itype.EQ.1 ) THEN
713 *
714 * Zero matrix
715 *
716  iinfo = 0
717 *
718  ELSE IF( itype.EQ.2 ) THEN
719 *
720 * Identity
721 *
722  DO 80 jcol = 1, mnmin
723  a( jcol, jcol ) = anorm
724  80 CONTINUE
725 *
726  ELSE IF( itype.EQ.4 ) THEN
727 *
728 * Diagonal Matrix, [Eigen]values Specified
729 *
730  CALL slatms( mnmin, mnmin, 'S', iseed, 'N', work, imode,
731  $ cond, anorm, 0, 0, 'N', a, lda,
732  $ work( mnmin+1 ), iinfo )
733 *
734  ELSE IF( itype.EQ.5 ) THEN
735 *
736 * Symmetric, eigenvalues specified
737 *
738  CALL slatms( mnmin, mnmin, 'S', iseed, 'S', work, imode,
739  $ cond, anorm, m, n, 'N', a, lda,
740  $ work( mnmin+1 ), iinfo )
741 *
742  ELSE IF( itype.EQ.6 ) THEN
743 *
744 * Nonsymmetric, singular values specified
745 *
746  CALL slatms( m, n, 'S', iseed, 'N', work, imode, cond,
747  $ anorm, m, n, 'N', a, lda, work( mnmin+1 ),
748  $ iinfo )
749 *
750  ELSE IF( itype.EQ.7 ) THEN
751 *
752 * Diagonal, random entries
753 *
754  CALL slatmr( mnmin, mnmin, 'S', iseed, 'N', work, 6, one,
755  $ one, 'T', 'N', work( mnmin+1 ), 1, one,
756  $ work( 2*mnmin+1 ), 1, one, 'N', iwork, 0, 0,
757  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
758 *
759  ELSE IF( itype.EQ.8 ) THEN
760 *
761 * Symmetric, random entries
762 *
763  CALL slatmr( mnmin, mnmin, 'S', iseed, 'S', work, 6, one,
764  $ one, 'T', 'N', work( mnmin+1 ), 1, one,
765  $ work( m+mnmin+1 ), 1, one, 'N', iwork, m, n,
766  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
767 *
768  ELSE IF( itype.EQ.9 ) THEN
769 *
770 * Nonsymmetric, random entries
771 *
772  CALL slatmr( m, n, 'S', iseed, 'N', work, 6, one, one,
773  $ 'T', 'N', work( mnmin+1 ), 1, one,
774  $ work( m+mnmin+1 ), 1, one, 'N', iwork, m, n,
775  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
776 *
777  ELSE IF( itype.EQ.10 ) THEN
778 *
779 * Bidiagonal, random entries
780 *
781  temp1 = -two*log( ulp )
782  DO 90 j = 1, mnmin
783  bd( j ) = exp( temp1*slarnd( 2, iseed ) )
784  IF( j.LT.mnmin )
785  $ be( j ) = exp( temp1*slarnd( 2, iseed ) )
786  90 CONTINUE
787 *
788  iinfo = 0
789  bidiag = .true.
790  IF( m.GE.n ) THEN
791  uplo = 'U'
792  ELSE
793  uplo = 'L'
794  END IF
795  ELSE
796  iinfo = 1
797  END IF
798 *
799  IF( iinfo.EQ.0 ) THEN
800 *
801 * Generate Right-Hand Side
802 *
803  IF( bidiag ) THEN
804  CALL slatmr( mnmin, nrhs, 'S', iseed, 'N', work, 6,
805  $ one, one, 'T', 'N', work( mnmin+1 ), 1,
806  $ one, work( 2*mnmin+1 ), 1, one, 'N',
807  $ iwork, mnmin, nrhs, zero, one, 'NO', y,
808  $ ldx, iwork, iinfo )
809  ELSE
810  CALL slatmr( m, nrhs, 'S', iseed, 'N', work, 6, one,
811  $ one, 'T', 'N', work( m+1 ), 1, one,
812  $ work( 2*m+1 ), 1, one, 'N', iwork, m,
813  $ nrhs, zero, one, 'NO', x, ldx, iwork,
814  $ iinfo )
815  END IF
816  END IF
817 *
818 * Error Exit
819 *
820  IF( iinfo.NE.0 ) THEN
821  WRITE( nout, fmt = 9998 )'Generator', iinfo, m, n,
822  $ jtype, ioldsd
823  info = abs( iinfo )
824  RETURN
825  END IF
826 *
827  100 CONTINUE
828 *
829 * Call SGEBRD and SORGBR to compute B, Q, and P, do tests.
830 *
831  IF( .NOT.bidiag ) THEN
832 *
833 * Compute transformations to reduce A to bidiagonal form:
834 * B := Q' * A * P.
835 *
836  CALL slacpy( ' ', m, n, a, lda, q, ldq )
837  CALL sgebrd( m, n, q, ldq, bd, be, work, work( mnmin+1 ),
838  $ work( 2*mnmin+1 ), lwork-2*mnmin, iinfo )
839 *
840 * Check error code from SGEBRD.
841 *
842  IF( iinfo.NE.0 ) THEN
843  WRITE( nout, fmt = 9998 )'SGEBRD', iinfo, m, n,
844  $ jtype, ioldsd
845  info = abs( iinfo )
846  RETURN
847  END IF
848 *
849  CALL slacpy( ' ', m, n, q, ldq, pt, ldpt )
850  IF( m.GE.n ) THEN
851  uplo = 'U'
852  ELSE
853  uplo = 'L'
854  END IF
855 *
856 * Generate Q
857 *
858  mq = m
859  IF( nrhs.LE.0 )
860  $ mq = mnmin
861  CALL sorgbr( 'Q', m, mq, n, q, ldq, work,
862  $ work( 2*mnmin+1 ), lwork-2*mnmin, iinfo )
863 *
864 * Check error code from SORGBR.
865 *
866  IF( iinfo.NE.0 ) THEN
867  WRITE( nout, fmt = 9998 )'SORGBR(Q)', iinfo, m, n,
868  $ jtype, ioldsd
869  info = abs( iinfo )
870  RETURN
871  END IF
872 *
873 * Generate P'
874 *
875  CALL sorgbr( 'P', mnmin, n, m, pt, ldpt, work( mnmin+1 ),
876  $ work( 2*mnmin+1 ), lwork-2*mnmin, iinfo )
877 *
878 * Check error code from SORGBR.
879 *
880  IF( iinfo.NE.0 ) THEN
881  WRITE( nout, fmt = 9998 )'SORGBR(P)', iinfo, m, n,
882  $ jtype, ioldsd
883  info = abs( iinfo )
884  RETURN
885  END IF
886 *
887 * Apply Q' to an M by NRHS matrix X: Y := Q' * X.
888 *
889  CALL sgemm( 'Transpose', 'No transpose', m, nrhs, m, one,
890  $ q, ldq, x, ldx, zero, y, ldx )
891 *
892 * Test 1: Check the decomposition A := Q * B * PT
893 * 2: Check the orthogonality of Q
894 * 3: Check the orthogonality of PT
895 *
896  CALL sbdt01( m, n, 1, a, lda, q, ldq, bd, be, pt, ldpt,
897  $ work, result( 1 ) )
898  CALL sort01( 'Columns', m, mq, q, ldq, work, lwork,
899  $ result( 2 ) )
900  CALL sort01( 'Rows', mnmin, n, pt, ldpt, work, lwork,
901  $ result( 3 ) )
902  END IF
903 *
904 * Use SBDSQR to form the SVD of the bidiagonal matrix B:
905 * B := U * S1 * VT, and compute Z = U' * Y.
906 *
907  CALL scopy( mnmin, bd, 1, s1, 1 )
908  IF( mnmin.GT.0 )
909  $ CALL scopy( mnmin-1, be, 1, work, 1 )
910  CALL slacpy( ' ', m, nrhs, y, ldx, z, ldx )
911  CALL slaset( 'Full', mnmin, mnmin, zero, one, u, ldpt )
912  CALL slaset( 'Full', mnmin, mnmin, zero, one, vt, ldpt )
913 *
914  CALL sbdsqr( uplo, mnmin, mnmin, mnmin, nrhs, s1, work, vt,
915  $ ldpt, u, ldpt, z, ldx, work( mnmin+1 ), iinfo )
916 *
917 * Check error code from SBDSQR.
918 *
919  IF( iinfo.NE.0 ) THEN
920  WRITE( nout, fmt = 9998 )'SBDSQR(vects)', iinfo, m, n,
921  $ jtype, ioldsd
922  info = abs( iinfo )
923  IF( iinfo.LT.0 ) THEN
924  RETURN
925  ELSE
926  result( 4 ) = ulpinv
927  GO TO 270
928  END IF
929  END IF
930 *
931 * Use SBDSQR to compute only the singular values of the
932 * bidiagonal matrix B; U, VT, and Z should not be modified.
933 *
934  CALL scopy( mnmin, bd, 1, s2, 1 )
935  IF( mnmin.GT.0 )
936  $ CALL scopy( mnmin-1, be, 1, work, 1 )
937 *
938  CALL sbdsqr( uplo, mnmin, 0, 0, 0, s2, work, vt, ldpt, u,
939  $ ldpt, z, ldx, work( mnmin+1 ), iinfo )
940 *
941 * Check error code from SBDSQR.
942 *
943  IF( iinfo.NE.0 ) THEN
944  WRITE( nout, fmt = 9998 )'SBDSQR(values)', iinfo, m, n,
945  $ jtype, ioldsd
946  info = abs( iinfo )
947  IF( iinfo.LT.0 ) THEN
948  RETURN
949  ELSE
950  result( 9 ) = ulpinv
951  GO TO 270
952  END IF
953  END IF
954 *
955 * Test 4: Check the decomposition B := U * S1 * VT
956 * 5: Check the computation Z := U' * Y
957 * 6: Check the orthogonality of U
958 * 7: Check the orthogonality of VT
959 *
960  CALL sbdt03( uplo, mnmin, 1, bd, be, u, ldpt, s1, vt, ldpt,
961  $ work, result( 4 ) )
962  CALL sbdt02( mnmin, nrhs, y, ldx, z, ldx, u, ldpt, work,
963  $ result( 5 ) )
964  CALL sort01( 'Columns', mnmin, mnmin, u, ldpt, work, lwork,
965  $ result( 6 ) )
966  CALL sort01( 'Rows', mnmin, mnmin, vt, ldpt, work, lwork,
967  $ result( 7 ) )
968 *
969 * Test 8: Check that the singular values are sorted in
970 * non-increasing order and are non-negative
971 *
972  result( 8 ) = zero
973  DO 110 i = 1, mnmin - 1
974  IF( s1( i ).LT.s1( i+1 ) )
975  $ result( 8 ) = ulpinv
976  IF( s1( i ).LT.zero )
977  $ result( 8 ) = ulpinv
978  110 CONTINUE
979  IF( mnmin.GE.1 ) THEN
980  IF( s1( mnmin ).LT.zero )
981  $ result( 8 ) = ulpinv
982  END IF
983 *
984 * Test 9: Compare SBDSQR with and without singular vectors
985 *
986  temp2 = zero
987 *
988  DO 120 j = 1, mnmin
989  temp1 = abs( s1( j )-s2( j ) ) /
990  $ max( sqrt( unfl )*max( s1( 1 ), one ),
991  $ ulp*max( abs( s1( j ) ), abs( s2( j ) ) ) )
992  temp2 = max( temp1, temp2 )
993  120 CONTINUE
994 *
995  result( 9 ) = temp2
996 *
997 * Test 10: Sturm sequence test of singular values
998 * Go up by factors of two until it succeeds
999 *
1000  temp1 = thresh*( half-ulp )
1001 *
1002  DO 130 j = 0, log2ui
1003 * CALL SSVDCH( MNMIN, BD, BE, S1, TEMP1, IINFO )
1004  IF( iinfo.EQ.0 )
1005  $ GO TO 140
1006  temp1 = temp1*two
1007  130 CONTINUE
1008 *
1009  140 CONTINUE
1010  result( 10 ) = temp1
1011 *
1012 * Use SBDSQR to form the decomposition A := (QU) S (VT PT)
1013 * from the bidiagonal form A := Q B PT.
1014 *
1015  IF( .NOT.bidiag ) THEN
1016  CALL scopy( mnmin, bd, 1, s2, 1 )
1017  IF( mnmin.GT.0 )
1018  $ CALL scopy( mnmin-1, be, 1, work, 1 )
1019 *
1020  CALL sbdsqr( uplo, mnmin, n, m, nrhs, s2, work, pt, ldpt,
1021  $ q, ldq, y, ldx, work( mnmin+1 ), iinfo )
1022 *
1023 * Test 11: Check the decomposition A := Q*U * S2 * VT*PT
1024 * 12: Check the computation Z := U' * Q' * X
1025 * 13: Check the orthogonality of Q*U
1026 * 14: Check the orthogonality of VT*PT
1027 *
1028  CALL sbdt01( m, n, 0, a, lda, q, ldq, s2, dumma, pt,
1029  $ ldpt, work, result( 11 ) )
1030  CALL sbdt02( m, nrhs, x, ldx, y, ldx, q, ldq, work,
1031  $ result( 12 ) )
1032  CALL sort01( 'Columns', m, mq, q, ldq, work, lwork,
1033  $ result( 13 ) )
1034  CALL sort01( 'Rows', mnmin, n, pt, ldpt, work, lwork,
1035  $ result( 14 ) )
1036  END IF
1037 *
1038 * Use SBDSDC to form the SVD of the bidiagonal matrix B:
1039 * B := U * S1 * VT
1040 *
1041  CALL scopy( mnmin, bd, 1, s1, 1 )
1042  IF( mnmin.GT.0 )
1043  $ CALL scopy( mnmin-1, be, 1, work, 1 )
1044  CALL slaset( 'Full', mnmin, mnmin, zero, one, u, ldpt )
1045  CALL slaset( 'Full', mnmin, mnmin, zero, one, vt, ldpt )
1046 *
1047  CALL sbdsdc( uplo, 'I', mnmin, s1, work, u, ldpt, vt, ldpt,
1048  $ dum, idum, work( mnmin+1 ), iwork, iinfo )
1049 *
1050 * Check error code from SBDSDC.
1051 *
1052  IF( iinfo.NE.0 ) THEN
1053  WRITE( nout, fmt = 9998 )'SBDSDC(vects)', iinfo, m, n,
1054  $ jtype, ioldsd
1055  info = abs( iinfo )
1056  IF( iinfo.LT.0 ) THEN
1057  RETURN
1058  ELSE
1059  result( 15 ) = ulpinv
1060  GO TO 270
1061  END IF
1062  END IF
1063 *
1064 * Use SBDSDC to compute only the singular values of the
1065 * bidiagonal matrix B; U and VT should not be modified.
1066 *
1067  CALL scopy( mnmin, bd, 1, s2, 1 )
1068  IF( mnmin.GT.0 )
1069  $ CALL scopy( mnmin-1, be, 1, work, 1 )
1070 *
1071  CALL sbdsdc( uplo, 'N', mnmin, s2, work, dum, 1, dum, 1,
1072  $ dum, idum, work( mnmin+1 ), iwork, iinfo )
1073 *
1074 * Check error code from SBDSDC.
1075 *
1076  IF( iinfo.NE.0 ) THEN
1077  WRITE( nout, fmt = 9998 )'SBDSDC(values)', iinfo, m, n,
1078  $ jtype, ioldsd
1079  info = abs( iinfo )
1080  IF( iinfo.LT.0 ) THEN
1081  RETURN
1082  ELSE
1083  result( 18 ) = ulpinv
1084  GO TO 270
1085  END IF
1086  END IF
1087 *
1088 * Test 15: Check the decomposition B := U * S1 * VT
1089 * 16: Check the orthogonality of U
1090 * 17: Check the orthogonality of VT
1091 *
1092  CALL sbdt03( uplo, mnmin, 1, bd, be, u, ldpt, s1, vt, ldpt,
1093  $ work, result( 15 ) )
1094  CALL sort01( 'Columns', mnmin, mnmin, u, ldpt, work, lwork,
1095  $ result( 16 ) )
1096  CALL sort01( 'Rows', mnmin, mnmin, vt, ldpt, work, lwork,
1097  $ result( 17 ) )
1098 *
1099 * Test 18: Check that the singular values are sorted in
1100 * non-increasing order and are non-negative
1101 *
1102  result( 18 ) = zero
1103  DO 150 i = 1, mnmin - 1
1104  IF( s1( i ).LT.s1( i+1 ) )
1105  $ result( 18 ) = ulpinv
1106  IF( s1( i ).LT.zero )
1107  $ result( 18 ) = ulpinv
1108  150 CONTINUE
1109  IF( mnmin.GE.1 ) THEN
1110  IF( s1( mnmin ).LT.zero )
1111  $ result( 18 ) = ulpinv
1112  END IF
1113 *
1114 * Test 19: Compare SBDSQR with and without singular vectors
1115 *
1116  temp2 = zero
1117 *
1118  DO 160 j = 1, mnmin
1119  temp1 = abs( s1( j )-s2( j ) ) /
1120  $ max( sqrt( unfl )*max( s1( 1 ), one ),
1121  $ ulp*max( abs( s1( 1 ) ), abs( s2( 1 ) ) ) )
1122  temp2 = max( temp1, temp2 )
1123  160 CONTINUE
1124 *
1125  result( 19 ) = temp2
1126 *
1127 *
1128 * Use SBDSVDX to compute the SVD of the bidiagonal matrix B:
1129 * B := U * S1 * VT
1130 *
1131  IF( jtype.EQ.10 .OR. jtype.EQ.16 ) THEN
1132 * =================================
1133 * Matrix types temporarily disabled
1134 * =================================
1135  result( 20:34 ) = zero
1136  GO TO 270
1137  END IF
1138 *
1139  iwbs = 1
1140  iwbd = iwbs + mnmin
1141  iwbe = iwbd + mnmin
1142  iwbz = iwbe + mnmin
1143  iwwork = iwbz + 2*mnmin*(mnmin+1)
1144  mnmin2 = max( 1,mnmin*2 )
1145 *
1146  CALL scopy( mnmin, bd, 1, work( iwbd ), 1 )
1147  IF( mnmin.GT.0 )
1148  $ CALL scopy( mnmin-1, be, 1, work( iwbe ), 1 )
1149 *
1150  CALL sbdsvdx( uplo, 'V', 'A', mnmin, work( iwbd ),
1151  $ work( iwbe ), zero, zero, 0, 0, ns1, s1,
1152  $ work( iwbz ), mnmin2, work( iwwork ),
1153  $ iwork, iinfo)
1154 *
1155 * Check error code from SBDSVDX.
1156 *
1157  IF( iinfo.NE.0 ) THEN
1158  WRITE( nout, fmt = 9998 )'SBDSVDX(vects,A)', iinfo, m, n,
1159  $ jtype, ioldsd
1160  info = abs( iinfo )
1161  IF( iinfo.LT.0 ) THEN
1162  RETURN
1163  ELSE
1164  result( 20 ) = ulpinv
1165  GO TO 270
1166  END IF
1167  END IF
1168 *
1169  j = iwbz
1170  DO 170 i = 1, ns1
1171  CALL scopy( mnmin, work( j ), 1, u( 1,i ), 1 )
1172  j = j + mnmin
1173  CALL scopy( mnmin, work( j ), 1, vt( i,1 ), ldpt )
1174  j = j + mnmin
1175  170 CONTINUE
1176 *
1177 * Use SBDSVDX to compute only the singular values of the
1178 * bidiagonal matrix B; U and VT should not be modified.
1179 *
1180  IF( jtype.EQ.9 ) THEN
1181 * =================================
1182 * Matrix types temporarily disabled
1183 * =================================
1184  result( 24 ) = zero
1185  GO TO 270
1186  END IF
1187 *
1188  CALL scopy( mnmin, bd, 1, work( iwbd ), 1 )
1189  IF( mnmin.GT.0 )
1190  $ CALL scopy( mnmin-1, be, 1, work( iwbe ), 1 )
1191 *
1192  CALL sbdsvdx( uplo, 'N', 'A', mnmin, work( iwbd ),
1193  $ work( iwbe ), zero, zero, 0, 0, ns2, s2,
1194  $ work( iwbz ), mnmin2, work( iwwork ),
1195  $ iwork, iinfo )
1196 *
1197 * Check error code from SBDSVDX.
1198 *
1199  IF( iinfo.NE.0 ) THEN
1200  WRITE( nout, fmt = 9998 )'SBDSVDX(values,A)', iinfo,
1201  $ m, n, jtype, ioldsd
1202  info = abs( iinfo )
1203  IF( iinfo.LT.0 ) THEN
1204  RETURN
1205  ELSE
1206  result( 24 ) = ulpinv
1207  GO TO 270
1208  END IF
1209  END IF
1210 *
1211 * Save S1 for tests 30-34.
1212 *
1213  CALL scopy( mnmin, s1, 1, work( iwbs ), 1 )
1214 *
1215 * Test 20: Check the decomposition B := U * S1 * VT
1216 * 21: Check the orthogonality of U
1217 * 22: Check the orthogonality of VT
1218 * 23: Check that the singular values are sorted in
1219 * non-increasing order and are non-negative
1220 * 24: Compare SBDSVDX with and without singular vectors
1221 *
1222  CALL sbdt03( uplo, mnmin, 1, bd, be, u, ldpt, s1, vt,
1223  $ ldpt, work( iwbs+mnmin ), result( 20 ) )
1224  CALL sort01( 'Columns', mnmin, mnmin, u, ldpt,
1225  $ work( iwbs+mnmin ), lwork-mnmin,
1226  $ result( 21 ) )
1227  CALL sort01( 'Rows', mnmin, mnmin, vt, ldpt,
1228  $ work( iwbs+mnmin ), lwork-mnmin,
1229  $ result( 22) )
1230 *
1231  result( 23 ) = zero
1232  DO 180 i = 1, mnmin - 1
1233  IF( s1( i ).LT.s1( i+1 ) )
1234  $ result( 23 ) = ulpinv
1235  IF( s1( i ).LT.zero )
1236  $ result( 23 ) = ulpinv
1237  180 CONTINUE
1238  IF( mnmin.GE.1 ) THEN
1239  IF( s1( mnmin ).LT.zero )
1240  $ result( 23 ) = ulpinv
1241  END IF
1242 *
1243  temp2 = zero
1244  DO 190 j = 1, mnmin
1245  temp1 = abs( s1( j )-s2( j ) ) /
1246  $ max( sqrt( unfl )*max( s1( 1 ), one ),
1247  $ ulp*max( abs( s1( 1 ) ), abs( s2( 1 ) ) ) )
1248  temp2 = max( temp1, temp2 )
1249  190 CONTINUE
1250  result( 24 ) = temp2
1251  anorm = s1( 1 )
1252 *
1253 * Use SBDSVDX with RANGE='I': choose random values for IL and
1254 * IU, and ask for the IL-th through IU-th singular values
1255 * and corresponding vectors.
1256 *
1257  DO 200 i = 1, 4
1258  iseed2( i ) = iseed( i )
1259  200 CONTINUE
1260  IF( mnmin.LE.1 ) THEN
1261  il = 1
1262  iu = mnmin
1263  ELSE
1264  il = 1 + int( ( mnmin-1 )*slarnd( 1, iseed2 ) )
1265  iu = 1 + int( ( mnmin-1 )*slarnd( 1, iseed2 ) )
1266  IF( iu.LT.il ) THEN
1267  itemp = iu
1268  iu = il
1269  il = itemp
1270  END IF
1271  END IF
1272 *
1273  CALL scopy( mnmin, bd, 1, work( iwbd ), 1 )
1274  IF( mnmin.GT.0 )
1275  $ CALL scopy( mnmin-1, be, 1, work( iwbe ), 1 )
1276 *
1277  CALL sbdsvdx( uplo, 'V', 'I', mnmin, work( iwbd ),
1278  $ work( iwbe ), zero, zero, il, iu, ns1, s1,
1279  $ work( iwbz ), mnmin2, work( iwwork ),
1280  $ iwork, iinfo)
1281 *
1282 * Check error code from SBDSVDX.
1283 *
1284  IF( iinfo.NE.0 ) THEN
1285  WRITE( nout, fmt = 9998 )'SBDSVDX(vects,I)', iinfo,
1286  $ m, n, jtype, ioldsd
1287  info = abs( iinfo )
1288  IF( iinfo.LT.0 ) THEN
1289  RETURN
1290  ELSE
1291  result( 25 ) = ulpinv
1292  GO TO 270
1293  END IF
1294  END IF
1295 *
1296  j = iwbz
1297  DO 210 i = 1, ns1
1298  CALL scopy( mnmin, work( j ), 1, u( 1,i ), 1 )
1299  j = j + mnmin
1300  CALL scopy( mnmin, work( j ), 1, vt( i,1 ), ldpt )
1301  j = j + mnmin
1302  210 CONTINUE
1303 *
1304 * Use SBDSVDX to compute only the singular values of the
1305 * bidiagonal matrix B; U and VT should not be modified.
1306 *
1307  CALL scopy( mnmin, bd, 1, work( iwbd ), 1 )
1308  IF( mnmin.GT.0 )
1309  $ CALL scopy( mnmin-1, be, 1, work( iwbe ), 1 )
1310 *
1311  CALL sbdsvdx( uplo, 'N', 'I', mnmin, work( iwbd ),
1312  $ work( iwbe ), zero, zero, il, iu, ns2, s2,
1313  $ work( iwbz ), mnmin2, work( iwwork ),
1314  $ iwork, iinfo )
1315 *
1316 * Check error code from SBDSVDX.
1317 *
1318  IF( iinfo.NE.0 ) THEN
1319  WRITE( nout, fmt = 9998 )'SBDSVDX(values,I)', iinfo,
1320  $ m, n, jtype, ioldsd
1321  info = abs( iinfo )
1322  IF( iinfo.LT.0 ) THEN
1323  RETURN
1324  ELSE
1325  result( 29 ) = ulpinv
1326  GO TO 270
1327  END IF
1328  END IF
1329 *
1330 * Test 25: Check S1 - U' * B * VT'
1331 * 26: Check the orthogonality of U
1332 * 27: Check the orthogonality of VT
1333 * 28: Check that the singular values are sorted in
1334 * non-increasing order and are non-negative
1335 * 29: Compare SBDSVDX with and without singular vectors
1336 *
1337  CALL sbdt04( uplo, mnmin, bd, be, s1, ns1, u,
1338  $ ldpt, vt, ldpt, work( iwbs+mnmin ),
1339  $ result( 25 ) )
1340  CALL sort01( 'Columns', mnmin, ns1, u, ldpt,
1341  $ work( iwbs+mnmin ), lwork-mnmin,
1342  $ result( 26 ) )
1343  CALL sort01( 'Rows', ns1, mnmin, vt, ldpt,
1344  $ work( iwbs+mnmin ), lwork-mnmin,
1345  $ result( 27 ) )
1346 *
1347  result( 28 ) = zero
1348  DO 220 i = 1, ns1 - 1
1349  IF( s1( i ).LT.s1( i+1 ) )
1350  $ result( 28 ) = ulpinv
1351  IF( s1( i ).LT.zero )
1352  $ result( 28 ) = ulpinv
1353  220 CONTINUE
1354  IF( ns1.GE.1 ) THEN
1355  IF( s1( ns1 ).LT.zero )
1356  $ result( 28 ) = ulpinv
1357  END IF
1358 *
1359  temp2 = zero
1360  DO 230 j = 1, ns1
1361  temp1 = abs( s1( j )-s2( j ) ) /
1362  $ max( sqrt( unfl )*max( s1( 1 ), one ),
1363  $ ulp*max( abs( s1( 1 ) ), abs( s2( 1 ) ) ) )
1364  temp2 = max( temp1, temp2 )
1365  230 CONTINUE
1366  result( 29 ) = temp2
1367 *
1368 * Use SBDSVDX with RANGE='V': determine the values VL and VU
1369 * of the IL-th and IU-th singular values and ask for all
1370 * singular values in this range.
1371 *
1372  CALL scopy( mnmin, work( iwbs ), 1, s1, 1 )
1373 *
1374  IF( mnmin.GT.0 ) THEN
1375  IF( il.NE.1 ) THEN
1376  vu = s1( il ) + max( half*abs( s1( il )-s1( il-1 ) ),
1377  $ ulp*anorm, two*rtunfl )
1378  ELSE
1379  vu = s1( 1 ) + max( half*abs( s1( mnmin )-s1( 1 ) ),
1380  $ ulp*anorm, two*rtunfl )
1381  END IF
1382  IF( iu.NE.ns1 ) THEN
1383  vl = s1( iu ) - max( ulp*anorm, two*rtunfl,
1384  $ half*abs( s1( iu+1 )-s1( iu ) ) )
1385  ELSE
1386  vl = s1( ns1 ) - max( ulp*anorm, two*rtunfl,
1387  $ half*abs( s1( mnmin )-s1( 1 ) ) )
1388  END IF
1389  vl = max( vl,zero )
1390  vu = max( vu,zero )
1391  IF( vl.GE.vu ) vu = max( vu*2, vu+vl+half )
1392  ELSE
1393  vl = zero
1394  vu = one
1395  END IF
1396 *
1397  CALL scopy( mnmin, bd, 1, work( iwbd ), 1 )
1398  IF( mnmin.GT.0 )
1399  $ CALL scopy( mnmin-1, be, 1, work( iwbe ), 1 )
1400 *
1401  CALL sbdsvdx( uplo, 'V', 'V', mnmin, work( iwbd ),
1402  $ work( iwbe ), vl, vu, 0, 0, ns1, s1,
1403  $ work( iwbz ), mnmin2, work( iwwork ),
1404  $ iwork, iinfo )
1405 *
1406 * Check error code from SBDSVDX.
1407 *
1408  IF( iinfo.NE.0 ) THEN
1409  WRITE( nout, fmt = 9998 )'SBDSVDX(vects,V)', iinfo,
1410  $ m, n, jtype, ioldsd
1411  info = abs( iinfo )
1412  IF( iinfo.LT.0 ) THEN
1413  RETURN
1414  ELSE
1415  result( 30 ) = ulpinv
1416  GO TO 270
1417  END IF
1418  END IF
1419 *
1420  j = iwbz
1421  DO 240 i = 1, ns1
1422  CALL scopy( mnmin, work( j ), 1, u( 1,i ), 1 )
1423  j = j + mnmin
1424  CALL scopy( mnmin, work( j ), 1, vt( i,1 ), ldpt )
1425  j = j + mnmin
1426  240 CONTINUE
1427 *
1428 * Use SBDSVDX to compute only the singular values of the
1429 * bidiagonal matrix B; U and VT should not be modified.
1430 *
1431  CALL scopy( mnmin, bd, 1, work( iwbd ), 1 )
1432  IF( mnmin.GT.0 )
1433  $ CALL scopy( mnmin-1, be, 1, work( iwbe ), 1 )
1434 *
1435  CALL sbdsvdx( uplo, 'N', 'V', mnmin, work( iwbd ),
1436  $ work( iwbe ), vl, vu, 0, 0, ns2, s2,
1437  $ work( iwbz ), mnmin2, work( iwwork ),
1438  $ iwork, iinfo )
1439 *
1440 * Check error code from SBDSVDX.
1441 *
1442  IF( iinfo.NE.0 ) THEN
1443  WRITE( nout, fmt = 9998 )'SBDSVDX(values,V)', iinfo,
1444  $ m, n, jtype, ioldsd
1445  info = abs( iinfo )
1446  IF( iinfo.LT.0 ) THEN
1447  RETURN
1448  ELSE
1449  result( 34 ) = ulpinv
1450  GO TO 270
1451  END IF
1452  END IF
1453 *
1454 * Test 30: Check S1 - U' * B * VT'
1455 * 31: Check the orthogonality of U
1456 * 32: Check the orthogonality of VT
1457 * 33: Check that the singular values are sorted in
1458 * non-increasing order and are non-negative
1459 * 34: Compare SBDSVDX with and without singular vectors
1460 *
1461  CALL sbdt04( uplo, mnmin, bd, be, s1, ns1, u,
1462  $ ldpt, vt, ldpt, work( iwbs+mnmin ),
1463  $ result( 30 ) )
1464  CALL sort01( 'Columns', mnmin, ns1, u, ldpt,
1465  $ work( iwbs+mnmin ), lwork-mnmin,
1466  $ result( 31 ) )
1467  CALL sort01( 'Rows', ns1, mnmin, vt, ldpt,
1468  $ work( iwbs+mnmin ), lwork-mnmin,
1469  $ result( 32 ) )
1470 *
1471  result( 33 ) = zero
1472  DO 250 i = 1, ns1 - 1
1473  IF( s1( i ).LT.s1( i+1 ) )
1474  $ result( 28 ) = ulpinv
1475  IF( s1( i ).LT.zero )
1476  $ result( 28 ) = ulpinv
1477  250 CONTINUE
1478  IF( ns1.GE.1 ) THEN
1479  IF( s1( ns1 ).LT.zero )
1480  $ result( 28 ) = ulpinv
1481  END IF
1482 *
1483  temp2 = zero
1484  DO 260 j = 1, ns1
1485  temp1 = abs( s1( j )-s2( j ) ) /
1486  $ max( sqrt( unfl )*max( s1( 1 ), one ),
1487  $ ulp*max( abs( s1( 1 ) ), abs( s2( 1 ) ) ) )
1488  temp2 = max( temp1, temp2 )
1489  260 CONTINUE
1490  result( 34 ) = temp2
1491 *
1492 * End of Loop -- Check for RESULT(j) > THRESH
1493 *
1494  270 CONTINUE
1495 *
1496  DO 280 j = 1, 34
1497  IF( result( j ).GE.thresh ) THEN
1498  IF( nfail.EQ.0 )
1499  $ CALL slahd2( nout, path )
1500  WRITE( nout, fmt = 9999 )m, n, jtype, ioldsd, j,
1501  $ result( j )
1502  nfail = nfail + 1
1503  END IF
1504  280 CONTINUE
1505  IF( .NOT.bidiag ) THEN
1506  ntest = ntest + 34
1507  ELSE
1508  ntest = ntest + 30
1509  END IF
1510 *
1511  290 CONTINUE
1512  300 CONTINUE
1513 *
1514 * Summary
1515 *
1516  CALL alasum( path, nout, nfail, ntest, 0 )
1517 *
1518  RETURN
1519 *
1520 * End of SCHKBD
1521 *
1522  9999 FORMAT( ' M=', i5, ', N=', i5, ', type ', i2, ', seed=',
1523  $ 4( i4, ',' ), ' test(', i2, ')=', g11.4 )
1524  9998 FORMAT( ' SCHKBD: ', a, ' returned INFO=', i6, '.', / 9x, 'M=',
1525  $ i6, ', N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ),
1526  $ i5, ')' )
1527 *
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:74
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition: slaset.f:110
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:103
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:73
subroutine sbdsqr(UPLO, N, NCVT, NRU, NCC, D, E, VT, LDVT, U, LDU, C, LDC, WORK, INFO)
SBDSQR
Definition: sbdsqr.f:240
subroutine sbdsdc(UPLO, COMPQ, N, D, E, U, LDU, VT, LDVT, Q, IQ, WORK, IWORK, INFO)
SBDSDC
Definition: sbdsdc.f:205
subroutine sbdt04(UPLO, N, D, E, S, NS, U, LDU, VT, LDVT, WORK, RESID)
SBDT04
Definition: sbdt04.f:131
real function slarnd(IDIST, ISEED)
SLARND
Definition: slarnd.f:73
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:321
subroutine slatmr(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER, CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM, PACK, A, LDA, IWORK, INFO)
SLATMR
Definition: slatmr.f:471
subroutine sorgbr(VECT, M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
SORGBR
Definition: sorgbr.f:157
subroutine sgebrd(M, N, A, LDA, D, E, TAUQ, TAUP, WORK, LWORK, INFO)
SGEBRD
Definition: sgebrd.f:205
subroutine sbdsvdx(UPLO, JOBZ, RANGE, N, D, E, VL, VU, IL, IU, NS, S, Z, LDZ, WORK, IWORK, INFO)
SBDSVDX
Definition: sbdsvdx.f:226
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:82
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:187
subroutine sort01(ROWCOL, M, N, U, LDU, WORK, LWORK, RESID)
SORT01
Definition: sort01.f:116
subroutine sbdt02(M, N, B, LDB, C, LDC, U, LDU, WORK, RESID)
SBDT02
Definition: sbdt02.f:112
subroutine sbdt01(M, N, KD, A, LDA, Q, LDQ, D, E, PT, LDPT, WORK, RESID)
SBDT01
Definition: sbdt01.f:141
real function ssxt1(IJOB, D1, N1, D2, N2, ABSTOL, ULP, UNFL)
SSXT1
Definition: ssxt1.f:106
subroutine sbdt03(UPLO, N, KD, D, E, U, LDU, S, VT, LDVT, WORK, RESID)
SBDT03
Definition: sbdt03.f:135
subroutine slahd2(IOUNIT, PATH)
SLAHD2
Definition: slahd2.f:65
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:68
Here is the call graph for this function:
Here is the caller graph for this function: