LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ cchkbd()

subroutine cchkbd ( integer  nsizes,
integer, dimension( * )  mval,
integer, dimension( * )  nval,
integer  ntypes,
logical, dimension( * )  dotype,
integer  nrhs,
integer, dimension( 4 )  iseed,
real  thresh,
complex, dimension( lda, * )  a,
integer  lda,
real, dimension( * )  bd,
real, dimension( * )  be,
real, dimension( * )  s1,
real, dimension( * )  s2,
complex, dimension( ldx, * )  x,
integer  ldx,
complex, dimension( ldx, * )  y,
complex, dimension( ldx, * )  z,
complex, dimension( ldq, * )  q,
integer  ldq,
complex, dimension( ldpt, * )  pt,
integer  ldpt,
complex, dimension( ldpt, * )  u,
complex, dimension( ldpt, * )  vt,
complex, dimension( * )  work,
integer  lwork,
real, dimension( * )  rwork,
integer  nout,
integer  info 
)

CCHKBD

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

 CGEBRD reduces a complex general m by n matrix A to real upper or
 lower bidiagonal form 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.

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

 CBDSQR 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, CBDSQR has an option to apply the left orthogonal matrix
 U to a matrix X, useful in least squares applications.

 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 CGEBRD and CUNGBR

 (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 CBDSQR 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)   0 if the true singular values of B are within THRESH of
       those in S1.  2*THRESH if they are not.  (Tested using
       SSVDCH)

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

 Test CBDSQR 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 )

 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) CGEBRD 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, CCHKBD
          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 CBDSQR.  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 CCHKBD 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 COMPLEX 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 COMPLEX 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 COMPLEX array, dimension (LDX,NRHS)
[out]Z
          Z is COMPLEX array, dimension (LDX,NRHS)
[out]Q
          Q is COMPLEX array, dimension (LDQ,MMAX)
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  LDQ >= max(1,MMAX).
[out]PT
          PT is COMPLEX 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 COMPLEX array, dimension
                      (LDPT,max(min(MVAL(j),NVAL(j))))
[out]VT
          VT is COMPLEX array, dimension
                      (LDPT,max(min(MVAL(j),NVAL(j))))
[out]WORK
          WORK is COMPLEX 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]RWORK
          RWORK is REAL array, dimension
                      (5*max(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: LDP < 1 or LDP < MNMAX.
          -27: LWORK too small.
          If  CLATMR, CLATMS, CGEBRD, CUNGBR, or CBDSQR,
              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 411 of file cchkbd.f.

415*
416* -- LAPACK test routine --
417* -- LAPACK is a software package provided by Univ. of Tennessee, --
418* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
419*
420* .. Scalar Arguments ..
421 INTEGER INFO, LDA, LDPT, LDQ, LDX, LWORK, NOUT, NRHS,
422 $ NSIZES, NTYPES
423 REAL THRESH
424* ..
425* .. Array Arguments ..
426 LOGICAL DOTYPE( * )
427 INTEGER ISEED( 4 ), MVAL( * ), NVAL( * )
428 REAL BD( * ), BE( * ), RWORK( * ), S1( * ), S2( * )
429 COMPLEX A( LDA, * ), PT( LDPT, * ), Q( LDQ, * ),
430 $ U( LDPT, * ), VT( LDPT, * ), WORK( * ),
431 $ X( LDX, * ), Y( LDX, * ), Z( LDX, * )
432* ..
433*
434* ======================================================================
435*
436* .. Parameters ..
437 REAL ZERO, ONE, TWO, HALF
438 parameter( zero = 0.0e0, one = 1.0e0, two = 2.0e0,
439 $ half = 0.5e0 )
440 COMPLEX CZERO, CONE
441 parameter( czero = ( 0.0e+0, 0.0e+0 ),
442 $ cone = ( 1.0e+0, 0.0e+0 ) )
443 INTEGER MAXTYP
444 parameter( maxtyp = 16 )
445* ..
446* .. Local Scalars ..
447 LOGICAL BADMM, BADNN, BIDIAG
448 CHARACTER UPLO
449 CHARACTER*3 PATH
450 INTEGER I, IINFO, IMODE, ITYPE, J, JCOL, JSIZE, JTYPE,
451 $ LOG2UI, M, MINWRK, MMAX, MNMAX, MNMIN, MQ,
452 $ MTYPES, N, NFAIL, NMAX, NTEST
453 REAL AMNINV, ANORM, COND, OVFL, RTOVFL, RTUNFL,
454 $ TEMP1, TEMP2, ULP, ULPINV, UNFL
455* ..
456* .. Local Arrays ..
457 INTEGER IOLDSD( 4 ), IWORK( 1 ), KMAGN( MAXTYP ),
458 $ KMODE( MAXTYP ), KTYPE( MAXTYP )
459 REAL DUMMA( 1 ), RESULT( 14 )
460* ..
461* .. External Functions ..
462 REAL SLAMCH, SLARND
463 EXTERNAL slamch, slarnd
464* ..
465* .. External Subroutines ..
466 EXTERNAL alasum, cbdsqr, cbdt01, cbdt02, cbdt03,
469 $ ssvdch, xerbla
470* ..
471* .. Intrinsic Functions ..
472 INTRINSIC abs, exp, int, log, max, min, sqrt
473* ..
474* .. Scalars in Common ..
475 LOGICAL LERR, OK
476 CHARACTER*32 SRNAMT
477 INTEGER INFOT, NUNIT
478* ..
479* .. Common blocks ..
480 COMMON / infoc / infot, nunit, ok, lerr
481 COMMON / srnamc / srnamt
482* ..
483* .. Data statements ..
484 DATA ktype / 1, 2, 5*4, 5*6, 3*9, 10 /
485 DATA kmagn / 2*1, 3*1, 2, 3, 3*1, 2, 3, 1, 2, 3, 0 /
486 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
487 $ 0, 0, 0 /
488* ..
489* .. Executable Statements ..
490*
491* Check for errors
492*
493 info = 0
494*
495 badmm = .false.
496 badnn = .false.
497 mmax = 1
498 nmax = 1
499 mnmax = 1
500 minwrk = 1
501 DO 10 j = 1, nsizes
502 mmax = max( mmax, mval( j ) )
503 IF( mval( j ).LT.0 )
504 $ badmm = .true.
505 nmax = max( nmax, nval( j ) )
506 IF( nval( j ).LT.0 )
507 $ badnn = .true.
508 mnmax = max( mnmax, min( mval( j ), nval( j ) ) )
509 minwrk = max( minwrk, 3*( mval( j )+nval( j ) ),
510 $ mval( j )*( mval( j )+max( mval( j ), nval( j ),
511 $ nrhs )+1 )+nval( j )*min( nval( j ), mval( j ) ) )
512 10 CONTINUE
513*
514* Check for errors
515*
516 IF( nsizes.LT.0 ) THEN
517 info = -1
518 ELSE IF( badmm ) THEN
519 info = -2
520 ELSE IF( badnn ) THEN
521 info = -3
522 ELSE IF( ntypes.LT.0 ) THEN
523 info = -4
524 ELSE IF( nrhs.LT.0 ) THEN
525 info = -6
526 ELSE IF( lda.LT.mmax ) THEN
527 info = -11
528 ELSE IF( ldx.LT.mmax ) THEN
529 info = -17
530 ELSE IF( ldq.LT.mmax ) THEN
531 info = -21
532 ELSE IF( ldpt.LT.mnmax ) THEN
533 info = -23
534 ELSE IF( minwrk.GT.lwork ) THEN
535 info = -27
536 END IF
537*
538 IF( info.NE.0 ) THEN
539 CALL xerbla( 'CCHKBD', -info )
540 RETURN
541 END IF
542*
543* Initialize constants
544*
545 path( 1: 1 ) = 'Complex precision'
546 path( 2: 3 ) = 'BD'
547 nfail = 0
548 ntest = 0
549 unfl = slamch( 'Safe minimum' )
550 ovfl = slamch( 'Overflow' )
551 ulp = slamch( 'Precision' )
552 ulpinv = one / ulp
553 log2ui = int( log( ulpinv ) / log( two ) )
554 rtunfl = sqrt( unfl )
555 rtovfl = sqrt( ovfl )
556 infot = 0
557*
558* Loop over sizes, types
559*
560 DO 180 jsize = 1, nsizes
561 m = mval( jsize )
562 n = nval( jsize )
563 mnmin = min( m, n )
564 amninv = one / max( m, n, 1 )
565*
566 IF( nsizes.NE.1 ) THEN
567 mtypes = min( maxtyp, ntypes )
568 ELSE
569 mtypes = min( maxtyp+1, ntypes )
570 END IF
571*
572 DO 170 jtype = 1, mtypes
573 IF( .NOT.dotype( jtype ) )
574 $ GO TO 170
575*
576 DO 20 j = 1, 4
577 ioldsd( j ) = iseed( j )
578 20 CONTINUE
579*
580 DO 30 j = 1, 14
581 result( j ) = -one
582 30 CONTINUE
583*
584 uplo = ' '
585*
586* Compute "A"
587*
588* Control parameters:
589*
590* KMAGN KMODE KTYPE
591* =1 O(1) clustered 1 zero
592* =2 large clustered 2 identity
593* =3 small exponential (none)
594* =4 arithmetic diagonal, (w/ eigenvalues)
595* =5 random symmetric, w/ eigenvalues
596* =6 nonsymmetric, w/ singular values
597* =7 random diagonal
598* =8 random symmetric
599* =9 random nonsymmetric
600* =10 random bidiagonal (log. distrib.)
601*
602 IF( mtypes.GT.maxtyp )
603 $ GO TO 100
604*
605 itype = ktype( jtype )
606 imode = kmode( jtype )
607*
608* Compute norm
609*
610 GO TO ( 40, 50, 60 )kmagn( jtype )
611*
612 40 CONTINUE
613 anorm = one
614 GO TO 70
615*
616 50 CONTINUE
617 anorm = ( rtovfl*ulp )*amninv
618 GO TO 70
619*
620 60 CONTINUE
621 anorm = rtunfl*max( m, n )*ulpinv
622 GO TO 70
623*
624 70 CONTINUE
625*
626 CALL claset( 'Full', lda, n, czero, czero, a, lda )
627 iinfo = 0
628 cond = ulpinv
629*
630 bidiag = .false.
631 IF( itype.EQ.1 ) THEN
632*
633* Zero matrix
634*
635 iinfo = 0
636*
637 ELSE IF( itype.EQ.2 ) THEN
638*
639* Identity
640*
641 DO 80 jcol = 1, mnmin
642 a( jcol, jcol ) = anorm
643 80 CONTINUE
644*
645 ELSE IF( itype.EQ.4 ) THEN
646*
647* Diagonal Matrix, [Eigen]values Specified
648*
649 CALL clatms( mnmin, mnmin, 'S', iseed, 'N', rwork, imode,
650 $ cond, anorm, 0, 0, 'N', a, lda, work,
651 $ iinfo )
652*
653 ELSE IF( itype.EQ.5 ) THEN
654*
655* Symmetric, eigenvalues specified
656*
657 CALL clatms( mnmin, mnmin, 'S', iseed, 'S', rwork, imode,
658 $ cond, anorm, m, n, 'N', a, lda, work,
659 $ iinfo )
660*
661 ELSE IF( itype.EQ.6 ) THEN
662*
663* Nonsymmetric, singular values specified
664*
665 CALL clatms( m, n, 'S', iseed, 'N', rwork, imode, cond,
666 $ anorm, m, n, 'N', a, lda, work, iinfo )
667*
668 ELSE IF( itype.EQ.7 ) THEN
669*
670* Diagonal, random entries
671*
672 CALL clatmr( mnmin, mnmin, 'S', iseed, 'N', work, 6, one,
673 $ cone, 'T', 'N', work( mnmin+1 ), 1, one,
674 $ work( 2*mnmin+1 ), 1, one, 'N', iwork, 0, 0,
675 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
676*
677 ELSE IF( itype.EQ.8 ) THEN
678*
679* Symmetric, random entries
680*
681 CALL clatmr( mnmin, mnmin, 'S', iseed, 'S', work, 6, one,
682 $ cone, 'T', 'N', work( mnmin+1 ), 1, one,
683 $ work( m+mnmin+1 ), 1, one, 'N', iwork, m, n,
684 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
685*
686 ELSE IF( itype.EQ.9 ) THEN
687*
688* Nonsymmetric, random entries
689*
690 CALL clatmr( m, n, 'S', iseed, 'N', work, 6, one, cone,
691 $ 'T', 'N', work( mnmin+1 ), 1, one,
692 $ work( m+mnmin+1 ), 1, one, 'N', iwork, m, n,
693 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
694*
695 ELSE IF( itype.EQ.10 ) THEN
696*
697* Bidiagonal, random entries
698*
699 temp1 = -two*log( ulp )
700 DO 90 j = 1, mnmin
701 bd( j ) = exp( temp1*slarnd( 2, iseed ) )
702 IF( j.LT.mnmin )
703 $ be( j ) = exp( temp1*slarnd( 2, iseed ) )
704 90 CONTINUE
705*
706 iinfo = 0
707 bidiag = .true.
708 IF( m.GE.n ) THEN
709 uplo = 'U'
710 ELSE
711 uplo = 'L'
712 END IF
713 ELSE
714 iinfo = 1
715 END IF
716*
717 IF( iinfo.EQ.0 ) THEN
718*
719* Generate Right-Hand Side
720*
721 IF( bidiag ) THEN
722 CALL clatmr( mnmin, nrhs, 'S', iseed, 'N', work, 6,
723 $ one, cone, 'T', 'N', work( mnmin+1 ), 1,
724 $ one, work( 2*mnmin+1 ), 1, one, 'N',
725 $ iwork, mnmin, nrhs, zero, one, 'NO', y,
726 $ ldx, iwork, iinfo )
727 ELSE
728 CALL clatmr( m, nrhs, 'S', iseed, 'N', work, 6, one,
729 $ cone, 'T', 'N', work( m+1 ), 1, one,
730 $ work( 2*m+1 ), 1, one, 'N', iwork, m,
731 $ nrhs, zero, one, 'NO', x, ldx, iwork,
732 $ iinfo )
733 END IF
734 END IF
735*
736* Error Exit
737*
738 IF( iinfo.NE.0 ) THEN
739 WRITE( nout, fmt = 9998 )'Generator', iinfo, m, n,
740 $ jtype, ioldsd
741 info = abs( iinfo )
742 RETURN
743 END IF
744*
745 100 CONTINUE
746*
747* Call CGEBRD and CUNGBR to compute B, Q, and P, do tests.
748*
749 IF( .NOT.bidiag ) THEN
750*
751* Compute transformations to reduce A to bidiagonal form:
752* B := Q' * A * P.
753*
754 CALL clacpy( ' ', m, n, a, lda, q, ldq )
755 CALL cgebrd( m, n, q, ldq, bd, be, work, work( mnmin+1 ),
756 $ work( 2*mnmin+1 ), lwork-2*mnmin, iinfo )
757*
758* Check error code from CGEBRD.
759*
760 IF( iinfo.NE.0 ) THEN
761 WRITE( nout, fmt = 9998 )'CGEBRD', iinfo, m, n,
762 $ jtype, ioldsd
763 info = abs( iinfo )
764 RETURN
765 END IF
766*
767 CALL clacpy( ' ', m, n, q, ldq, pt, ldpt )
768 IF( m.GE.n ) THEN
769 uplo = 'U'
770 ELSE
771 uplo = 'L'
772 END IF
773*
774* Generate Q
775*
776 mq = m
777 IF( nrhs.LE.0 )
778 $ mq = mnmin
779 CALL cungbr( 'Q', m, mq, n, q, ldq, work,
780 $ work( 2*mnmin+1 ), lwork-2*mnmin, iinfo )
781*
782* Check error code from CUNGBR.
783*
784 IF( iinfo.NE.0 ) THEN
785 WRITE( nout, fmt = 9998 )'CUNGBR(Q)', iinfo, m, n,
786 $ jtype, ioldsd
787 info = abs( iinfo )
788 RETURN
789 END IF
790*
791* Generate P'
792*
793 CALL cungbr( 'P', mnmin, n, m, pt, ldpt, work( mnmin+1 ),
794 $ work( 2*mnmin+1 ), lwork-2*mnmin, iinfo )
795*
796* Check error code from CUNGBR.
797*
798 IF( iinfo.NE.0 ) THEN
799 WRITE( nout, fmt = 9998 )'CUNGBR(P)', iinfo, m, n,
800 $ jtype, ioldsd
801 info = abs( iinfo )
802 RETURN
803 END IF
804*
805* Apply Q' to an M by NRHS matrix X: Y := Q' * X.
806*
807 CALL cgemm( 'Conjugate transpose', 'No transpose', m,
808 $ nrhs, m, cone, q, ldq, x, ldx, czero, y,
809 $ ldx )
810*
811* Test 1: Check the decomposition A := Q * B * PT
812* 2: Check the orthogonality of Q
813* 3: Check the orthogonality of PT
814*
815 CALL cbdt01( m, n, 1, a, lda, q, ldq, bd, be, pt, ldpt,
816 $ work, rwork, result( 1 ) )
817 CALL cunt01( 'Columns', m, mq, q, ldq, work, lwork,
818 $ rwork, result( 2 ) )
819 CALL cunt01( 'Rows', mnmin, n, pt, ldpt, work, lwork,
820 $ rwork, result( 3 ) )
821 END IF
822*
823* Use CBDSQR to form the SVD of the bidiagonal matrix B:
824* B := U * S1 * VT, and compute Z = U' * Y.
825*
826 CALL scopy( mnmin, bd, 1, s1, 1 )
827 IF( mnmin.GT.0 )
828 $ CALL scopy( mnmin-1, be, 1, rwork, 1 )
829 CALL clacpy( ' ', m, nrhs, y, ldx, z, ldx )
830 CALL claset( 'Full', mnmin, mnmin, czero, cone, u, ldpt )
831 CALL claset( 'Full', mnmin, mnmin, czero, cone, vt, ldpt )
832*
833 CALL cbdsqr( uplo, mnmin, mnmin, mnmin, nrhs, s1, rwork, vt,
834 $ ldpt, u, ldpt, z, ldx, rwork( mnmin+1 ),
835 $ iinfo )
836*
837* Check error code from CBDSQR.
838*
839 IF( iinfo.NE.0 ) THEN
840 WRITE( nout, fmt = 9998 )'CBDSQR(vects)', iinfo, m, n,
841 $ jtype, ioldsd
842 info = abs( iinfo )
843 IF( iinfo.LT.0 ) THEN
844 RETURN
845 ELSE
846 result( 4 ) = ulpinv
847 GO TO 150
848 END IF
849 END IF
850*
851* Use CBDSQR to compute only the singular values of the
852* bidiagonal matrix B; U, VT, and Z should not be modified.
853*
854 CALL scopy( mnmin, bd, 1, s2, 1 )
855 IF( mnmin.GT.0 )
856 $ CALL scopy( mnmin-1, be, 1, rwork, 1 )
857*
858 CALL cbdsqr( uplo, mnmin, 0, 0, 0, s2, rwork, vt, ldpt, u,
859 $ ldpt, z, ldx, rwork( mnmin+1 ), iinfo )
860*
861* Check error code from CBDSQR.
862*
863 IF( iinfo.NE.0 ) THEN
864 WRITE( nout, fmt = 9998 )'CBDSQR(values)', iinfo, m, n,
865 $ jtype, ioldsd
866 info = abs( iinfo )
867 IF( iinfo.LT.0 ) THEN
868 RETURN
869 ELSE
870 result( 9 ) = ulpinv
871 GO TO 150
872 END IF
873 END IF
874*
875* Test 4: Check the decomposition B := U * S1 * VT
876* 5: Check the computation Z := U' * Y
877* 6: Check the orthogonality of U
878* 7: Check the orthogonality of VT
879*
880 CALL cbdt03( uplo, mnmin, 1, bd, be, u, ldpt, s1, vt, ldpt,
881 $ work, result( 4 ) )
882 CALL cbdt02( mnmin, nrhs, y, ldx, z, ldx, u, ldpt, work,
883 $ rwork, result( 5 ) )
884 CALL cunt01( 'Columns', mnmin, mnmin, u, ldpt, work, lwork,
885 $ rwork, result( 6 ) )
886 CALL cunt01( 'Rows', mnmin, mnmin, vt, ldpt, work, lwork,
887 $ rwork, result( 7 ) )
888*
889* Test 8: Check that the singular values are sorted in
890* non-increasing order and are non-negative
891*
892 result( 8 ) = zero
893 DO 110 i = 1, mnmin - 1
894 IF( s1( i ).LT.s1( i+1 ) )
895 $ result( 8 ) = ulpinv
896 IF( s1( i ).LT.zero )
897 $ result( 8 ) = ulpinv
898 110 CONTINUE
899 IF( mnmin.GE.1 ) THEN
900 IF( s1( mnmin ).LT.zero )
901 $ result( 8 ) = ulpinv
902 END IF
903*
904* Test 9: Compare CBDSQR with and without singular vectors
905*
906 temp2 = zero
907*
908 DO 120 j = 1, mnmin
909 temp1 = abs( s1( j )-s2( j ) ) /
910 $ max( sqrt( unfl )*max( s1( 1 ), one ),
911 $ ulp*max( abs( s1( j ) ), abs( s2( j ) ) ) )
912 temp2 = max( temp1, temp2 )
913 120 CONTINUE
914*
915 result( 9 ) = temp2
916*
917* Test 10: Sturm sequence test of singular values
918* Go up by factors of two until it succeeds
919*
920 temp1 = thresh*( half-ulp )
921*
922 DO 130 j = 0, log2ui
923 CALL ssvdch( mnmin, bd, be, s1, temp1, iinfo )
924 IF( iinfo.EQ.0 )
925 $ GO TO 140
926 temp1 = temp1*two
927 130 CONTINUE
928*
929 140 CONTINUE
930 result( 10 ) = temp1
931*
932* Use CBDSQR to form the decomposition A := (QU) S (VT PT)
933* from the bidiagonal form A := Q B PT.
934*
935 IF( .NOT.bidiag ) THEN
936 CALL scopy( mnmin, bd, 1, s2, 1 )
937 IF( mnmin.GT.0 )
938 $ CALL scopy( mnmin-1, be, 1, rwork, 1 )
939*
940 CALL cbdsqr( uplo, mnmin, n, m, nrhs, s2, rwork, pt,
941 $ ldpt, q, ldq, y, ldx, rwork( mnmin+1 ),
942 $ iinfo )
943*
944* Test 11: Check the decomposition A := Q*U * S2 * VT*PT
945* 12: Check the computation Z := U' * Q' * X
946* 13: Check the orthogonality of Q*U
947* 14: Check the orthogonality of VT*PT
948*
949 CALL cbdt01( m, n, 0, a, lda, q, ldq, s2, dumma, pt,
950 $ ldpt, work, rwork, result( 11 ) )
951 CALL cbdt02( m, nrhs, x, ldx, y, ldx, q, ldq, work,
952 $ rwork, result( 12 ) )
953 CALL cunt01( 'Columns', m, mq, q, ldq, work, lwork,
954 $ rwork, result( 13 ) )
955 CALL cunt01( 'Rows', mnmin, n, pt, ldpt, work, lwork,
956 $ rwork, result( 14 ) )
957 END IF
958*
959* End of Loop -- Check for RESULT(j) > THRESH
960*
961 150 CONTINUE
962 DO 160 j = 1, 14
963 IF( result( j ).GE.thresh ) THEN
964 IF( nfail.EQ.0 )
965 $ CALL slahd2( nout, path )
966 WRITE( nout, fmt = 9999 )m, n, jtype, ioldsd, j,
967 $ result( j )
968 nfail = nfail + 1
969 END IF
970 160 CONTINUE
971 IF( .NOT.bidiag ) THEN
972 ntest = ntest + 14
973 ELSE
974 ntest = ntest + 5
975 END IF
976*
977 170 CONTINUE
978 180 CONTINUE
979*
980* Summary
981*
982 CALL alasum( path, nout, nfail, ntest, 0 )
983*
984 RETURN
985*
986* End of CCHKBD
987*
988 9999 FORMAT( ' M=', i5, ', N=', i5, ', type ', i2, ', seed=',
989 $ 4( i4, ',' ), ' test(', i2, ')=', g11.4 )
990 9998 FORMAT( ' CCHKBD: ', a, ' returned INFO=', i6, '.', / 9x, 'M=',
991 $ i6, ', N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ),
992 $ i5, ')' )
993*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine cbdt01(m, n, kd, a, lda, q, ldq, d, e, pt, ldpt, work, rwork, resid)
CBDT01
Definition cbdt01.f:147
subroutine cbdt02(m, n, b, ldb, c, ldc, u, ldu, work, rwork, resid)
CBDT02
Definition cbdt02.f:120
subroutine cbdt03(uplo, n, kd, d, e, u, ldu, s, vt, ldvt, work, resid)
CBDT03
Definition cbdt03.f:135
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine clatmr(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)
CLATMR
Definition clatmr.f:490
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine cunt01(rowcol, m, n, u, ldu, work, lwork, rwork, resid)
CUNT01
Definition cunt01.f:126
subroutine cbdsqr(uplo, n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c, ldc, rwork, info)
CBDSQR
Definition cbdsqr.f:233
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine cgebrd(m, n, a, lda, d, e, tauq, taup, work, lwork, info)
CGEBRD
Definition cgebrd.f:206
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:106
subroutine cungbr(vect, m, n, k, a, lda, tau, work, lwork, info)
CUNGBR
Definition cungbr.f:157
subroutine slahd2(iounit, path)
SLAHD2
Definition slahd2.f:65
real function slarnd(idist, iseed)
SLARND
Definition slarnd.f:73
subroutine ssvdch(n, s, e, svd, tol, info)
SSVDCH
Definition ssvdch.f:97
Here is the call graph for this function:
Here is the caller graph for this function: