LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zdrvbd ( integer  NSIZES,
integer, dimension( * )  MM,
integer, dimension( * )  NN,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer, dimension( 4 )  ISEED,
double precision  THRESH,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( ldvt, * )  VT,
integer  LDVT,
complex*16, dimension( lda, * )  ASAV,
complex*16, dimension( ldu, * )  USAV,
complex*16, dimension( ldvt, * )  VTSAV,
double precision, dimension( * )  S,
double precision, dimension( * )  SSAV,
double precision, dimension( * )  E,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUNIT,
integer  INFO 
)

ZDRVBD

Purpose:
 ZDRVBD checks the singular value decomposition (SVD) driver ZGESVD
 and ZGESDD.
 ZGESVD and ZGESDD factors A = U diag(S) VT, where U and VT are
 unitary and diag(S) is diagonal with the entries of the array S on
 its diagonal. The entries of S are the singular values, nonnegative
 and stored in decreasing order.  U and VT can be optionally not
 computed, overwritten on A, or computed partially.

 A is M by N. Let MNMIN = min( M, N ). S has dimension MNMIN.
 U can be M by M or M by MNMIN. VT can be N by N or MNMIN by N.

 When ZDRVBD is called, a number of matrix "sizes" (M's and N's)
 and a number of matrix "types" are specified.  For each size (M,N)
 and each type of matrix, and for the minimal workspace as well as
 workspace adequate to permit blocking, an  M x N  matrix "A" will be
 generated and used to test the SVD routines.  For each matrix, A will
 be factored as A = U diag(S) VT and the following 12 tests computed:

 Test for ZGESVD:

 (1)   | A - U diag(S) VT | / ( |A| max(M,N) ulp )

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

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

 (4)   S contains MNMIN nonnegative values in decreasing order.
       (Return 0 if true, 1/ULP if false.)

 (5)   | U - Upartial | / ( M ulp ) where Upartial is a partially
       computed U.

 (6)   | VT - VTpartial | / ( N ulp ) where VTpartial is a partially
       computed VT.

 (7)   | S - Spartial | / ( MNMIN ulp |S| ) where Spartial is the
       vector of singular values from the partial SVD

 Test for ZGESDD:

 (1)   | A - U diag(S) VT | / ( |A| max(M,N) ulp )

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

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

 (4)   S contains MNMIN nonnegative values in decreasing order.
       (Return 0 if true, 1/ULP if false.)

 (5)   | U - Upartial | / ( M ulp ) where Upartial is a partially
       computed U.

 (6)   | VT - VTpartial | / ( N ulp ) where VTpartial is a partially
       computed VT.

 (7)   | S - Spartial | / ( MNMIN ulp |S| ) where Spartial is the
       vector of singular values from the partial SVD

 Test for ZGESVJ:

 (1)   | A - U diag(S) VT | / ( |A| max(M,N) ulp )

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

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

 (4)   S contains MNMIN nonnegative values in decreasing order.
       (Return 0 if true, 1/ULP if false.)

 Test for ZGEJSV:

 (1)   | A - U diag(S) VT | / ( |A| max(M,N) ulp )

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

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

 (4)   S contains MNMIN nonnegative values in decreasing order.
        (Return 0 if true, 1/ULP if false.)

 Test for ZGESVDX( 'V', 'V', 'A' )/ZGESVDX( 'N', 'N', 'A' )

 (1)   | A - U diag(S) VT | / ( |A| max(M,N) ulp )

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

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

 (4)   S contains MNMIN nonnegative values in decreasing order.
       (Return 0 if true, 1/ULP if false.)

 (5)   | U - Upartial | / ( M ulp ) where Upartial is a partially
       computed U.

 (6)   | VT - VTpartial | / ( N ulp ) where VTpartial is a partially
       computed VT.

 (7)   | S - Spartial | / ( MNMIN ulp |S| ) where Spartial is the
       vector of singular values from the partial SVD
      
 Test for ZGESVDX( 'V', 'V', 'I' )

 (8)   | U' A VT''' - diag(S) | / ( |A| max(M,N) ulp )

 (9)   | I - U'U | / ( M ulp )

 (10)  | I - VT VT' | / ( N ulp )
      
 Test for ZGESVDX( 'V', 'V', 'V' )

 (11)   | U' A VT''' - diag(S) | / ( |A| max(M,N) ulp )

 (12)   | I - U'U | / ( M ulp )

 (13)   | I - VT VT' | / ( N ulp )

 The "sizes" are specified by the arrays MM(1:NSIZES) and
 NN(1:NSIZES); the value of each element pair (MM(j),NN(j))
 specifies one size.  The "types" are specified by a logical array
 DOTYPE( 1:NTYPES ); if DOTYPE(j) is .TRUE., then matrix type "j"
 will be generated.
 Currently, the list of possible types is:

 (1)  The zero matrix.
 (2)  The identity matrix.
 (3)  A matrix of the form  U D V, where U and V are unitary and
      D has evenly spaced entries 1, ..., ULP with random signs
      on the diagonal.
 (4)  Same as (3), but multiplied by the underflow-threshold / ULP.
 (5)  Same as (3), but multiplied by the overflow-threshold * ULP.
Parameters
[in]NSIZES
          NSIZES is INTEGER
          The number of sizes of matrices to use.  If it is zero,
          ZDRVBD does nothing.  It must be at least zero.
[in]MM
          MM is INTEGER array, dimension (NSIZES)
          An array containing the matrix "heights" to be used.  For
          each j=1,...,NSIZES, if MM(j) is zero, then MM(j) and NN(j)
          will be ignored.  The MM(j) values must be at least zero.
[in]NN
          NN is INTEGER array, dimension (NSIZES)
          An array containing the matrix "widths" to be used.  For
          each j=1,...,NSIZES, if NN(j) is zero, then MM(j) and NN(j)
          will be ignored.  The NN(j) values must be at least zero.
[in]NTYPES
          NTYPES is INTEGER
          The number of elements in DOTYPE.   If it is zero, ZDRVBD
          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,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 random number generator uses a linear
          congruential sequence limited to small integers, and so
          should produce machine independent random numbers. The
          values of ISEED are changed on exit, and can be used in the
          next call to ZDRVBD to continue the same random number
          sequence.
[in]THRESH
          THRESH is DOUBLE PRECISION
          A test will count as "failed" if the "error", computed as
          described above, exceeds THRESH.  Note that the error
          is scaled to be O(1), so THRESH should be a reasonably
          small multiple of 1, e.g., 10 or 100.  In particular,
          it should not depend on the precision (single vs. double)
          or the size of the matrix.  It must be at least zero.
[out]A
          A is COMPLEX*16 array, dimension (LDA,max(NN))
          Used to hold the matrix whose singular values are to be
          computed.  On exit, A contains the last matrix actually
          used.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  It must be at
          least 1 and at least max( MM ).
[out]U
          U is COMPLEX*16 array, dimension (LDU,max(MM))
          Used to hold the computed matrix of right singular vectors.
          On exit, U contains the last such vectors actually computed.
[in]LDU
          LDU is INTEGER
          The leading dimension of U.  It must be at
          least 1 and at least max( MM ).
[out]VT
          VT is COMPLEX*16 array, dimension (LDVT,max(NN))
          Used to hold the computed matrix of left singular vectors.
          On exit, VT contains the last such vectors actually computed.
[in]LDVT
          LDVT is INTEGER
          The leading dimension of VT.  It must be at
          least 1 and at least max( NN ).
[out]ASAV
          ASAV is COMPLEX*16 array, dimension (LDA,max(NN))
          Used to hold a different copy of the matrix whose singular
          values are to be computed.  On exit, A contains the last
          matrix actually used.
[out]USAV
          USAV is COMPLEX*16 array, dimension (LDU,max(MM))
          Used to hold a different copy of the computed matrix of
          right singular vectors. On exit, USAV contains the last such
          vectors actually computed.
[out]VTSAV
          VTSAV is COMPLEX*16 array, dimension (LDVT,max(NN))
          Used to hold a different copy of the computed matrix of
          left singular vectors. On exit, VTSAV contains the last such
          vectors actually computed.
[out]S
          S is DOUBLE PRECISION array, dimension (max(min(MM,NN)))
          Contains the computed singular values.
[out]SSAV
          SSAV is DOUBLE PRECISION array, dimension (max(min(MM,NN)))
          Contains another copy of the computed singular values.
[out]E
          E is DOUBLE PRECISION array, dimension (max(min(MM,NN)))
          Workspace for ZGESVD.
[out]WORK
          WORK is COMPLEX*16 array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The number of entries in WORK.  This must be at least
          MAX(3*MIN(M,N)+MAX(M,N)**2,5*MIN(M,N),3*MAX(M,N)) for all
          pairs  (M,N)=(MM(j),NN(j))
[out]RWORK
          RWORK is DOUBLE PRECISION array,
                      dimension ( 5*max(max(MM,NN)) )
[out]IWORK
          IWORK is INTEGER array, dimension at least 8*min(M,N)
[in]NOUNIT
          NOUNIT 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
           -7: THRESH < 0
          -10: LDA < 1 or LDA < MMAX, where MMAX is max( MM(j) ).
          -12: LDU < 1 or LDU < MMAX.
          -14: LDVT < 1 or LDVT < NMAX, where NMAX is max( NN(j) ).
          -21: LWORK too small.
          If  ZLATMS, or ZGESVD returns an error code, the
              absolute value of it is returned.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
June 2016

Definition at line 391 of file zdrvbd.f.

391 *
392 * -- LAPACK test routine (version 3.6.1) --
393 * -- LAPACK is a software package provided by Univ. of Tennessee, --
394 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
395 * June 2016
396 *
397 * .. Scalar Arguments ..
398  INTEGER info, lda, ldu, ldvt, lwork, nounit, nsizes,
399  $ ntypes
400  DOUBLE PRECISION thresh
401 * ..
402 * .. Array Arguments ..
403  LOGICAL dotype( * )
404  INTEGER iseed( 4 ), iwork( * ), mm( * ), nn( * )
405  DOUBLE PRECISION e( * ), rwork( * ), s( * ), ssav( * )
406  COMPLEX*16 a( lda, * ), asav( lda, * ), u( ldu, * ),
407  $ usav( ldu, * ), vt( ldvt, * ),
408  $ vtsav( ldvt, * ), work( * )
409 * ..
410 *
411 * =====================================================================
412 *
413 * .. Parameters ..
414  DOUBLE PRECISION zero, one, two, half
415  parameter ( zero = 0.0d0, one = 1.0d0, two = 2.0d0,
416  $ half = 0.5d0 )
417  COMPLEX*16 czero, cone
418  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
419  $ cone = ( 1.0d+0, 0.0d+0 ) )
420  INTEGER maxtyp
421  parameter ( maxtyp = 5 )
422 * ..
423 * .. Local Scalars ..
424  LOGICAL badmm, badnn
425  CHARACTER jobq, jobu, jobvt, range
426  INTEGER i, iinfo, ijq, iju, ijvt, il, iu, itemp,
427  $ iwspc, iwtmp, j, jsize, jtype, lswork, m,
428  $ minwrk, mmax, mnmax, mnmin, mtypes, n,
429  $ nerrs, nfail, nmax, ns, nsi, nsv, ntest,
430  $ ntestf, ntestt, lrwork
431  DOUBLE PRECISION anorm, dif, div, ovfl, rtunfl, ulp, ulpinv,
432  $ unfl, vl, vu
433 * ..
434 * .. Local Arrays ..
435  CHARACTER cjob( 4 ), cjobr( 3 ), cjobv( 2 )
436  INTEGER ioldsd( 4 ), iseed2( 4 )
437  DOUBLE PRECISION result( 35 )
438 * ..
439 * .. External Functions ..
440  DOUBLE PRECISION dlamch, dlarnd
441  EXTERNAL dlamch, dlarnd
442 * ..
443 * .. External Subroutines ..
444  EXTERNAL alasvm, xerbla, zbdt01, zbdt05, zgesdd,
447 * ..
448 * .. Intrinsic Functions ..
449  INTRINSIC abs, dble, max, min
450 * ..
451 * .. Scalars in Common ..
452  CHARACTER*32 srnamt
453 * ..
454 * .. Common blocks ..
455  COMMON / srnamc / srnamt
456 * ..
457 * .. Data statements ..
458  DATA cjob / 'N', 'O', 'S', 'A' /
459  DATA cjobr / 'A', 'V', 'I' /
460  DATA cjobv / 'N', 'V' /
461 * ..
462 * .. Executable Statements ..
463 *
464 * Check for errors
465 *
466  info = 0
467 *
468 * Important constants
469 *
470  nerrs = 0
471  ntestt = 0
472  ntestf = 0
473  badmm = .false.
474  badnn = .false.
475  mmax = 1
476  nmax = 1
477  mnmax = 1
478  minwrk = 1
479  DO 10 j = 1, nsizes
480  mmax = max( mmax, mm( j ) )
481  IF( mm( j ).LT.0 )
482  $ badmm = .true.
483  nmax = max( nmax, nn( j ) )
484  IF( nn( j ).LT.0 )
485  $ badnn = .true.
486  mnmax = max( mnmax, min( mm( j ), nn( j ) ) )
487  minwrk = max( minwrk, max( 3*min( mm( j ),
488  $ nn( j ) )+max( mm( j ), nn( j ) )**2, 5*min( mm( j ),
489  $ nn( j ) ), 3*max( mm( j ), nn( j ) ) ) )
490  10 CONTINUE
491 *
492 * Check for errors
493 *
494  IF( nsizes.LT.0 ) THEN
495  info = -1
496  ELSE IF( badmm ) THEN
497  info = -2
498  ELSE IF( badnn ) THEN
499  info = -3
500  ELSE IF( ntypes.LT.0 ) THEN
501  info = -4
502  ELSE IF( lda.LT.max( 1, mmax ) ) THEN
503  info = -10
504  ELSE IF( ldu.LT.max( 1, mmax ) ) THEN
505  info = -12
506  ELSE IF( ldvt.LT.max( 1, nmax ) ) THEN
507  info = -14
508  ELSE IF( minwrk.GT.lwork ) THEN
509  info = -21
510  END IF
511 *
512  IF( info.NE.0 ) THEN
513  CALL xerbla( 'ZDRVBD', -info )
514  RETURN
515  END IF
516 *
517 * Quick return if nothing to do
518 *
519  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
520  $ RETURN
521 *
522 * More Important constants
523 *
524  unfl = dlamch( 'S' )
525  ovfl = one / unfl
526  ulp = dlamch( 'E' )
527  ulpinv = one / ulp
528  rtunfl = sqrt( unfl )
529 *
530 * Loop over sizes, types
531 *
532  nerrs = 0
533 *
534  DO 230 jsize = 1, nsizes
535  m = mm( jsize )
536  n = nn( jsize )
537  mnmin = min( m, n )
538 *
539  IF( nsizes.NE.1 ) THEN
540  mtypes = min( maxtyp, ntypes )
541  ELSE
542  mtypes = min( maxtyp+1, ntypes )
543  END IF
544 *
545  DO 220 jtype = 1, mtypes
546  IF( .NOT.dotype( jtype ) )
547  $ GO TO 220
548  ntest = 0
549 *
550  DO 20 j = 1, 4
551  ioldsd( j ) = iseed( j )
552  20 CONTINUE
553 *
554 * Compute "A"
555 *
556  IF( mtypes.GT.maxtyp )
557  $ GO TO 50
558 *
559  IF( jtype.EQ.1 ) THEN
560 *
561 * Zero matrix
562 *
563  CALL zlaset( 'Full', m, n, czero, czero, a, lda )
564  DO 30 i = 1, min( m, n )
565  s( i ) = zero
566  30 CONTINUE
567 *
568  ELSE IF( jtype.EQ.2 ) THEN
569 *
570 * Identity matrix
571 *
572  CALL zlaset( 'Full', m, n, czero, cone, a, lda )
573  DO 40 i = 1, min( m, n )
574  s( i ) = one
575  40 CONTINUE
576 *
577  ELSE
578 *
579 * (Scaled) random matrix
580 *
581  IF( jtype.EQ.3 )
582  $ anorm = one
583  IF( jtype.EQ.4 )
584  $ anorm = unfl / ulp
585  IF( jtype.EQ.5 )
586  $ anorm = ovfl*ulp
587  CALL zlatms( m, n, 'U', iseed, 'N', s, 4, dble( mnmin ),
588  $ anorm, m-1, n-1, 'N', a, lda, work, iinfo )
589  IF( iinfo.NE.0 ) THEN
590  WRITE( nounit, fmt = 9996 )'Generator', iinfo, m, n,
591  $ jtype, ioldsd
592  info = abs( iinfo )
593  RETURN
594  END IF
595  END IF
596 *
597  50 CONTINUE
598  CALL zlacpy( 'F', m, n, a, lda, asav, lda )
599 *
600 * Do for minimal and adequate (for blocking) workspace
601 *
602  DO 210 iwspc = 1, 4
603 *
604 * Test for ZGESVD
605 *
606  iwtmp = 2*min( m, n )+max( m, n )
607  lswork = iwtmp + ( iwspc-1 )*( lwork-iwtmp ) / 3
608  lswork = min( lswork, lwork )
609  lswork = max( lswork, 1 )
610  IF( iwspc.EQ.4 )
611  $ lswork = lwork
612 *
613  DO 60 j = 1, 35
614  result( j ) = -one
615  60 CONTINUE
616 *
617 * Factorize A
618 *
619  IF( iwspc.GT.1 )
620  $ CALL zlacpy( 'F', m, n, asav, lda, a, lda )
621  srnamt = 'ZGESVD'
622  CALL zgesvd( 'A', 'A', m, n, a, lda, ssav, usav, ldu,
623  $ vtsav, ldvt, work, lswork, rwork, iinfo )
624  IF( iinfo.NE.0 ) THEN
625  WRITE( nounit, fmt = 9995 )'GESVD', iinfo, m, n,
626  $ jtype, lswork, ioldsd
627  info = abs( iinfo )
628  RETURN
629  END IF
630 *
631 * Do tests 1--4
632 *
633  CALL zbdt01( m, n, 0, asav, lda, usav, ldu, ssav, e,
634  $ vtsav, ldvt, work, rwork, result( 1 ) )
635  IF( m.NE.0 .AND. n.NE.0 ) THEN
636  CALL zunt01( 'Columns', mnmin, m, usav, ldu, work,
637  $ lwork, rwork, result( 2 ) )
638  CALL zunt01( 'Rows', mnmin, n, vtsav, ldvt, work,
639  $ lwork, rwork, result( 3 ) )
640  END IF
641  result( 4 ) = 0
642  DO 70 i = 1, mnmin - 1
643  IF( ssav( i ).LT.ssav( i+1 ) )
644  $ result( 4 ) = ulpinv
645  IF( ssav( i ).LT.zero )
646  $ result( 4 ) = ulpinv
647  70 CONTINUE
648  IF( mnmin.GE.1 ) THEN
649  IF( ssav( mnmin ).LT.zero )
650  $ result( 4 ) = ulpinv
651  END IF
652 *
653 * Do partial SVDs, comparing to SSAV, USAV, and VTSAV
654 *
655  result( 5 ) = zero
656  result( 6 ) = zero
657  result( 7 ) = zero
658  DO 100 iju = 0, 3
659  DO 90 ijvt = 0, 3
660  IF( ( iju.EQ.3 .AND. ijvt.EQ.3 ) .OR.
661  $ ( iju.EQ.1 .AND. ijvt.EQ.1 ) )GO TO 90
662  jobu = cjob( iju+1 )
663  jobvt = cjob( ijvt+1 )
664  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
665  srnamt = 'ZGESVD'
666  CALL zgesvd( jobu, jobvt, m, n, a, lda, s, u, ldu,
667  $ vt, ldvt, work, lswork, rwork, iinfo )
668 *
669 * Compare U
670 *
671  dif = zero
672  IF( m.GT.0 .AND. n.GT.0 ) THEN
673  IF( iju.EQ.1 ) THEN
674  CALL zunt03( 'C', m, mnmin, m, mnmin, usav,
675  $ ldu, a, lda, work, lwork, rwork,
676  $ dif, iinfo )
677  ELSE IF( iju.EQ.2 ) THEN
678  CALL zunt03( 'C', m, mnmin, m, mnmin, usav,
679  $ ldu, u, ldu, work, lwork, rwork,
680  $ dif, iinfo )
681  ELSE IF( iju.EQ.3 ) THEN
682  CALL zunt03( 'C', m, m, m, mnmin, usav, ldu,
683  $ u, ldu, work, lwork, rwork, dif,
684  $ iinfo )
685  END IF
686  END IF
687  result( 5 ) = max( result( 5 ), dif )
688 *
689 * Compare VT
690 *
691  dif = zero
692  IF( m.GT.0 .AND. n.GT.0 ) THEN
693  IF( ijvt.EQ.1 ) THEN
694  CALL zunt03( 'R', n, mnmin, n, mnmin, vtsav,
695  $ ldvt, a, lda, work, lwork,
696  $ rwork, dif, iinfo )
697  ELSE IF( ijvt.EQ.2 ) THEN
698  CALL zunt03( 'R', n, mnmin, n, mnmin, vtsav,
699  $ ldvt, vt, ldvt, work, lwork,
700  $ rwork, dif, iinfo )
701  ELSE IF( ijvt.EQ.3 ) THEN
702  CALL zunt03( 'R', n, n, n, mnmin, vtsav,
703  $ ldvt, vt, ldvt, work, lwork,
704  $ rwork, dif, iinfo )
705  END IF
706  END IF
707  result( 6 ) = max( result( 6 ), dif )
708 *
709 * Compare S
710 *
711  dif = zero
712  div = max( dble( mnmin )*ulp*s( 1 ),
713  $ dlamch( 'Safe minimum' ) )
714  DO 80 i = 1, mnmin - 1
715  IF( ssav( i ).LT.ssav( i+1 ) )
716  $ dif = ulpinv
717  IF( ssav( i ).LT.zero )
718  $ dif = ulpinv
719  dif = max( dif, abs( ssav( i )-s( i ) ) / div )
720  80 CONTINUE
721  result( 7 ) = max( result( 7 ), dif )
722  90 CONTINUE
723  100 CONTINUE
724 *
725 * Test for ZGESDD
726 *
727  iwtmp = 2*mnmin*mnmin + 2*mnmin + max( m, n )
728  lswork = iwtmp + ( iwspc-1 )*( lwork-iwtmp ) / 3
729  lswork = min( lswork, lwork )
730  lswork = max( lswork, 1 )
731  IF( iwspc.EQ.4 )
732  $ lswork = lwork
733 *
734 * Factorize A
735 *
736  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
737  srnamt = 'ZGESDD'
738  CALL zgesdd( 'A', m, n, a, lda, ssav, usav, ldu, vtsav,
739  $ ldvt, work, lswork, rwork, iwork, iinfo )
740  IF( iinfo.NE.0 ) THEN
741  WRITE( nounit, fmt = 9995 )'GESDD', iinfo, m, n,
742  $ jtype, lswork, ioldsd
743  info = abs( iinfo )
744  RETURN
745  END IF
746 *
747 * Do tests 1--4
748 *
749  CALL zbdt01( m, n, 0, asav, lda, usav, ldu, ssav, e,
750  $ vtsav, ldvt, work, rwork, result( 8 ) )
751  IF( m.NE.0 .AND. n.NE.0 ) THEN
752  CALL zunt01( 'Columns', mnmin, m, usav, ldu, work,
753  $ lwork, rwork, result( 9 ) )
754  CALL zunt01( 'Rows', mnmin, n, vtsav, ldvt, work,
755  $ lwork, rwork, result( 10 ) )
756  END IF
757  result( 11 ) = 0
758  DO 110 i = 1, mnmin - 1
759  IF( ssav( i ).LT.ssav( i+1 ) )
760  $ result( 11 ) = ulpinv
761  IF( ssav( i ).LT.zero )
762  $ result( 11 ) = ulpinv
763  110 CONTINUE
764  IF( mnmin.GE.1 ) THEN
765  IF( ssav( mnmin ).LT.zero )
766  $ result( 11 ) = ulpinv
767  END IF
768 *
769 * Do partial SVDs, comparing to SSAV, USAV, and VTSAV
770 *
771  result( 12 ) = zero
772  result( 13 ) = zero
773  result( 14 ) = zero
774  DO 130 ijq = 0, 2
775  jobq = cjob( ijq+1 )
776  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
777  srnamt = 'ZGESDD'
778  CALL zgesdd( jobq, m, n, a, lda, s, u, ldu, vt, ldvt,
779  $ work, lswork, rwork, iwork, iinfo )
780 *
781 * Compare U
782 *
783  dif = zero
784  IF( m.GT.0 .AND. n.GT.0 ) THEN
785  IF( ijq.EQ.1 ) THEN
786  IF( m.GE.n ) THEN
787  CALL zunt03( 'C', m, mnmin, m, mnmin, usav,
788  $ ldu, a, lda, work, lwork, rwork,
789  $ dif, iinfo )
790  ELSE
791  CALL zunt03( 'C', m, mnmin, m, mnmin, usav,
792  $ ldu, u, ldu, work, lwork, rwork,
793  $ dif, iinfo )
794  END IF
795  ELSE IF( ijq.EQ.2 ) THEN
796  CALL zunt03( 'C', m, mnmin, m, mnmin, usav, ldu,
797  $ u, ldu, work, lwork, rwork, dif,
798  $ iinfo )
799  END IF
800  END IF
801  result( 12 ) = max( result( 12 ), dif )
802 *
803 * Compare VT
804 *
805  dif = zero
806  IF( m.GT.0 .AND. n.GT.0 ) THEN
807  IF( ijq.EQ.1 ) THEN
808  IF( m.GE.n ) THEN
809  CALL zunt03( 'R', n, mnmin, n, mnmin, vtsav,
810  $ ldvt, vt, ldvt, work, lwork,
811  $ rwork, dif, iinfo )
812  ELSE
813  CALL zunt03( 'R', n, mnmin, n, mnmin, vtsav,
814  $ ldvt, a, lda, work, lwork,
815  $ rwork, dif, iinfo )
816  END IF
817  ELSE IF( ijq.EQ.2 ) THEN
818  CALL zunt03( 'R', n, mnmin, n, mnmin, vtsav,
819  $ ldvt, vt, ldvt, work, lwork, rwork,
820  $ dif, iinfo )
821  END IF
822  END IF
823  result( 13 ) = max( result( 13 ), dif )
824 *
825 * Compare S
826 *
827  dif = zero
828  div = max( dble( mnmin )*ulp*s( 1 ),
829  $ dlamch( 'Safe minimum' ) )
830  DO 120 i = 1, mnmin - 1
831  IF( ssav( i ).LT.ssav( i+1 ) )
832  $ dif = ulpinv
833  IF( ssav( i ).LT.zero )
834  $ dif = ulpinv
835  dif = max( dif, abs( ssav( i )-s( i ) ) / div )
836  120 CONTINUE
837  result( 14 ) = max( result( 14 ), dif )
838  130 CONTINUE
839 
840 *
841 * Test ZGESVJ: Factorize A
842 * Note: ZGESVJ does not work for M < N
843 *
844  result( 15 ) = zero
845  result( 16 ) = zero
846  result( 17 ) = zero
847  result( 18 ) = zero
848 *
849  IF( m.GE.n ) THEN
850  iwtmp = 2*mnmin*mnmin + 2*mnmin + max( m, n )
851  lswork = iwtmp + ( iwspc-1 )*( lwork-iwtmp ) / 3
852  lswork = min( lswork, lwork )
853  lswork = max( lswork, 1 )
854  lrwork = max(6,n)
855  IF( iwspc.EQ.4 )
856  $ lswork = lwork
857 *
858  CALL zlacpy( 'F', m, n, asav, lda, usav, lda )
859  srnamt = 'ZGESVJ'
860  CALL zgesvj( 'G', 'U', 'V', m, n, usav, lda, ssav,
861  & 0, a, ldvt, work, lwork, rwork,
862  & lrwork, iinfo )
863 *
864 * ZGESVJ retuns V not VT, so we transpose to use the same
865 * test suite.
866 *
867  DO j=1,n
868  DO i=1,n
869  vtsav(j,i) = conjg(a(i,j))
870  END DO
871  END DO
872 *
873  IF( iinfo.NE.0 ) THEN
874  WRITE( nounit, fmt = 9995 )'GESVJ', iinfo, m, n,
875  $ jtype, lswork, ioldsd
876  info = abs( iinfo )
877  RETURN
878  END IF
879 *
880 * Do tests 15--18
881 *
882  CALL zbdt01( m, n, 0, asav, lda, usav, ldu, ssav, e,
883  $ vtsav, ldvt, work, rwork, result( 15 ) )
884  IF( m.NE.0 .AND. n.NE.0 ) THEN
885  CALL zunt01( 'Columns', m, m, usav, ldu, work,
886  $ lwork, rwork, result( 16 ) )
887  CALL zunt01( 'Rows', n, n, vtsav, ldvt, work,
888  $ lwork, rwork, result( 17 ) )
889  END IF
890  result( 18 ) = zero
891  DO 131 i = 1, mnmin - 1
892  IF( ssav( i ).LT.ssav( i+1 ) )
893  $ result( 18 ) = ulpinv
894  IF( ssav( i ).LT.zero )
895  $ result( 18 ) = ulpinv
896  131 CONTINUE
897  IF( mnmin.GE.1 ) THEN
898  IF( ssav( mnmin ).LT.zero )
899  $ result( 18 ) = ulpinv
900  END IF
901  END IF
902 *
903 * Test ZGEJSV: Factorize A
904 * Note: ZGEJSV does not work for M < N
905 *
906  result( 19 ) = zero
907  result( 20 ) = zero
908  result( 21 ) = zero
909  result( 22 ) = zero
910  IF( m.GE.n ) THEN
911  iwtmp = 2*mnmin*mnmin + 2*mnmin + max( m, n )
912  lswork = iwtmp + ( iwspc-1 )*( lwork-iwtmp ) / 3
913  lswork = min( lswork, lwork )
914  lswork = max( lswork, 1 )
915  IF( iwspc.EQ.4 )
916  $ lswork = lwork
917  lrwork = max( 7, n + 2*m)
918 *
919  CALL zlacpy( 'F', m, n, asav, lda, vtsav, lda )
920  srnamt = 'ZGEJSV'
921  CALL zgejsv( 'G', 'U', 'V', 'R', 'N', 'N',
922  & m, n, vtsav, lda, ssav, usav, ldu, a, ldvt,
923  & work, lwork, rwork,
924  & lrwork, iwork, iinfo )
925 *
926 * ZGEJSV retuns V not VT, so we transpose to use the same
927 * test suite.
928 *
929  DO 133 j=1,n
930  DO 132 i=1,n
931  vtsav(j,i) = conjg(a(i,j))
932  132 END DO
933  133 END DO
934 *
935  IF( iinfo.NE.0 ) THEN
936  WRITE( nounit, fmt = 9995 )'GESVJ', iinfo, m, n,
937  $ jtype, lswork, ioldsd
938  info = abs( iinfo )
939  RETURN
940  END IF
941 *
942 * Do tests 19--22
943 *
944  CALL zbdt01( m, n, 0, asav, lda, usav, ldu, ssav, e,
945  $ vtsav, ldvt, work, rwork, result( 19 ) )
946  IF( m.NE.0 .AND. n.NE.0 ) THEN
947  CALL zunt01( 'Columns', m, m, usav, ldu, work,
948  $ lwork, rwork, result( 20 ) )
949  CALL zunt01( 'Rows', n, n, vtsav, ldvt, work,
950  $ lwork, rwork, result( 21 ) )
951  END IF
952  result( 22 ) = zero
953  DO 134 i = 1, mnmin - 1
954  IF( ssav( i ).LT.ssav( i+1 ) )
955  $ result( 22 ) = ulpinv
956  IF( ssav( i ).LT.zero )
957  $ result( 22 ) = ulpinv
958  134 CONTINUE
959  IF( mnmin.GE.1 ) THEN
960  IF( ssav( mnmin ).LT.zero )
961  $ result( 22 ) = ulpinv
962  END IF
963  END IF
964 *
965 * Test ZGESVDX
966 *
967 * Factorize A
968 *
969  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
970  srnamt = 'ZGESVDX'
971  CALL zgesvdx( 'V', 'V', 'A', m, n, a, lda,
972  $ vl, vu, il, iu, ns, ssav, usav, ldu,
973  $ vtsav, ldvt, work, lwork, rwork,
974  $ iwork, iinfo )
975  IF( iinfo.NE.0 ) THEN
976  WRITE( nounit, fmt = 9995 )'GESVDX', iinfo, m, n,
977  $ jtype, lswork, ioldsd
978  info = abs( iinfo )
979  RETURN
980  END IF
981 *
982 * Do tests 1--4
983 *
984  result( 23 ) = zero
985  result( 24 ) = zero
986  result( 25 ) = zero
987  CALL zbdt01( m, n, 0, asav, lda, usav, ldu, ssav, e,
988  $ vtsav, ldvt, work, rwork, result( 23 ) )
989  IF( m.NE.0 .AND. n.NE.0 ) THEN
990  CALL zunt01( 'Columns', mnmin, m, usav, ldu, work,
991  $ lwork, rwork, result( 24 ) )
992  CALL zunt01( 'Rows', mnmin, n, vtsav, ldvt, work,
993  $ lwork, rwork, result( 25 ) )
994  END IF
995  result( 26 ) = zero
996  DO 140 i = 1, mnmin - 1
997  IF( ssav( i ).LT.ssav( i+1 ) )
998  $ result( 26 ) = ulpinv
999  IF( ssav( i ).LT.zero )
1000  $ result( 26 ) = ulpinv
1001  140 CONTINUE
1002  IF( mnmin.GE.1 ) THEN
1003  IF( ssav( mnmin ).LT.zero )
1004  $ result( 26 ) = ulpinv
1005  END IF
1006 *
1007 * Do partial SVDs, comparing to SSAV, USAV, and VTSAV
1008 *
1009  result( 27 ) = zero
1010  result( 28 ) = zero
1011  result( 29 ) = zero
1012  DO 170 iju = 0, 1
1013  DO 160 ijvt = 0, 1
1014  IF( ( iju.EQ.0 .AND. ijvt.EQ.0 ) .OR.
1015  $ ( iju.EQ.1 .AND. ijvt.EQ.1 ) ) GO TO 160
1016  jobu = cjobv( iju+1 )
1017  jobvt = cjobv( ijvt+1 )
1018  range = cjobr( 1 )
1019  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
1020  srnamt = 'ZGESVDX'
1021  CALL zgesvdx( jobu, jobvt, 'A', m, n, a, lda,
1022  $ vl, vu, il, iu, ns, ssav, u, ldu,
1023  $ vt, ldvt, work, lwork, rwork,
1024  $ iwork, iinfo )
1025 *
1026 * Compare U
1027 *
1028  dif = zero
1029  IF( m.GT.0 .AND. n.GT.0 ) THEN
1030  IF( iju.EQ.1 ) THEN
1031  CALL zunt03( 'C', m, mnmin, m, mnmin, usav,
1032  $ ldu, u, ldu, work, lwork, rwork,
1033  $ dif, iinfo )
1034  END IF
1035  END IF
1036  result( 27 ) = max( result( 27 ), dif )
1037 *
1038 * Compare VT
1039 *
1040  dif = zero
1041  IF( m.GT.0 .AND. n.GT.0 ) THEN
1042  IF( ijvt.EQ.1 ) THEN
1043  CALL zunt03( 'R', n, mnmin, n, mnmin, vtsav,
1044  $ ldvt, vt, ldvt, work, lwork,
1045  $ rwork, dif, iinfo )
1046  END IF
1047  END IF
1048  result( 28 ) = max( result( 28 ), dif )
1049 *
1050 * Compare S
1051 *
1052  dif = zero
1053  div = max( dble( mnmin )*ulp*s( 1 ),
1054  $ dlamch( 'Safe minimum' ) )
1055  DO 150 i = 1, mnmin - 1
1056  IF( ssav( i ).LT.ssav( i+1 ) )
1057  $ dif = ulpinv
1058  IF( ssav( i ).LT.zero )
1059  $ dif = ulpinv
1060  dif = max( dif, abs( ssav( i )-s( i ) ) / div )
1061  150 CONTINUE
1062  result( 29) = max( result( 29 ), dif )
1063  160 CONTINUE
1064  170 CONTINUE
1065 *
1066 * Do tests 8--10
1067 *
1068  DO 180 i = 1, 4
1069  iseed2( i ) = iseed( i )
1070  180 CONTINUE
1071  IF( mnmin.LE.1 ) THEN
1072  il = 1
1073  iu = max( 1, mnmin )
1074  ELSE
1075  il = 1 + int( ( mnmin-1 )*dlarnd( 1, iseed2 ) )
1076  iu = 1 + int( ( mnmin-1 )*dlarnd( 1, iseed2 ) )
1077  IF( iu.LT.il ) THEN
1078  itemp = iu
1079  iu = il
1080  il = itemp
1081  END IF
1082  END IF
1083  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
1084  srnamt = 'ZGESVDX'
1085  CALL zgesvdx( 'V', 'V', 'I', m, n, a, lda,
1086  $ vl, vu, il, iu, nsi, s, u, ldu,
1087  $ vt, ldvt, work, lwork, rwork,
1088  $ iwork, iinfo )
1089  IF( iinfo.NE.0 ) THEN
1090  WRITE( nounit, fmt = 9995 )'GESVDX', iinfo, m, n,
1091  $ jtype, lswork, ioldsd
1092  info = abs( iinfo )
1093  RETURN
1094  END IF
1095 *
1096  result( 30 ) = zero
1097  result( 31 ) = zero
1098  result( 32 ) = zero
1099  CALL zbdt05( m, n, asav, lda, s, nsi, u, ldu,
1100  $ vt, ldvt, work, result( 30 ) )
1101  IF( m.NE.0 .AND. n.NE.0 ) THEN
1102  CALL zunt01( 'Columns', m, nsi, u, ldu, work,
1103  $ lwork, rwork, result( 31 ) )
1104  CALL zunt01( 'Rows', nsi, n, vt, ldvt, work,
1105  $ lwork, rwork, result( 32 ) )
1106  END IF
1107 *
1108 * Do tests 11--13
1109 *
1110  IF( mnmin.GT.0 .AND. nsi.GT.1 ) THEN
1111  IF( il.NE.1 ) THEN
1112  vu = ssav( il ) +
1113  $ max( half*abs( ssav( il )-ssav( il-1 ) ),
1114  $ ulp*anorm, two*rtunfl )
1115  ELSE
1116  vu = ssav( 1 ) +
1117  $ max( half*abs( ssav( ns )-ssav( 1 ) ),
1118  $ ulp*anorm, two*rtunfl )
1119  END IF
1120  IF( iu.NE.ns ) THEN
1121  vl = ssav( iu ) - max( ulp*anorm, two*rtunfl,
1122  $ half*abs( ssav( iu+1 )-ssav( iu ) ) )
1123  ELSE
1124  vl = ssav( ns ) - max( ulp*anorm, two*rtunfl,
1125  $ half*abs( ssav( ns )-ssav( 1 ) ) )
1126  END IF
1127  vl = max( vl,zero )
1128  vu = max( vu,zero )
1129  IF( vl.GE.vu ) vu = max( vu*2, vu+vl+half )
1130  ELSE
1131  vl = zero
1132  vu = one
1133  END IF
1134  CALL zlacpy( 'F', m, n, asav, lda, a, lda )
1135  srnamt = 'ZGESVDX'
1136  CALL zgesvdx( 'V', 'V', 'V', m, n, a, lda,
1137  $ vl, vu, il, iu, nsv, s, u, ldu,
1138  $ vt, ldvt, work, lwork, rwork,
1139  $ iwork, iinfo )
1140  IF( iinfo.NE.0 ) THEN
1141  WRITE( nounit, fmt = 9995 )'GESVDX', iinfo, m, n,
1142  $ jtype, lswork, ioldsd
1143  info = abs( iinfo )
1144  RETURN
1145  END IF
1146 *
1147  result( 33 ) = zero
1148  result( 34 ) = zero
1149  result( 35 ) = zero
1150  CALL zbdt05( m, n, asav, lda, s, nsv, u, ldu,
1151  $ vt, ldvt, work, result( 33 ) )
1152  IF( m.NE.0 .AND. n.NE.0 ) THEN
1153  CALL zunt01( 'Columns', m, nsv, u, ldu, work,
1154  $ lwork, rwork, result( 34 ) )
1155  CALL zunt01( 'Rows', nsv, n, vt, ldvt, work,
1156  $ lwork, rwork, result( 35 ) )
1157  END IF
1158 *
1159 * End of Loop -- Check for RESULT(j) > THRESH
1160 *
1161  ntest = 0
1162  nfail = 0
1163  DO 190 j = 1, 35
1164  IF( result( j ).GE.zero )
1165  $ ntest = ntest + 1
1166  IF( result( j ).GE.thresh )
1167  $ nfail = nfail + 1
1168  190 CONTINUE
1169 *
1170  IF( nfail.GT.0 )
1171  $ ntestf = ntestf + 1
1172  IF( ntestf.EQ.1 ) THEN
1173  WRITE( nounit, fmt = 9999 )
1174  WRITE( nounit, fmt = 9998 )thresh
1175  ntestf = 2
1176  END IF
1177 *
1178  DO 200 j = 1, 35
1179  IF( result( j ).GE.thresh ) THEN
1180  WRITE( nounit, fmt = 9997 )m, n, jtype, iwspc,
1181  $ ioldsd, j, result( j )
1182  END IF
1183  200 CONTINUE
1184 *
1185  nerrs = nerrs + nfail
1186  ntestt = ntestt + ntest
1187 *
1188  210 CONTINUE
1189 *
1190  220 CONTINUE
1191  230 CONTINUE
1192 *
1193 * Summary
1194 *
1195  CALL alasvm( 'ZBD', nounit, nerrs, ntestt, 0 )
1196 *
1197  9999 FORMAT( ' SVD -- Complex Singular Value Decomposition Driver ',
1198  $ / ' Matrix types (see ZDRVBD for details):',
1199  $ / / ' 1 = Zero matrix', / ' 2 = Identity matrix',
1200  $ / ' 3 = Evenly spaced singular values near 1',
1201  $ / ' 4 = Evenly spaced singular values near underflow',
1202  $ / ' 5 = Evenly spaced singular values near overflow',
1203  $ / / ' Tests performed: ( A is dense, U and V are unitary,',
1204  $ / 19x, ' S is an array, and Upartial, VTpartial, and',
1205  $ / 19x, ' Spartial are partially computed U, VT and S),', / )
1206  9998 FORMAT( ' Tests performed with Test Threshold = ', f8.2,
1207  $ / ' ZGESVD: ', /
1208  $ ' 1 = | A - U diag(S) VT | / ( |A| max(M,N) ulp ) ',
1209  $ / ' 2 = | I - U**T U | / ( M ulp ) ',
1210  $ / ' 3 = | I - VT VT**T | / ( N ulp ) ',
1211  $ / ' 4 = 0 if S contains min(M,N) nonnegative values in',
1212  $ ' decreasing order, else 1/ulp',
1213  $ / ' 5 = | U - Upartial | / ( M ulp )',
1214  $ / ' 6 = | VT - VTpartial | / ( N ulp )',
1215  $ / ' 7 = | S - Spartial | / ( min(M,N) ulp |S| )',
1216  $ / ' ZGESDD: ', /
1217  $ ' 8 = | A - U diag(S) VT | / ( |A| max(M,N) ulp ) ',
1218  $ / ' 9 = | I - U**T U | / ( M ulp ) ',
1219  $ / '10 = | I - VT VT**T | / ( N ulp ) ',
1220  $ / '11 = 0 if S contains min(M,N) nonnegative values in',
1221  $ ' decreasing order, else 1/ulp',
1222  $ / '12 = | U - Upartial | / ( M ulp )',
1223  $ / '13 = | VT - VTpartial | / ( N ulp )',
1224  $ / '14 = | S - Spartial | / ( min(M,N) ulp |S| )',
1225  $ / ' ZGESVJ: ', /
1226  $ / '15 = | A - U diag(S) VT | / ( |A| max(M,N) ulp ) ',
1227  $ / '16 = | I - U**T U | / ( M ulp ) ',
1228  $ / '17 = | I - VT VT**T | / ( N ulp ) ',
1229  $ / '18 = 0 if S contains min(M,N) nonnegative values in',
1230  $ ' decreasing order, else 1/ulp',
1231  $ / ' ZGESJV: ', /
1232  $ / '19 = | A - U diag(S) VT | / ( |A| max(M,N) ulp )',
1233  $ / '20 = | I - U**T U | / ( M ulp ) ',
1234  $ / '21 = | I - VT VT**T | / ( N ulp ) ',
1235  $ / '22 = 0 if S contains min(M,N) nonnegative values in',
1236  $ ' decreasing order, else 1/ulp',
1237  $ / ' ZGESVDX(V,V,A): ', /
1238  $ '23 = | A - U diag(S) VT | / ( |A| max(M,N) ulp ) ',
1239  $ / '24 = | I - U**T U | / ( M ulp ) ',
1240  $ / '25 = | I - VT VT**T | / ( N ulp ) ',
1241  $ / '26 = 0 if S contains min(M,N) nonnegative values in',
1242  $ ' decreasing order, else 1/ulp',
1243  $ / '27 = | U - Upartial | / ( M ulp )',
1244  $ / '28 = | VT - VTpartial | / ( N ulp )',
1245  $ / '29 = | S - Spartial | / ( min(M,N) ulp |S| )',
1246  $ / ' ZGESVDX(V,V,I): ',
1247  $ / '30 = | U**T A VT**T - diag(S) | / ( |A| max(M,N) ulp )',
1248  $ / '31 = | I - U**T U | / ( M ulp ) ',
1249  $ / '32 = | I - VT VT**T | / ( N ulp ) ',
1250  $ / ' ZGESVDX(V,V,V) ',
1251  $ / '33 = | U**T A VT**T - diag(S) | / ( |A| max(M,N) ulp )',
1252  $ / '34 = | I - U**T U | / ( M ulp ) ',
1253  $ / '35 = | I - VT VT**T | / ( N ulp ) ',
1254  $ / / )
1255  9997 FORMAT( ' M=', i5, ', N=', i5, ', type ', i1, ', IWS=', i1,
1256  $ ', seed=', 4( i4, ',' ), ' test(', i2, ')=', g11.4 )
1257  9996 FORMAT( ' ZDRVBD: ', a, ' returned INFO=', i6, '.', / 9x, 'M=',
1258  $ i6, ', N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ),
1259  $ i5, ')' )
1260  9995 FORMAT( ' ZDRVBD: ', a, ' returned INFO=', i6, '.', / 9x, 'M=',
1261  $ i6, ', N=', i6, ', JTYPE=', i6, ', LSWORK=', i6, / 9x,
1262  $ 'ISEED=(', 3( i5, ',' ), i5, ')' )
1263 *
1264  RETURN
1265 *
1266 * End of ZDRVBD
1267 *
subroutine zgesdd(JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, IWORK, INFO)
ZGESDD
Definition: zgesdd.f:229
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zgesvd(JOBU, JOBVT, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, INFO)
ZGESVD computes the singular value decomposition (SVD) for GE matrices
Definition: zgesvd.f:216
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zgesvdx(JOBU, JOBVT, RANGE, M, N, A, LDA, VL, VU, IL, IU, NS, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, IWORK, INFO)
ZGESVDX computes the singular value decomposition (SVD) for GE matrices
Definition: zgesvdx.f:272
subroutine zunt01(ROWCOL, M, N, U, LDU, WORK, LWORK, RWORK, RESID)
ZUNT01
Definition: zunt01.f:128
subroutine zbdt01(M, N, KD, A, LDA, Q, LDQ, D, E, PT, LDPT, WORK, RWORK, RESID)
ZBDT01
Definition: zbdt01.f:148
double precision function dlarnd(IDIST, ISEED)
DLARND
Definition: dlarnd.f:75
subroutine zgesvj(JOBA, JOBU, JOBV, M, N, A, LDA, SVA, MV, V, LDV, CWORK, LWORK, RWORK, LRWORK, INFO)
ZGESVJ
Definition: zgesvj.f:344
subroutine zgejsv(JOBA, JOBU, JOBV, JOBR, JOBT, JOBP, M, N, A, LDA, SVA, U, LDU, V, LDV, CWORK, LWORK, RWORK, LRWORK, IWORK, INFO)
ZGEJSV
Definition: zgejsv.f:519
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zbdt05(M, N, A, LDA, S, NS, U, LDU, VT, LDVT, WORK, RESID)
Definition: zbdt05.f:125
subroutine zunt03(RC, MU, MV, N, K, U, LDU, V, LDV, WORK, LWORK, RWORK, RESULT, INFO)
ZUNT03
Definition: zunt03.f:164

Here is the call graph for this function:

Here is the caller graph for this function: