LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dchkbb ( integer  NSIZES,
integer, dimension( * )  MVAL,
integer, dimension( * )  NVAL,
integer  NWDTHS,
integer, dimension( * )  KK,
integer  NTYPES,
logical, dimension( * )  DOTYPE,
integer  NRHS,
integer, dimension( 4 )  ISEED,
double precision  THRESH,
integer  NOUNIT,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldab, * )  AB,
integer  LDAB,
double precision, dimension( * )  BD,
double precision, dimension( * )  BE,
double precision, dimension( ldq, * )  Q,
integer  LDQ,
double precision, dimension( ldp, * )  P,
integer  LDP,
double precision, dimension( ldc, * )  C,
integer  LDC,
double precision, dimension( ldc, * )  CC,
double precision, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RESULT,
integer  INFO 
)

DCHKBB

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

 DGBBRD factors a general band matrix A as  Q B P* , where * means
 transpose, B is upper bidiagonal, and Q and P are orthogonal;
 DGBBRD 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, DCHKBB 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,
          DCHKBB 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, DCHKBB
          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 DCHKBB 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.
[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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (max(NN))
          Used to hold the diagonal of the bidiagonal matrix computed
          by DGBBRD.
[out]BE
          BE is DOUBLE PRECISION array, dimension (max(NN))
          Used to hold the off-diagonal of the bidiagonal matrix
          computed by DGBBRD.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDQ, max(NN))
          Used to hold the orthogonal matrix Q computed by DGBBRD.
[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 DOUBLE PRECISION array, dimension (LDP, max(NN))
          Used to hold the orthogonal matrix P computed by DGBBRD.
[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 DOUBLE PRECISION array, dimension (LDC, max(NN))
          Used to hold the matrix C updated by DGBBRD.
[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 DOUBLE PRECISION array, dimension (LDC, max(NN))
          Used to hold a copy of the matrix C.
[out]WORK
          WORK is DOUBLE PRECISION 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 DOUBLE PRECISION 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.
Date
November 2011

Definition at line 357 of file dchkbb.f.

357 *
358 * -- LAPACK test routine (input) --
359 * -- LAPACK is a software package provided by Univ. of Tennessee, --
360 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
361 * November 2011
362 *
363 * .. Scalar Arguments ..
364  INTEGER info, lda, ldab, ldc, ldp, ldq, lwork, nounit,
365  $ nrhs, nsizes, ntypes, nwdths
366  DOUBLE PRECISION thresh
367 * ..
368 * .. Array Arguments ..
369  LOGICAL dotype( * )
370  INTEGER iseed( 4 ), kk( * ), mval( * ), nval( * )
371  DOUBLE PRECISION a( lda, * ), ab( ldab, * ), bd( * ), be( * ),
372  $ c( ldc, * ), cc( ldc, * ), p( ldp, * ),
373  $ q( ldq, * ), result( * ), work( * )
374 * ..
375 *
376 * =====================================================================
377 *
378 * .. Parameters ..
379  DOUBLE PRECISION zero, one
380  parameter ( zero = 0.0d0, one = 1.0d0 )
381  INTEGER maxtyp
382  parameter ( maxtyp = 15 )
383 * ..
384 * .. Local Scalars ..
385  LOGICAL badmm, badnn, badnnb
386  INTEGER i, iinfo, imode, itype, j, jcol, jr, jsize,
387  $ jtype, jwidth, k, kl, kmax, ku, m, mmax, mnmax,
388  $ mnmin, mtypes, n, nerrs, nmats, nmax, ntest,
389  $ ntestt
390  DOUBLE PRECISION amninv, anorm, cond, ovfl, rtovfl, rtunfl, ulp,
391  $ ulpinv, unfl
392 * ..
393 * .. Local Arrays ..
394  INTEGER idumma( 1 ), ioldsd( 4 ), kmagn( maxtyp ),
395  $ kmode( maxtyp ), ktype( maxtyp )
396 * ..
397 * .. External Functions ..
398  DOUBLE PRECISION dlamch
399  EXTERNAL dlamch
400 * ..
401 * .. External Subroutines ..
402  EXTERNAL dbdt01, dbdt02, dgbbrd, dlacpy, dlahd2, dlaset,
404 * ..
405 * .. Intrinsic Functions ..
406  INTRINSIC abs, dble, max, min, sqrt
407 * ..
408 * .. Data statements ..
409  DATA ktype / 1, 2, 5*4, 5*6, 3*9 /
410  DATA kmagn / 2*1, 3*1, 2, 3, 3*1, 2, 3, 1, 2, 3 /
411  DATA kmode / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
412  $ 0, 0 /
413 * ..
414 * .. Executable Statements ..
415 *
416 * Check for errors
417 *
418  ntestt = 0
419  info = 0
420 *
421 * Important constants
422 *
423  badmm = .false.
424  badnn = .false.
425  mmax = 1
426  nmax = 1
427  mnmax = 1
428  DO 10 j = 1, nsizes
429  mmax = max( mmax, mval( j ) )
430  IF( mval( j ).LT.0 )
431  $ badmm = .true.
432  nmax = max( nmax, nval( j ) )
433  IF( nval( j ).LT.0 )
434  $ badnn = .true.
435  mnmax = max( mnmax, min( mval( j ), nval( j ) ) )
436  10 CONTINUE
437 *
438  badnnb = .false.
439  kmax = 0
440  DO 20 j = 1, nwdths
441  kmax = max( kmax, kk( j ) )
442  IF( kk( j ).LT.0 )
443  $ badnnb = .true.
444  20 CONTINUE
445 *
446 * Check for errors
447 *
448  IF( nsizes.LT.0 ) THEN
449  info = -1
450  ELSE IF( badmm ) THEN
451  info = -2
452  ELSE IF( badnn ) THEN
453  info = -3
454  ELSE IF( nwdths.LT.0 ) THEN
455  info = -4
456  ELSE IF( badnnb ) THEN
457  info = -5
458  ELSE IF( ntypes.LT.0 ) THEN
459  info = -6
460  ELSE IF( nrhs.LT.0 ) THEN
461  info = -8
462  ELSE IF( lda.LT.nmax ) THEN
463  info = -13
464  ELSE IF( ldab.LT.2*kmax+1 ) THEN
465  info = -15
466  ELSE IF( ldq.LT.nmax ) THEN
467  info = -19
468  ELSE IF( ldp.LT.nmax ) THEN
469  info = -21
470  ELSE IF( ldc.LT.nmax ) THEN
471  info = -23
472  ELSE IF( ( max( lda, nmax )+1 )*nmax.GT.lwork ) THEN
473  info = -26
474  END IF
475 *
476  IF( info.NE.0 ) THEN
477  CALL xerbla( 'DCHKBB', -info )
478  RETURN
479  END IF
480 *
481 * Quick return if possible
482 *
483  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 .OR. nwdths.EQ.0 )
484  $ RETURN
485 *
486 * More Important constants
487 *
488  unfl = dlamch( 'Safe minimum' )
489  ovfl = one / unfl
490  ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
491  ulpinv = one / ulp
492  rtunfl = sqrt( unfl )
493  rtovfl = sqrt( ovfl )
494 *
495 * Loop over sizes, widths, types
496 *
497  nerrs = 0
498  nmats = 0
499 *
500  DO 160 jsize = 1, nsizes
501  m = mval( jsize )
502  n = nval( jsize )
503  mnmin = min( m, n )
504  amninv = one / dble( max( 1, m, n ) )
505 *
506  DO 150 jwidth = 1, nwdths
507  k = kk( jwidth )
508  IF( k.GE.m .AND. k.GE.n )
509  $ GO TO 150
510  kl = max( 0, min( m-1, k ) )
511  ku = max( 0, min( n-1, k ) )
512 *
513  IF( nsizes.NE.1 ) THEN
514  mtypes = min( maxtyp, ntypes )
515  ELSE
516  mtypes = min( maxtyp+1, ntypes )
517  END IF
518 *
519  DO 140 jtype = 1, mtypes
520  IF( .NOT.dotype( jtype ) )
521  $ GO TO 140
522  nmats = nmats + 1
523  ntest = 0
524 *
525  DO 30 j = 1, 4
526  ioldsd( j ) = iseed( j )
527  30 CONTINUE
528 *
529 * Compute "A".
530 *
531 * Control parameters:
532 *
533 * KMAGN KMODE KTYPE
534 * =1 O(1) clustered 1 zero
535 * =2 large clustered 2 identity
536 * =3 small exponential (none)
537 * =4 arithmetic diagonal, (w/ singular values)
538 * =5 random log (none)
539 * =6 random nonhermitian, w/ singular values
540 * =7 (none)
541 * =8 (none)
542 * =9 random nonhermitian
543 *
544  IF( mtypes.GT.maxtyp )
545  $ GO TO 90
546 *
547  itype = ktype( jtype )
548  imode = kmode( jtype )
549 *
550 * Compute norm
551 *
552  GO TO ( 40, 50, 60 )kmagn( jtype )
553 *
554  40 CONTINUE
555  anorm = one
556  GO TO 70
557 *
558  50 CONTINUE
559  anorm = ( rtovfl*ulp )*amninv
560  GO TO 70
561 *
562  60 CONTINUE
563  anorm = rtunfl*max( m, n )*ulpinv
564  GO TO 70
565 *
566  70 CONTINUE
567 *
568  CALL dlaset( 'Full', lda, n, zero, zero, a, lda )
569  CALL dlaset( 'Full', ldab, n, zero, zero, ab, ldab )
570  iinfo = 0
571  cond = ulpinv
572 *
573 * Special Matrices -- Identity & Jordan block
574 *
575 * Zero
576 *
577  IF( itype.EQ.1 ) THEN
578  iinfo = 0
579 *
580  ELSE IF( itype.EQ.2 ) THEN
581 *
582 * Identity
583 *
584  DO 80 jcol = 1, n
585  a( jcol, jcol ) = anorm
586  80 CONTINUE
587 *
588  ELSE IF( itype.EQ.4 ) THEN
589 *
590 * Diagonal Matrix, singular values specified
591 *
592  CALL dlatms( m, n, 'S', iseed, 'N', work, imode, cond,
593  $ anorm, 0, 0, 'N', a, lda, work( m+1 ),
594  $ iinfo )
595 *
596  ELSE IF( itype.EQ.6 ) THEN
597 *
598 * Nonhermitian, singular values specified
599 *
600  CALL dlatms( m, n, 'S', iseed, 'N', work, imode, cond,
601  $ anorm, kl, ku, 'N', a, lda, work( m+1 ),
602  $ iinfo )
603 *
604  ELSE IF( itype.EQ.9 ) THEN
605 *
606 * Nonhermitian, random entries
607 *
608  CALL dlatmr( m, n, 'S', iseed, 'N', work, 6, one, one,
609  $ 'T', 'N', work( n+1 ), 1, one,
610  $ work( 2*n+1 ), 1, one, 'N', idumma, kl,
611  $ ku, zero, anorm, 'N', a, lda, idumma,
612  $ iinfo )
613 *
614  ELSE
615 *
616  iinfo = 1
617  END IF
618 *
619 * Generate Right-Hand Side
620 *
621  CALL dlatmr( m, nrhs, 'S', iseed, 'N', work, 6, one, one,
622  $ 'T', 'N', work( m+1 ), 1, one,
623  $ work( 2*m+1 ), 1, one, 'N', idumma, m, nrhs,
624  $ zero, one, 'NO', c, ldc, idumma, iinfo )
625 *
626  IF( iinfo.NE.0 ) THEN
627  WRITE( nounit, fmt = 9999 )'Generator', iinfo, n,
628  $ jtype, ioldsd
629  info = abs( iinfo )
630  RETURN
631  END IF
632 *
633  90 CONTINUE
634 *
635 * Copy A to band storage.
636 *
637  DO 110 j = 1, n
638  DO 100 i = max( 1, j-ku ), min( m, j+kl )
639  ab( ku+1+i-j, j ) = a( i, j )
640  100 CONTINUE
641  110 CONTINUE
642 *
643 * Copy C
644 *
645  CALL dlacpy( 'Full', m, nrhs, c, ldc, cc, ldc )
646 *
647 * Call DGBBRD to compute B, Q and P, and to update C.
648 *
649  CALL dgbbrd( 'B', m, n, nrhs, kl, ku, ab, ldab, bd, be,
650  $ q, ldq, p, ldp, cc, ldc, work, iinfo )
651 *
652  IF( iinfo.NE.0 ) THEN
653  WRITE( nounit, fmt = 9999 )'DGBBRD', iinfo, n, jtype,
654  $ ioldsd
655  info = abs( iinfo )
656  IF( iinfo.LT.0 ) THEN
657  RETURN
658  ELSE
659  result( 1 ) = ulpinv
660  GO TO 120
661  END IF
662  END IF
663 *
664 * Test 1: Check the decomposition A := Q * B * P'
665 * 2: Check the orthogonality of Q
666 * 3: Check the orthogonality of P
667 * 4: Check the computation of Q' * C
668 *
669  CALL dbdt01( m, n, -1, a, lda, q, ldq, bd, be, p, ldp,
670  $ work, result( 1 ) )
671  CALL dort01( 'Columns', m, m, q, ldq, work, lwork,
672  $ result( 2 ) )
673  CALL dort01( 'Rows', n, n, p, ldp, work, lwork,
674  $ result( 3 ) )
675  CALL dbdt02( m, nrhs, c, ldc, cc, ldc, q, ldq, work,
676  $ result( 4 ) )
677 *
678 * End of Loop -- Check for RESULT(j) > THRESH
679 *
680  ntest = 4
681  120 CONTINUE
682  ntestt = ntestt + ntest
683 *
684 * Print out tests which fail.
685 *
686  DO 130 jr = 1, ntest
687  IF( result( jr ).GE.thresh ) THEN
688  IF( nerrs.EQ.0 )
689  $ CALL dlahd2( nounit, 'DBB' )
690  nerrs = nerrs + 1
691  WRITE( nounit, fmt = 9998 )m, n, k, ioldsd, jtype,
692  $ jr, result( jr )
693  END IF
694  130 CONTINUE
695 *
696  140 CONTINUE
697  150 CONTINUE
698  160 CONTINUE
699 *
700 * Summary
701 *
702  CALL dlasum( 'DBB', nounit, nerrs, ntestt )
703  RETURN
704 *
705  9999 FORMAT( ' DCHKBB: ', a, ' returned INFO=', i5, '.', / 9x, 'M=',
706  $ i5, ' N=', i5, ' K=', i5, ', JTYPE=', i5, ', ISEED=(',
707  $ 3( i5, ',' ), i5, ')' )
708  9998 FORMAT( ' M =', i4, ' N=', i4, ', K=', i3, ', seed=',
709  $ 4( i4, ',' ), ' type ', i2, ', test(', i2, ')=', g10.3 )
710 *
711 * End of DCHKBB
712 *
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dlatmr(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER, CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM, PACK, A, LDA, IWORK, INFO)
DLATMR
Definition: dlatmr.f:473
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dbdt01(M, N, KD, A, LDA, Q, LDQ, D, E, PT, LDPT, WORK, RESID)
DBDT01
Definition: dbdt01.f:142
subroutine dgbbrd(VECT, M, N, NCC, KL, KU, AB, LDAB, D, E, Q, LDQ, PT, LDPT, C, LDC, WORK, INFO)
DGBBRD
Definition: dgbbrd.f:189
subroutine dbdt02(M, N, B, LDB, C, LDC, U, LDU, WORK, RESID)
DBDT02
Definition: dbdt02.f:113
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:45
subroutine dlahd2(IOUNIT, PATH)
DLAHD2
Definition: dlahd2.f:67
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dort01(ROWCOL, M, N, U, LDU, WORK, LWORK, RESID)
DORT01
Definition: dort01.f:118

Here is the call graph for this function:

Here is the caller graph for this function: