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

◆ zchkbb()

subroutine zchkbb ( 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,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( ldab, * )  ab,
integer  ldab,
double precision, dimension( * )  bd,
double precision, dimension( * )  be,
complex*16, dimension( ldq, * )  q,
integer  ldq,
complex*16, dimension( ldp, * )  p,
integer  ldp,
complex*16, dimension( ldc, * )  c,
integer  ldc,
complex*16, dimension( ldc, * )  cc,
complex*16, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  rwork,
double precision, dimension( * )  result,
integer  info 
)

ZCHKBB

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

 ZGBBRD factors a general band matrix A as  Q B P* , where * means
 conjugate transpose, B is upper bidiagonal, and Q and P are unitary;
 ZGBBRD 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, ZCHKBB 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,
          ZCHKBB 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, ZCHKBB
          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 ZCHKBB 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 ZGBBRD.
[out]BE
          BE is DOUBLE PRECISION array, dimension (max(NN))
          Used to hold the off-diagonal of the bidiagonal matrix
          computed by ZGBBRD.
[out]Q
          Q is COMPLEX*16 array, dimension (LDQ, max(NN))
          Used to hold the unitary matrix Q computed by ZGBBRD.
[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 COMPLEX*16 array, dimension (LDP, max(NN))
          Used to hold the unitary matrix P computed by ZGBBRD.
[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 COMPLEX*16 array, dimension (LDC, max(NN))
          Used to hold the matrix C updated by ZGBBRD.
[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 COMPLEX*16 array, dimension (LDC, max(NN))
          Used to hold a copy of the matrix C.
[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( LDA+1, max(NN)+1 )*max(NN).
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (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.

Definition at line 357 of file zchkbb.f.

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