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

◆ schkbb()

subroutine schkbb ( integer  nsizes,
integer, dimension( * )  mval,
integer, dimension( * )  nval,
integer  nwdths,
integer, dimension( * )  kk,
integer  ntypes,
logical, dimension( * )  dotype,
integer  nrhs,
integer, dimension( 4 )  iseed,
real  thresh,
integer  nounit,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( ldab, * )  ab,
integer  ldab,
real, dimension( * )  bd,
real, dimension( * )  be,
real, dimension( ldq, * )  q,
integer  ldq,
real, dimension( ldp, * )  p,
integer  ldp,
real, dimension( ldc, * )  c,
integer  ldc,
real, dimension( ldc, * )  cc,
real, dimension( * )  work,
integer  lwork,
real, dimension( * )  result,
integer  info 
)

SCHKBB

Purpose:
 SCHKBB tests the reduction of a general real rectangular band
 matrix to bidiagonal form.

 SGBBRD factors a general band matrix A as  Q B P* , where * means
 transpose, B is upper bidiagonal, and Q and P are orthogonal;
 SGBBRD can also overwrite a given matrix C with Q* C .

 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 C are generated.
 The problem dimensions are as follows
    A:          M x N
    Q:          M x M
    P:          N x N
    B:          min(M,N) x min(M,N)
    C:          M x NRHS

 For each generated matrix, 4 tests are performed:

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

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

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

 (4)   | Y - Q' C | / ( |Y| max(M,NRHS) ulp ), where Y = Q' C.

 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:

 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 )
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).
          If NSIZES is zero, SCHKBB does nothing.  NSIZES must be at
          least zero.
[in]MVAL
          MVAL is INTEGER array, dimension (NSIZES)
          The values of the matrix row dimension M.
[in]NVAL
          NVAL is INTEGER array, dimension (NSIZES)
          The values of the matrix column dimension N.
[in]NWDTHS
          NWDTHS is INTEGER
          The number of bandwidths to use.  If it is zero,
          SCHKBB does nothing.  It must be at least zero.
[in]KK
          KK is INTEGER array, dimension (NWDTHS)
          An array containing the bandwidths to be used for the band
          matrices.  The values must be at least zero.
[in]NTYPES
          NTYPES is INTEGER
          The number of elements in DOTYPE.   If it is zero, SCHKBB
          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 matrix is in A.  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 in NN a
          matrix of that size and 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" matrix C.
          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 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 SCHKBB to continue the same random number
          sequence.
[in]THRESH
          THRESH is REAL
          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.
[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.)
[in,out]A
          A is REAL array, dimension
                            (LDA, max(NN))
          Used to hold the matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  It must be at least 1
          and at least max( NN ).
[out]AB
          AB is REAL array, dimension (LDAB, max(NN))
          Used to hold A in band storage format.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of AB.  It must be at least 2 (not 1!)
          and at least max( KK )+1.
[out]BD
          BD is REAL array, dimension (max(NN))
          Used to hold the diagonal of the bidiagonal matrix computed
          by SGBBRD.
[out]BE
          BE is REAL array, dimension (max(NN))
          Used to hold the off-diagonal of the bidiagonal matrix
          computed by SGBBRD.
[out]Q
          Q is REAL array, dimension (LDQ, max(NN))
          Used to hold the orthogonal matrix Q computed by SGBBRD.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of Q.  It must be at least 1
          and at least max( NN ).
[out]P
          P is REAL array, dimension (LDP, max(NN))
          Used to hold the orthogonal matrix P computed by SGBBRD.
[in]LDP
          LDP is INTEGER
          The leading dimension of P.  It must be at least 1
          and at least max( NN ).
[out]C
          C is REAL array, dimension (LDC, max(NN))
          Used to hold the matrix C updated by SGBBRD.
[in]LDC
          LDC is INTEGER
          The leading dimension of U.  It must be at least 1
          and at least max( NN ).
[out]CC
          CC is REAL array, dimension (LDC, max(NN))
          Used to hold a copy of the matrix C.
[out]WORK
          WORK is REAL array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The number of entries in WORK.  This must be at least
          max( LDA+1, max(NN)+1 )*max(NN).
[out]RESULT
          RESULT is REAL array, dimension (4)
          The values computed by the tests described above.
          The values are currently limited to 1/ulp, to avoid
          overflow.
[out]INFO
          INFO is INTEGER
          If 0, then everything ran OK.

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

       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.
       NTESTT          The total number of tests performed so far.
       NMAX            Largest value in NN.
       NMATS           The number of matrices generated so far.
       NERRS           The number of tests which have exceeded THRESH
                       so far.
       COND, IMODE     Values to be passed to the matrix generators.
       ANORM           Norm of A; passed to matrix generators.

       OVFL, UNFL      Overflow and underflow thresholds.
       ULP, ULPINV     Finest relative precision and its inverse.
       RTOVFL, RTUNFL  Square roots of the previous 2 values.
               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 351 of file schkbb.f.

355*
356* -- LAPACK test routine (input) --
357* -- LAPACK is a software package provided by Univ. of Tennessee, --
358* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
359*
360* .. Scalar Arguments ..
361 INTEGER INFO, LDA, LDAB, LDC, LDP, LDQ, LWORK, NOUNIT,
362 $ NRHS, NSIZES, NTYPES, NWDTHS
363 REAL THRESH
364* ..
365* .. Array Arguments ..
366 LOGICAL DOTYPE( * )
367 INTEGER ISEED( 4 ), KK( * ), MVAL( * ), NVAL( * )
368 REAL A( LDA, * ), AB( LDAB, * ), BD( * ), BE( * ),
369 $ C( LDC, * ), CC( LDC, * ), P( LDP, * ),
370 $ Q( LDQ, * ), RESULT( * ), WORK( * )
371* ..
372*
373* =====================================================================
374*
375* .. Parameters ..
376 REAL ZERO, ONE
377 parameter( zero = 0.0e0, one = 1.0e0 )
378 INTEGER MAXTYP
379 parameter( maxtyp = 15 )
380* ..
381* .. Local Scalars ..
382 LOGICAL BADMM, BADNN, BADNNB
383 INTEGER I, IINFO, IMODE, ITYPE, J, JCOL, JR, JSIZE,
384 $ JTYPE, JWIDTH, K, KL, KMAX, KU, M, MMAX, MNMAX,
385 $ MNMIN, MTYPES, N, NERRS, NMATS, NMAX, NTEST,
386 $ NTESTT
387 REAL AMNINV, ANORM, COND, OVFL, RTOVFL, RTUNFL, ULP,
388 $ ULPINV, UNFL
389* ..
390* .. Local Arrays ..
391 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KMAGN( MAXTYP ),
392 $ KMODE( MAXTYP ), KTYPE( MAXTYP )
393* ..
394* .. External Functions ..
395 REAL SLAMCH
396 EXTERNAL slamch
397* ..
398* .. External Subroutines ..
399 EXTERNAL sbdt01, sbdt02, sgbbrd, slacpy, slahd2, slaset,
401* ..
402* .. Intrinsic Functions ..
403 INTRINSIC abs, max, min, real, sqrt
404* ..
405* .. Data statements ..
406 DATA ktype / 1, 2, 5*4, 5*6, 3*9 /
407 DATA kmagn / 2*1, 3*1, 2, 3, 3*1, 2, 3, 1, 2, 3 /
408 DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
409 $ 0, 0 /
410* ..
411* .. Executable Statements ..
412*
413* Check for errors
414*
415 ntestt = 0
416 info = 0
417*
418* Important constants
419*
420 badmm = .false.
421 badnn = .false.
422 mmax = 1
423 nmax = 1
424 mnmax = 1
425 DO 10 j = 1, nsizes
426 mmax = max( mmax, mval( j ) )
427 IF( mval( j ).LT.0 )
428 $ badmm = .true.
429 nmax = max( nmax, nval( j ) )
430 IF( nval( j ).LT.0 )
431 $ badnn = .true.
432 mnmax = max( mnmax, min( mval( j ), nval( j ) ) )
433 10 CONTINUE
434*
435 badnnb = .false.
436 kmax = 0
437 DO 20 j = 1, nwdths
438 kmax = max( kmax, kk( j ) )
439 IF( kk( j ).LT.0 )
440 $ badnnb = .true.
441 20 CONTINUE
442*
443* Check for errors
444*
445 IF( nsizes.LT.0 ) THEN
446 info = -1
447 ELSE IF( badmm ) THEN
448 info = -2
449 ELSE IF( badnn ) THEN
450 info = -3
451 ELSE IF( nwdths.LT.0 ) THEN
452 info = -4
453 ELSE IF( badnnb ) THEN
454 info = -5
455 ELSE IF( ntypes.LT.0 ) THEN
456 info = -6
457 ELSE IF( nrhs.LT.0 ) THEN
458 info = -8
459 ELSE IF( lda.LT.nmax ) THEN
460 info = -13
461 ELSE IF( ldab.LT.2*kmax+1 ) THEN
462 info = -15
463 ELSE IF( ldq.LT.nmax ) THEN
464 info = -19
465 ELSE IF( ldp.LT.nmax ) THEN
466 info = -21
467 ELSE IF( ldc.LT.nmax ) THEN
468 info = -23
469 ELSE IF( ( max( lda, nmax )+1 )*nmax.GT.lwork ) THEN
470 info = -26
471 END IF
472*
473 IF( info.NE.0 ) THEN
474 CALL xerbla( 'SCHKBB', -info )
475 RETURN
476 END IF
477*
478* Quick return if possible
479*
480 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 .OR. nwdths.EQ.0 )
481 $ RETURN
482*
483* More Important constants
484*
485 unfl = slamch( 'Safe minimum' )
486 ovfl = one / unfl
487 ulp = slamch( 'Epsilon' )*slamch( 'Base' )
488 ulpinv = one / ulp
489 rtunfl = sqrt( unfl )
490 rtovfl = sqrt( ovfl )
491*
492* Loop over sizes, widths, types
493*
494 nerrs = 0
495 nmats = 0
496*
497 DO 160 jsize = 1, nsizes
498 m = mval( jsize )
499 n = nval( jsize )
500 mnmin = min( m, n )
501 amninv = one / real( max( 1, m, n ) )
502*
503 DO 150 jwidth = 1, nwdths
504 k = kk( jwidth )
505 IF( k.GE.m .AND. k.GE.n )
506 $ GO TO 150
507 kl = max( 0, min( m-1, k ) )
508 ku = max( 0, min( n-1, k ) )
509*
510 IF( nsizes.NE.1 ) THEN
511 mtypes = min( maxtyp, ntypes )
512 ELSE
513 mtypes = min( maxtyp+1, ntypes )
514 END IF
515*
516 DO 140 jtype = 1, mtypes
517 IF( .NOT.dotype( jtype ) )
518 $ GO TO 140
519 nmats = nmats + 1
520 ntest = 0
521*
522 DO 30 j = 1, 4
523 ioldsd( j ) = iseed( j )
524 30 CONTINUE
525*
526* Compute "A".
527*
528* Control parameters:
529*
530* KMAGN KMODE KTYPE
531* =1 O(1) clustered 1 zero
532* =2 large clustered 2 identity
533* =3 small exponential (none)
534* =4 arithmetic diagonal, (w/ singular values)
535* =5 random log (none)
536* =6 random nonhermitian, w/ singular values
537* =7 (none)
538* =8 (none)
539* =9 random nonhermitian
540*
541 IF( mtypes.GT.maxtyp )
542 $ GO TO 90
543*
544 itype = ktype( jtype )
545 imode = kmode( jtype )
546*
547* Compute norm
548*
549 GO TO ( 40, 50, 60 )kmagn( jtype )
550*
551 40 CONTINUE
552 anorm = one
553 GO TO 70
554*
555 50 CONTINUE
556 anorm = ( rtovfl*ulp )*amninv
557 GO TO 70
558*
559 60 CONTINUE
560 anorm = rtunfl*max( m, n )*ulpinv
561 GO TO 70
562*
563 70 CONTINUE
564*
565 CALL slaset( 'Full', lda, n, zero, zero, a, lda )
566 CALL slaset( 'Full', ldab, n, zero, zero, ab, ldab )
567 iinfo = 0
568 cond = ulpinv
569*
570* Special Matrices -- Identity & Jordan block
571*
572* Zero
573*
574 IF( itype.EQ.1 ) THEN
575 iinfo = 0
576*
577 ELSE IF( itype.EQ.2 ) THEN
578*
579* Identity
580*
581 DO 80 jcol = 1, n
582 a( jcol, jcol ) = anorm
583 80 CONTINUE
584*
585 ELSE IF( itype.EQ.4 ) THEN
586*
587* Diagonal Matrix, singular values specified
588*
589 CALL slatms( m, n, 'S', iseed, 'N', work, imode, cond,
590 $ anorm, 0, 0, 'N', a, lda, work( m+1 ),
591 $ iinfo )
592*
593 ELSE IF( itype.EQ.6 ) THEN
594*
595* Nonhermitian, singular values specified
596*
597 CALL slatms( m, n, 'S', iseed, 'N', work, imode, cond,
598 $ anorm, kl, ku, 'N', a, lda, work( m+1 ),
599 $ iinfo )
600*
601 ELSE IF( itype.EQ.9 ) THEN
602*
603* Nonhermitian, random entries
604*
605 CALL slatmr( m, n, 'S', iseed, 'N', work, 6, one, one,
606 $ 'T', 'N', work( n+1 ), 1, one,
607 $ work( 2*n+1 ), 1, one, 'N', idumma, kl,
608 $ ku, zero, anorm, 'N', a, lda, idumma,
609 $ iinfo )
610*
611 ELSE
612*
613 iinfo = 1
614 END IF
615*
616* Generate Right-Hand Side
617*
618 CALL slatmr( m, nrhs, 'S', iseed, 'N', work, 6, one, one,
619 $ 'T', 'N', work( m+1 ), 1, one,
620 $ work( 2*m+1 ), 1, one, 'N', idumma, m, nrhs,
621 $ zero, one, 'NO', c, ldc, idumma, iinfo )
622*
623 IF( iinfo.NE.0 ) THEN
624 WRITE( nounit, fmt = 9999 )'Generator', iinfo, n,
625 $ jtype, ioldsd
626 info = abs( iinfo )
627 RETURN
628 END IF
629*
630 90 CONTINUE
631*
632* Copy A to band storage.
633*
634 DO 110 j = 1, n
635 DO 100 i = max( 1, j-ku ), min( m, j+kl )
636 ab( ku+1+i-j, j ) = a( i, j )
637 100 CONTINUE
638 110 CONTINUE
639*
640* Copy C
641*
642 CALL slacpy( 'Full', m, nrhs, c, ldc, cc, ldc )
643*
644* Call SGBBRD to compute B, Q and P, and to update C.
645*
646 CALL sgbbrd( 'B', m, n, nrhs, kl, ku, ab, ldab, bd, be,
647 $ q, ldq, p, ldp, cc, ldc, work, iinfo )
648*
649 IF( iinfo.NE.0 ) THEN
650 WRITE( nounit, fmt = 9999 )'SGBBRD', iinfo, n, jtype,
651 $ ioldsd
652 info = abs( iinfo )
653 IF( iinfo.LT.0 ) THEN
654 RETURN
655 ELSE
656 result( 1 ) = ulpinv
657 GO TO 120
658 END IF
659 END IF
660*
661* Test 1: Check the decomposition A := Q * B * P'
662* 2: Check the orthogonality of Q
663* 3: Check the orthogonality of P
664* 4: Check the computation of Q' * C
665*
666 CALL sbdt01( m, n, -1, a, lda, q, ldq, bd, be, p, ldp,
667 $ work, result( 1 ) )
668 CALL sort01( 'Columns', m, m, q, ldq, work, lwork,
669 $ result( 2 ) )
670 CALL sort01( 'Rows', n, n, p, ldp, work, lwork,
671 $ result( 3 ) )
672 CALL sbdt02( m, nrhs, c, ldc, cc, ldc, q, ldq, work,
673 $ result( 4 ) )
674*
675* End of Loop -- Check for RESULT(j) > THRESH
676*
677 ntest = 4
678 120 CONTINUE
679 ntestt = ntestt + ntest
680*
681* Print out tests which fail.
682*
683 DO 130 jr = 1, ntest
684 IF( result( jr ).GE.thresh ) THEN
685 IF( nerrs.EQ.0 )
686 $ CALL slahd2( nounit, 'SBB' )
687 nerrs = nerrs + 1
688 WRITE( nounit, fmt = 9998 )m, n, k, ioldsd, jtype,
689 $ jr, result( jr )
690 END IF
691 130 CONTINUE
692*
693 140 CONTINUE
694 150 CONTINUE
695 160 CONTINUE
696*
697* Summary
698*
699 CALL slasum( 'SBB', nounit, nerrs, ntestt )
700 RETURN
701*
702 9999 FORMAT( ' SCHKBB: ', a, ' returned INFO=', i5, '.', / 9x, 'M=',
703 $ i5, ' N=', i5, ' K=', i5, ', JTYPE=', i5, ', ISEED=(',
704 $ 3( i5, ',' ), i5, ')' )
705 9998 FORMAT( ' M =', i4, ' N=', i4, ', K=', i3, ', seed=',
706 $ 4( i4, ',' ), ' type ', i2, ', test(', i2, ')=', g10.3 )
707*
708* End of SCHKBB
709*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgbbrd(vect, m, n, ncc, kl, ku, ab, ldab, d, e, q, ldq, pt, ldpt, c, ldc, work, info)
SGBBRD
Definition sgbbrd.f:187
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
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
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 sbdt01(m, n, kd, a, lda, q, ldq, d, e, pt, ldpt, work, resid)
SBDT01
Definition sbdt01.f:141
subroutine sbdt02(m, n, b, ldb, c, ldc, u, ldu, work, resid)
SBDT02
Definition sbdt02.f:112
subroutine slahd2(iounit, path)
SLAHD2
Definition slahd2.f:65
subroutine slasum(type, iounit, ie, nrun)
SLASUM
Definition slasum.f:41
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 slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321
subroutine sort01(rowcol, m, n, u, ldu, work, lwork, resid)
SORT01
Definition sort01.f:116
Here is the call graph for this function:
Here is the caller graph for this function: