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

◆ slatmr()

subroutine slatmr ( integer  m,
integer  n,
character  dist,
integer, dimension( 4 )  iseed,
character  sym,
real, dimension( * )  d,
integer  mode,
real  cond,
real  dmax,
character  rsign,
character  grade,
real, dimension( * )  dl,
integer  model,
real  condl,
real, dimension( * )  dr,
integer  moder,
real  condr,
character  pivtng,
integer, dimension( * )  ipivot,
integer  kl,
integer  ku,
real  sparse,
real  anorm,
character  pack,
real, dimension( lda, * )  a,
integer  lda,
integer, dimension( * )  iwork,
integer  info 
)

SLATMR

Purpose:
    SLATMR generates random matrices of various types for testing
    LAPACK programs.

    SLATMR operates by applying the following sequence of
    operations:

      Generate a matrix A with random entries of distribution DIST
         which is symmetric if SYM='S', and nonsymmetric
         if SYM='N'.

      Set the diagonal to D, where D may be input or
         computed according to MODE, COND, DMAX and RSIGN
         as described below.

      Grade the matrix, if desired, from the left and/or right
         as specified by GRADE. The inputs DL, MODEL, CONDL, DR,
         MODER and CONDR also determine the grading as described
         below.

      Permute, if desired, the rows and/or columns as specified by
         PIVTNG and IPIVOT.

      Set random entries to zero, if desired, to get a random sparse
         matrix as specified by SPARSE.

      Make A a band matrix, if desired, by zeroing out the matrix
         outside a band of lower bandwidth KL and upper bandwidth KU.

      Scale A, if desired, to have maximum entry ANORM.

      Pack the matrix if desired. Options specified by PACK are:
         no packing
         zero out upper half (if symmetric)
         zero out lower half (if symmetric)
         store the upper half columnwise (if symmetric or
             square upper triangular)
         store the lower half columnwise (if symmetric or
             square lower triangular)
             same as upper half rowwise if symmetric
         store the lower triangle in banded format (if symmetric)
         store the upper triangle in banded format (if symmetric)
         store the entire matrix in banded format

    Note: If two calls to SLATMR differ only in the PACK parameter,
          they will generate mathematically equivalent matrices.

          If two calls to SLATMR both have full bandwidth (KL = M-1
          and KU = N-1), and differ only in the PIVTNG and PACK
          parameters, then the matrices generated will differ only
          in the order of the rows and/or columns, and otherwise
          contain the same data. This consistency cannot be and
          is not maintained with less than full bandwidth.
Parameters
[in]M
          M is INTEGER
           Number of rows of A. Not modified.
[in]N
          N is INTEGER
           Number of columns of A. Not modified.
[in]DIST
          DIST is CHARACTER*1
           On entry, DIST specifies the type of distribution to be used
           to generate a random matrix .
           'U' => UNIFORM( 0, 1 )  ( 'U' for uniform )
           'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
           'N' => NORMAL( 0, 1 )   ( 'N' for normal )
           Not modified.
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
           On entry ISEED specifies the seed of the random number
           generator. They should lie between 0 and 4095 inclusive,
           and ISEED(4) should 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 SLATMR
           to continue the same random number sequence.
           Changed on exit.
[in]SYM
          SYM is CHARACTER*1
           If SYM='S' or 'H', generated matrix is symmetric.
           If SYM='N', generated matrix is nonsymmetric.
           Not modified.
[in]D
          D is REAL array, dimension (min(M,N))
           On entry this array specifies the diagonal entries
           of the diagonal of A.  D may either be specified
           on entry, or set according to MODE and COND as described
           below. May be changed on exit if MODE is nonzero.
[in]MODE
          MODE is INTEGER
           On entry describes how D is to be used:
           MODE = 0 means use D as input
           MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
           MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
           MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
           MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
           MODE = 5 sets D to random numbers in the range
                    ( 1/COND , 1 ) such that their logarithms
                    are uniformly distributed.
           MODE = 6 set D to random numbers from same distribution
                    as the rest of the matrix.
           MODE < 0 has the same meaning as ABS(MODE), except that
              the order of the elements of D is reversed.
           Thus if MODE is positive, D has entries ranging from
              1 to 1/COND, if negative, from 1/COND to 1,
           Not modified.
[in]COND
          COND is REAL
           On entry, used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]DMAX
          DMAX is REAL
           If MODE neither -6, 0 nor 6, the diagonal is scaled by
           DMAX / max(abs(D(i))), so that maximum absolute entry
           of diagonal is abs(DMAX). If DMAX is negative (or zero),
           diagonal will be scaled by a negative number (or zero).
[in]RSIGN
          RSIGN is CHARACTER*1
           If MODE neither -6, 0 nor 6, specifies sign of diagonal
           as follows:
           'T' => diagonal entries are multiplied by 1 or -1
                  with probability .5
           'F' => diagonal unchanged
           Not modified.
[in]GRADE
          GRADE is CHARACTER*1
           Specifies grading of matrix as follows:
           'N'  => no grading
           'L'  => matrix premultiplied by diag( DL )
                   (only if matrix nonsymmetric)
           'R'  => matrix postmultiplied by diag( DR )
                   (only if matrix nonsymmetric)
           'B'  => matrix premultiplied by diag( DL ) and
                         postmultiplied by diag( DR )
                   (only if matrix nonsymmetric)
           'S' or 'H'  => matrix premultiplied by diag( DL ) and
                          postmultiplied by diag( DL )
                          ('S' for symmetric, or 'H' for Hermitian)
           'E'  => matrix premultiplied by diag( DL ) and
                         postmultiplied by inv( diag( DL ) )
                         ( 'E' for eigenvalue invariance)
                   (only if matrix nonsymmetric)
                   Note: if GRADE='E', then M must equal N.
           Not modified.
[in,out]DL
          DL is REAL array, dimension (M)
           If MODEL=0, then on entry this array specifies the diagonal
           entries of a diagonal matrix used as described under GRADE
           above. If MODEL is not zero, then DL will be set according
           to MODEL and CONDL, analogous to the way D is set according
           to MODE and COND (except there is no DMAX parameter for DL).
           If GRADE='E', then DL cannot have zero entries.
           Not referenced if GRADE = 'N' or 'R'. Changed on exit.
[in]MODEL
          MODEL is INTEGER
           This specifies how the diagonal array DL is to be computed,
           just as MODE specifies how D is to be computed.
           Not modified.
[in]CONDL
          CONDL is REAL
           When MODEL is not zero, this specifies the condition number
           of the computed DL.  Not modified.
[in,out]DR
          DR is REAL array, dimension (N)
           If MODER=0, then on entry this array specifies the diagonal
           entries of a diagonal matrix used as described under GRADE
           above. If MODER is not zero, then DR will be set according
           to MODER and CONDR, analogous to the way D is set according
           to MODE and COND (except there is no DMAX parameter for DR).
           Not referenced if GRADE = 'N', 'L', 'H', 'S' or 'E'.
           Changed on exit.
[in]MODER
          MODER is INTEGER
           This specifies how the diagonal array DR is to be computed,
           just as MODE specifies how D is to be computed.
           Not modified.
[in]CONDR
          CONDR is REAL
           When MODER is not zero, this specifies the condition number
           of the computed DR.  Not modified.
[in]PIVTNG
          PIVTNG is CHARACTER*1
           On entry specifies pivoting permutations as follows:
           'N' or ' ' => none.
           'L' => left or row pivoting (matrix must be nonsymmetric).
           'R' => right or column pivoting (matrix must be
                  nonsymmetric).
           'B' or 'F' => both or full pivoting, i.e., on both sides.
                         In this case, M must equal N

           If two calls to SLATMR both have full bandwidth (KL = M-1
           and KU = N-1), and differ only in the PIVTNG and PACK
           parameters, then the matrices generated will differ only
           in the order of the rows and/or columns, and otherwise
           contain the same data. This consistency cannot be
           maintained with less than full bandwidth.
[in]IPIVOT
          IPIVOT is INTEGER array, dimension (N or M)
           This array specifies the permutation used.  After the
           basic matrix is generated, the rows, columns, or both
           are permuted.   If, say, row pivoting is selected, SLATMR
           starts with the *last* row and interchanges the M-th and
           IPIVOT(M)-th rows, then moves to the next-to-last row,
           interchanging the (M-1)-th and the IPIVOT(M-1)-th rows,
           and so on.  In terms of "2-cycles", the permutation is
           (1 IPIVOT(1)) (2 IPIVOT(2)) ... (M IPIVOT(M))
           where the rightmost cycle is applied first.  This is the
           *inverse* of the effect of pivoting in LINPACK.  The idea
           is that factoring (with pivoting) an identity matrix
           which has been inverse-pivoted in this way should
           result in a pivot vector identical to IPIVOT.
           Not referenced if PIVTNG = 'N'. Not modified.
[in]KL
          KL is INTEGER
           On entry specifies the lower bandwidth of the  matrix. For
           example, KL=0 implies upper triangular, KL=1 implies upper
           Hessenberg, and KL at least M-1 implies the matrix is not
           banded. Must equal KU if matrix is symmetric.
           Not modified.
[in]KU
          KU is INTEGER
           On entry specifies the upper bandwidth of the  matrix. For
           example, KU=0 implies lower triangular, KU=1 implies lower
           Hessenberg, and KU at least N-1 implies the matrix is not
           banded. Must equal KL if matrix is symmetric.
           Not modified.
[in]SPARSE
          SPARSE is REAL
           On entry specifies the sparsity of the matrix if a sparse
           matrix is to be generated. SPARSE should lie between
           0 and 1. To generate a sparse matrix, for each matrix entry
           a uniform ( 0, 1 ) random number x is generated and
           compared to SPARSE; if x is larger the matrix entry
           is unchanged and if x is smaller the entry is set
           to zero. Thus on the average a fraction SPARSE of the
           entries will be set to zero.
           Not modified.
[in]ANORM
          ANORM is REAL
           On entry specifies maximum entry of output matrix
           (output matrix will by multiplied by a constant so that
           its largest absolute entry equal ANORM)
           if ANORM is nonnegative. If ANORM is negative no scaling
           is done. Not modified.
[in]PACK
          PACK is CHARACTER*1
           On entry specifies packing of matrix as follows:
           'N' => no packing
           'U' => zero out all subdiagonal entries (if symmetric)
           'L' => zero out all superdiagonal entries (if symmetric)
           'C' => store the upper triangle columnwise
                  (only if matrix symmetric or square upper triangular)
           'R' => store the lower triangle columnwise
                  (only if matrix symmetric or square lower triangular)
                  (same as upper half rowwise if symmetric)
           'B' => store the lower triangle in band storage scheme
                  (only if matrix symmetric)
           'Q' => store the upper triangle in band storage scheme
                  (only if matrix symmetric)
           'Z' => store the entire matrix in band storage scheme
                      (pivoting can be provided for by using this
                      option to store A in the trailing rows of
                      the allocated storage)

           Using these options, the various LAPACK packed and banded
           storage schemes can be obtained:
           GB               - use 'Z'
           PB, SB or TB     - use 'B' or 'Q'
           PP, SP or TP     - use 'C' or 'R'

           If two calls to SLATMR differ only in the PACK parameter,
           they will generate mathematically equivalent matrices.
           Not modified.
[in,out]A
          A is REAL array, dimension (LDA,N)
           On exit A is the desired test matrix. Only those
           entries of A which are significant on output
           will be referenced (even if A is in packed or band
           storage format). The 'unoccupied corners' of A in
           band format will be zeroed out.
[in]LDA
          LDA is INTEGER
           on entry LDA specifies the first dimension of A as
           declared in the calling program.
           If PACK='N', 'U' or 'L', LDA must be at least max ( 1, M ).
           If PACK='C' or 'R', LDA must be at least 1.
           If PACK='B', or 'Q', LDA must be MIN ( KU+1, N )
           If PACK='Z', LDA must be at least KUU+KLL+1, where
           KUU = MIN ( KU, N-1 ) and KLL = MIN ( KL, M-1 )
           Not modified.
[out]IWORK
          IWORK is INTEGER array, dimension ( N or M)
           Workspace. Not referenced if PIVTNG = 'N'. Changed on exit.
[out]INFO
          INFO is INTEGER
           Error parameter on exit:
             0 => normal return
            -1 => M negative or unequal to N and SYM='S' or 'H'
            -2 => N negative
            -3 => DIST illegal string
            -5 => SYM illegal string
            -7 => MODE not in range -6 to 6
            -8 => COND less than 1.0, and MODE neither -6, 0 nor 6
           -10 => MODE neither -6, 0 nor 6 and RSIGN illegal string
           -11 => GRADE illegal string, or GRADE='E' and
                  M not equal to N, or GRADE='L', 'R', 'B' or 'E' and
                  SYM = 'S' or 'H'
           -12 => GRADE = 'E' and DL contains zero
           -13 => MODEL not in range -6 to 6 and GRADE= 'L', 'B', 'H',
                  'S' or 'E'
           -14 => CONDL less than 1.0, GRADE='L', 'B', 'H', 'S' or 'E',
                  and MODEL neither -6, 0 nor 6
           -16 => MODER not in range -6 to 6 and GRADE= 'R' or 'B'
           -17 => CONDR less than 1.0, GRADE='R' or 'B', and
                  MODER neither -6, 0 nor 6
           -18 => PIVTNG illegal string, or PIVTNG='B' or 'F' and
                  M not equal to N, or PIVTNG='L' or 'R' and SYM='S'
                  or 'H'
           -19 => IPIVOT contains out of range number and
                  PIVTNG not equal to 'N'
           -20 => KL negative
           -21 => KU negative, or SYM='S' or 'H' and KU not equal to KL
           -22 => SPARSE not in range 0. to 1.
           -24 => PACK illegal string, or PACK='U', 'L', 'B' or 'Q'
                  and SYM='N', or PACK='C' and SYM='N' and either KL
                  not equal to 0 or N not equal to M, or PACK='R' and
                  SYM='N', and either KU not equal to 0 or N not equal
                  to M
           -26 => LDA too small
             1 => Error return from SLATM1 (computing D)
             2 => Cannot scale diagonal to DMAX (max. entry is 0)
             3 => Error return from SLATM1 (computing DL)
             4 => Error return from SLATM1 (computing DR)
             5 => ANORM is positive, but matrix constructed prior to
                  attempting to scale it to have norm ANORM, is zero
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 467 of file slatmr.f.

471*
472* -- LAPACK computational routine --
473* -- LAPACK is a software package provided by Univ. of Tennessee, --
474* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
475*
476* .. Scalar Arguments ..
477 CHARACTER DIST, GRADE, PACK, PIVTNG, RSIGN, SYM
478 INTEGER INFO, KL, KU, LDA, M, MODE, MODEL, MODER, N
479 REAL ANORM, COND, CONDL, CONDR, DMAX, SPARSE
480* ..
481* .. Array Arguments ..
482 INTEGER IPIVOT( * ), ISEED( 4 ), IWORK( * )
483 REAL A( LDA, * ), D( * ), DL( * ), DR( * )
484* ..
485*
486* =====================================================================
487*
488* .. Parameters ..
489 REAL ZERO
490 parameter( zero = 0.0e0 )
491 REAL ONE
492 parameter( one = 1.0e0 )
493* ..
494* .. Local Scalars ..
495 LOGICAL BADPVT, DZERO, FULBND
496 INTEGER I, IDIST, IGRADE, IISUB, IPACK, IPVTNG, IRSIGN,
497 $ ISUB, ISYM, J, JJSUB, JSUB, K, KLL, KUU, MNMIN,
498 $ MNSUB, MXSUB, NPVTS
499 REAL ALPHA, ONORM, TEMP
500* ..
501* .. Local Arrays ..
502 REAL TEMPA( 1 )
503* ..
504* .. External Functions ..
505 LOGICAL LSAME
506 REAL SLANGB, SLANGE, SLANSB, SLANSP, SLANSY, SLATM2,
507 $ SLATM3
508 EXTERNAL lsame, slangb, slange, slansb, slansp, slansy,
509 $ slatm2, slatm3
510* ..
511* .. External Subroutines ..
512 EXTERNAL slatm1, sscal, xerbla
513* ..
514* .. Intrinsic Functions ..
515 INTRINSIC abs, max, min, mod
516* ..
517* .. Executable Statements ..
518*
519* 1) Decode and Test the input parameters.
520* Initialize flags & seed.
521*
522 info = 0
523*
524* Quick return if possible
525*
526 IF( m.EQ.0 .OR. n.EQ.0 )
527 $ RETURN
528*
529* Decode DIST
530*
531 IF( lsame( dist, 'U' ) ) THEN
532 idist = 1
533 ELSE IF( lsame( dist, 'S' ) ) THEN
534 idist = 2
535 ELSE IF( lsame( dist, 'N' ) ) THEN
536 idist = 3
537 ELSE
538 idist = -1
539 END IF
540*
541* Decode SYM
542*
543 IF( lsame( sym, 'S' ) ) THEN
544 isym = 0
545 ELSE IF( lsame( sym, 'N' ) ) THEN
546 isym = 1
547 ELSE IF( lsame( sym, 'H' ) ) THEN
548 isym = 0
549 ELSE
550 isym = -1
551 END IF
552*
553* Decode RSIGN
554*
555 IF( lsame( rsign, 'F' ) ) THEN
556 irsign = 0
557 ELSE IF( lsame( rsign, 'T' ) ) THEN
558 irsign = 1
559 ELSE
560 irsign = -1
561 END IF
562*
563* Decode PIVTNG
564*
565 IF( lsame( pivtng, 'N' ) ) THEN
566 ipvtng = 0
567 ELSE IF( lsame( pivtng, ' ' ) ) THEN
568 ipvtng = 0
569 ELSE IF( lsame( pivtng, 'L' ) ) THEN
570 ipvtng = 1
571 npvts = m
572 ELSE IF( lsame( pivtng, 'R' ) ) THEN
573 ipvtng = 2
574 npvts = n
575 ELSE IF( lsame( pivtng, 'B' ) ) THEN
576 ipvtng = 3
577 npvts = min( n, m )
578 ELSE IF( lsame( pivtng, 'F' ) ) THEN
579 ipvtng = 3
580 npvts = min( n, m )
581 ELSE
582 ipvtng = -1
583 END IF
584*
585* Decode GRADE
586*
587 IF( lsame( grade, 'N' ) ) THEN
588 igrade = 0
589 ELSE IF( lsame( grade, 'L' ) ) THEN
590 igrade = 1
591 ELSE IF( lsame( grade, 'R' ) ) THEN
592 igrade = 2
593 ELSE IF( lsame( grade, 'B' ) ) THEN
594 igrade = 3
595 ELSE IF( lsame( grade, 'E' ) ) THEN
596 igrade = 4
597 ELSE IF( lsame( grade, 'H' ) .OR. lsame( grade, 'S' ) ) THEN
598 igrade = 5
599 ELSE
600 igrade = -1
601 END IF
602*
603* Decode PACK
604*
605 IF( lsame( pack, 'N' ) ) THEN
606 ipack = 0
607 ELSE IF( lsame( pack, 'U' ) ) THEN
608 ipack = 1
609 ELSE IF( lsame( pack, 'L' ) ) THEN
610 ipack = 2
611 ELSE IF( lsame( pack, 'C' ) ) THEN
612 ipack = 3
613 ELSE IF( lsame( pack, 'R' ) ) THEN
614 ipack = 4
615 ELSE IF( lsame( pack, 'B' ) ) THEN
616 ipack = 5
617 ELSE IF( lsame( pack, 'Q' ) ) THEN
618 ipack = 6
619 ELSE IF( lsame( pack, 'Z' ) ) THEN
620 ipack = 7
621 ELSE
622 ipack = -1
623 END IF
624*
625* Set certain internal parameters
626*
627 mnmin = min( m, n )
628 kll = min( kl, m-1 )
629 kuu = min( ku, n-1 )
630*
631* If inv(DL) is used, check to see if DL has a zero entry.
632*
633 dzero = .false.
634 IF( igrade.EQ.4 .AND. model.EQ.0 ) THEN
635 DO 10 i = 1, m
636 IF( dl( i ).EQ.zero )
637 $ dzero = .true.
638 10 CONTINUE
639 END IF
640*
641* Check values in IPIVOT
642*
643 badpvt = .false.
644 IF( ipvtng.GT.0 ) THEN
645 DO 20 j = 1, npvts
646 IF( ipivot( j ).LE.0 .OR. ipivot( j ).GT.npvts )
647 $ badpvt = .true.
648 20 CONTINUE
649 END IF
650*
651* Set INFO if an error
652*
653 IF( m.LT.0 ) THEN
654 info = -1
655 ELSE IF( m.NE.n .AND. isym.EQ.0 ) THEN
656 info = -1
657 ELSE IF( n.LT.0 ) THEN
658 info = -2
659 ELSE IF( idist.EQ.-1 ) THEN
660 info = -3
661 ELSE IF( isym.EQ.-1 ) THEN
662 info = -5
663 ELSE IF( mode.LT.-6 .OR. mode.GT.6 ) THEN
664 info = -7
665 ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
666 $ cond.LT.one ) THEN
667 info = -8
668 ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
669 $ irsign.EQ.-1 ) THEN
670 info = -10
671 ELSE IF( igrade.EQ.-1 .OR. ( igrade.EQ.4 .AND. m.NE.n ) .OR.
672 $ ( ( igrade.GE.1 .AND. igrade.LE.4 ) .AND. isym.EQ.0 ) )
673 $ THEN
674 info = -11
675 ELSE IF( igrade.EQ.4 .AND. dzero ) THEN
676 info = -12
677 ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
678 $ igrade.EQ.5 ) .AND. ( model.LT.-6 .OR. model.GT.6 ) )
679 $ THEN
680 info = -13
681 ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
682 $ igrade.EQ.5 ) .AND. ( model.NE.-6 .AND. model.NE.0 .AND.
683 $ model.NE.6 ) .AND. condl.LT.one ) THEN
684 info = -14
685 ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
686 $ ( moder.LT.-6 .OR. moder.GT.6 ) ) THEN
687 info = -16
688 ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
689 $ ( moder.NE.-6 .AND. moder.NE.0 .AND. moder.NE.6 ) .AND.
690 $ condr.LT.one ) THEN
691 info = -17
692 ELSE IF( ipvtng.EQ.-1 .OR. ( ipvtng.EQ.3 .AND. m.NE.n ) .OR.
693 $ ( ( ipvtng.EQ.1 .OR. ipvtng.EQ.2 ) .AND. isym.EQ.0 ) )
694 $ THEN
695 info = -18
696 ELSE IF( ipvtng.NE.0 .AND. badpvt ) THEN
697 info = -19
698 ELSE IF( kl.LT.0 ) THEN
699 info = -20
700 ELSE IF( ku.LT.0 .OR. ( isym.EQ.0 .AND. kl.NE.ku ) ) THEN
701 info = -21
702 ELSE IF( sparse.LT.zero .OR. sparse.GT.one ) THEN
703 info = -22
704 ELSE IF( ipack.EQ.-1 .OR. ( ( ipack.EQ.1 .OR. ipack.EQ.2 .OR.
705 $ ipack.EQ.5 .OR. ipack.EQ.6 ) .AND. isym.EQ.1 ) .OR.
706 $ ( ipack.EQ.3 .AND. isym.EQ.1 .AND. ( kl.NE.0 .OR. m.NE.
707 $ n ) ) .OR. ( ipack.EQ.4 .AND. isym.EQ.1 .AND. ( ku.NE.
708 $ 0 .OR. m.NE.n ) ) ) THEN
709 info = -24
710 ELSE IF( ( ( ipack.EQ.0 .OR. ipack.EQ.1 .OR. ipack.EQ.2 ) .AND.
711 $ lda.LT.max( 1, m ) ) .OR. ( ( ipack.EQ.3 .OR. ipack.EQ.
712 $ 4 ) .AND. lda.LT.1 ) .OR. ( ( ipack.EQ.5 .OR. ipack.EQ.
713 $ 6 ) .AND. lda.LT.kuu+1 ) .OR.
714 $ ( ipack.EQ.7 .AND. lda.LT.kll+kuu+1 ) ) THEN
715 info = -26
716 END IF
717*
718 IF( info.NE.0 ) THEN
719 CALL xerbla( 'SLATMR', -info )
720 RETURN
721 END IF
722*
723* Decide if we can pivot consistently
724*
725 fulbnd = .false.
726 IF( kuu.EQ.n-1 .AND. kll.EQ.m-1 )
727 $ fulbnd = .true.
728*
729* Initialize random number generator
730*
731 DO 30 i = 1, 4
732 iseed( i ) = mod( abs( iseed( i ) ), 4096 )
733 30 CONTINUE
734*
735 iseed( 4 ) = 2*( iseed( 4 ) / 2 ) + 1
736*
737* 2) Set up D, DL, and DR, if indicated.
738*
739* Compute D according to COND and MODE
740*
741 CALL slatm1( mode, cond, irsign, idist, iseed, d, mnmin, info )
742 IF( info.NE.0 ) THEN
743 info = 1
744 RETURN
745 END IF
746 IF( mode.NE.0 .AND. mode.NE.-6 .AND. mode.NE.6 ) THEN
747*
748* Scale by DMAX
749*
750 temp = abs( d( 1 ) )
751 DO 40 i = 2, mnmin
752 temp = max( temp, abs( d( i ) ) )
753 40 CONTINUE
754 IF( temp.EQ.zero .AND. dmax.NE.zero ) THEN
755 info = 2
756 RETURN
757 END IF
758 IF( temp.NE.zero ) THEN
759 alpha = dmax / temp
760 ELSE
761 alpha = one
762 END IF
763 DO 50 i = 1, mnmin
764 d( i ) = alpha*d( i )
765 50 CONTINUE
766*
767 END IF
768*
769* Compute DL if grading set
770*
771 IF( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR. igrade.EQ.
772 $ 5 ) THEN
773 CALL slatm1( model, condl, 0, idist, iseed, dl, m, info )
774 IF( info.NE.0 ) THEN
775 info = 3
776 RETURN
777 END IF
778 END IF
779*
780* Compute DR if grading set
781*
782 IF( igrade.EQ.2 .OR. igrade.EQ.3 ) THEN
783 CALL slatm1( moder, condr, 0, idist, iseed, dr, n, info )
784 IF( info.NE.0 ) THEN
785 info = 4
786 RETURN
787 END IF
788 END IF
789*
790* 3) Generate IWORK if pivoting
791*
792 IF( ipvtng.GT.0 ) THEN
793 DO 60 i = 1, npvts
794 iwork( i ) = i
795 60 CONTINUE
796 IF( fulbnd ) THEN
797 DO 70 i = 1, npvts
798 k = ipivot( i )
799 j = iwork( i )
800 iwork( i ) = iwork( k )
801 iwork( k ) = j
802 70 CONTINUE
803 ELSE
804 DO 80 i = npvts, 1, -1
805 k = ipivot( i )
806 j = iwork( i )
807 iwork( i ) = iwork( k )
808 iwork( k ) = j
809 80 CONTINUE
810 END IF
811 END IF
812*
813* 4) Generate matrices for each kind of PACKing
814* Always sweep matrix columnwise (if symmetric, upper
815* half only) so that matrix generated does not depend
816* on PACK
817*
818 IF( fulbnd ) THEN
819*
820* Use SLATM3 so matrices generated with differing PIVOTing only
821* differ only in the order of their rows and/or columns.
822*
823 IF( ipack.EQ.0 ) THEN
824 IF( isym.EQ.0 ) THEN
825 DO 100 j = 1, n
826 DO 90 i = 1, j
827 temp = slatm3( m, n, i, j, isub, jsub, kl, ku,
828 $ idist, iseed, d, igrade, dl, dr, ipvtng,
829 $ iwork, sparse )
830 a( isub, jsub ) = temp
831 a( jsub, isub ) = temp
832 90 CONTINUE
833 100 CONTINUE
834 ELSE IF( isym.EQ.1 ) THEN
835 DO 120 j = 1, n
836 DO 110 i = 1, m
837 temp = slatm3( m, n, i, j, isub, jsub, kl, ku,
838 $ idist, iseed, d, igrade, dl, dr, ipvtng,
839 $ iwork, sparse )
840 a( isub, jsub ) = temp
841 110 CONTINUE
842 120 CONTINUE
843 END IF
844*
845 ELSE IF( ipack.EQ.1 ) THEN
846*
847 DO 140 j = 1, n
848 DO 130 i = 1, j
849 temp = slatm3( m, n, i, j, isub, jsub, kl, ku, idist,
850 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
851 $ sparse )
852 mnsub = min( isub, jsub )
853 mxsub = max( isub, jsub )
854 a( mnsub, mxsub ) = temp
855 IF( mnsub.NE.mxsub )
856 $ a( mxsub, mnsub ) = zero
857 130 CONTINUE
858 140 CONTINUE
859*
860 ELSE IF( ipack.EQ.2 ) THEN
861*
862 DO 160 j = 1, n
863 DO 150 i = 1, j
864 temp = slatm3( m, n, i, j, isub, jsub, kl, ku, idist,
865 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
866 $ sparse )
867 mnsub = min( isub, jsub )
868 mxsub = max( isub, jsub )
869 a( mxsub, mnsub ) = temp
870 IF( mnsub.NE.mxsub )
871 $ a( mnsub, mxsub ) = zero
872 150 CONTINUE
873 160 CONTINUE
874*
875 ELSE IF( ipack.EQ.3 ) THEN
876*
877 DO 180 j = 1, n
878 DO 170 i = 1, j
879 temp = slatm3( m, n, i, j, isub, jsub, kl, ku, idist,
880 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
881 $ sparse )
882*
883* Compute K = location of (ISUB,JSUB) entry in packed
884* array
885*
886 mnsub = min( isub, jsub )
887 mxsub = max( isub, jsub )
888 k = mxsub*( mxsub-1 ) / 2 + mnsub
889*
890* Convert K to (IISUB,JJSUB) location
891*
892 jjsub = ( k-1 ) / lda + 1
893 iisub = k - lda*( jjsub-1 )
894*
895 a( iisub, jjsub ) = temp
896 170 CONTINUE
897 180 CONTINUE
898*
899 ELSE IF( ipack.EQ.4 ) THEN
900*
901 DO 200 j = 1, n
902 DO 190 i = 1, j
903 temp = slatm3( m, n, i, j, isub, jsub, kl, ku, idist,
904 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
905 $ sparse )
906*
907* Compute K = location of (I,J) entry in packed array
908*
909 mnsub = min( isub, jsub )
910 mxsub = max( isub, jsub )
911 IF( mnsub.EQ.1 ) THEN
912 k = mxsub
913 ELSE
914 k = n*( n+1 ) / 2 - ( n-mnsub+1 )*( n-mnsub+2 ) /
915 $ 2 + mxsub - mnsub + 1
916 END IF
917*
918* Convert K to (IISUB,JJSUB) location
919*
920 jjsub = ( k-1 ) / lda + 1
921 iisub = k - lda*( jjsub-1 )
922*
923 a( iisub, jjsub ) = temp
924 190 CONTINUE
925 200 CONTINUE
926*
927 ELSE IF( ipack.EQ.5 ) THEN
928*
929 DO 220 j = 1, n
930 DO 210 i = j - kuu, j
931 IF( i.LT.1 ) THEN
932 a( j-i+1, i+n ) = zero
933 ELSE
934 temp = slatm3( m, n, i, j, isub, jsub, kl, ku,
935 $ idist, iseed, d, igrade, dl, dr, ipvtng,
936 $ iwork, sparse )
937 mnsub = min( isub, jsub )
938 mxsub = max( isub, jsub )
939 a( mxsub-mnsub+1, mnsub ) = temp
940 END IF
941 210 CONTINUE
942 220 CONTINUE
943*
944 ELSE IF( ipack.EQ.6 ) THEN
945*
946 DO 240 j = 1, n
947 DO 230 i = j - kuu, j
948 temp = slatm3( m, n, i, j, isub, jsub, kl, ku, idist,
949 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
950 $ sparse )
951 mnsub = min( isub, jsub )
952 mxsub = max( isub, jsub )
953 a( mnsub-mxsub+kuu+1, mxsub ) = temp
954 230 CONTINUE
955 240 CONTINUE
956*
957 ELSE IF( ipack.EQ.7 ) THEN
958*
959 IF( isym.EQ.0 ) THEN
960 DO 260 j = 1, n
961 DO 250 i = j - kuu, j
962 temp = slatm3( m, n, i, j, isub, jsub, kl, ku,
963 $ idist, iseed, d, igrade, dl, dr, ipvtng,
964 $ iwork, sparse )
965 mnsub = min( isub, jsub )
966 mxsub = max( isub, jsub )
967 a( mnsub-mxsub+kuu+1, mxsub ) = temp
968 IF( i.LT.1 )
969 $ a( j-i+1+kuu, i+n ) = zero
970 IF( i.GE.1 .AND. mnsub.NE.mxsub )
971 $ a( mxsub-mnsub+1+kuu, mnsub ) = temp
972 250 CONTINUE
973 260 CONTINUE
974 ELSE IF( isym.EQ.1 ) THEN
975 DO 280 j = 1, n
976 DO 270 i = j - kuu, j + kll
977 temp = slatm3( m, n, i, j, isub, jsub, kl, ku,
978 $ idist, iseed, d, igrade, dl, dr, ipvtng,
979 $ iwork, sparse )
980 a( isub-jsub+kuu+1, jsub ) = temp
981 270 CONTINUE
982 280 CONTINUE
983 END IF
984*
985 END IF
986*
987 ELSE
988*
989* Use SLATM2
990*
991 IF( ipack.EQ.0 ) THEN
992 IF( isym.EQ.0 ) THEN
993 DO 300 j = 1, n
994 DO 290 i = 1, j
995 a( i, j ) = slatm2( m, n, i, j, kl, ku, idist,
996 $ iseed, d, igrade, dl, dr, ipvtng,
997 $ iwork, sparse )
998 a( j, i ) = a( i, j )
999 290 CONTINUE
1000 300 CONTINUE
1001 ELSE IF( isym.EQ.1 ) THEN
1002 DO 320 j = 1, n
1003 DO 310 i = 1, m
1004 a( i, j ) = slatm2( m, n, i, j, kl, ku, idist,
1005 $ iseed, d, igrade, dl, dr, ipvtng,
1006 $ iwork, sparse )
1007 310 CONTINUE
1008 320 CONTINUE
1009 END IF
1010*
1011 ELSE IF( ipack.EQ.1 ) THEN
1012*
1013 DO 340 j = 1, n
1014 DO 330 i = 1, j
1015 a( i, j ) = slatm2( m, n, i, j, kl, ku, idist, iseed,
1016 $ d, igrade, dl, dr, ipvtng, iwork, sparse )
1017 IF( i.NE.j )
1018 $ a( j, i ) = zero
1019 330 CONTINUE
1020 340 CONTINUE
1021*
1022 ELSE IF( ipack.EQ.2 ) THEN
1023*
1024 DO 360 j = 1, n
1025 DO 350 i = 1, j
1026 a( j, i ) = slatm2( m, n, i, j, kl, ku, idist, iseed,
1027 $ d, igrade, dl, dr, ipvtng, iwork, sparse )
1028 IF( i.NE.j )
1029 $ a( i, j ) = zero
1030 350 CONTINUE
1031 360 CONTINUE
1032*
1033 ELSE IF( ipack.EQ.3 ) THEN
1034*
1035 isub = 0
1036 jsub = 1
1037 DO 380 j = 1, n
1038 DO 370 i = 1, j
1039 isub = isub + 1
1040 IF( isub.GT.lda ) THEN
1041 isub = 1
1042 jsub = jsub + 1
1043 END IF
1044 a( isub, jsub ) = slatm2( m, n, i, j, kl, ku, idist,
1045 $ iseed, d, igrade, dl, dr, ipvtng,
1046 $ iwork, sparse )
1047 370 CONTINUE
1048 380 CONTINUE
1049*
1050 ELSE IF( ipack.EQ.4 ) THEN
1051*
1052 IF( isym.EQ.0 ) THEN
1053 DO 400 j = 1, n
1054 DO 390 i = 1, j
1055*
1056* Compute K = location of (I,J) entry in packed array
1057*
1058 IF( i.EQ.1 ) THEN
1059 k = j
1060 ELSE
1061 k = n*( n+1 ) / 2 - ( n-i+1 )*( n-i+2 ) / 2 +
1062 $ j - i + 1
1063 END IF
1064*
1065* Convert K to (ISUB,JSUB) location
1066*
1067 jsub = ( k-1 ) / lda + 1
1068 isub = k - lda*( jsub-1 )
1069*
1070 a( isub, jsub ) = slatm2( m, n, i, j, kl, ku,
1071 $ idist, iseed, d, igrade, dl, dr,
1072 $ ipvtng, iwork, sparse )
1073 390 CONTINUE
1074 400 CONTINUE
1075 ELSE
1076 isub = 0
1077 jsub = 1
1078 DO 420 j = 1, n
1079 DO 410 i = j, m
1080 isub = isub + 1
1081 IF( isub.GT.lda ) THEN
1082 isub = 1
1083 jsub = jsub + 1
1084 END IF
1085 a( isub, jsub ) = slatm2( m, n, i, j, kl, ku,
1086 $ idist, iseed, d, igrade, dl, dr,
1087 $ ipvtng, iwork, sparse )
1088 410 CONTINUE
1089 420 CONTINUE
1090 END IF
1091*
1092 ELSE IF( ipack.EQ.5 ) THEN
1093*
1094 DO 440 j = 1, n
1095 DO 430 i = j - kuu, j
1096 IF( i.LT.1 ) THEN
1097 a( j-i+1, i+n ) = zero
1098 ELSE
1099 a( j-i+1, i ) = slatm2( m, n, i, j, kl, ku, idist,
1100 $ iseed, d, igrade, dl, dr, ipvtng,
1101 $ iwork, sparse )
1102 END IF
1103 430 CONTINUE
1104 440 CONTINUE
1105*
1106 ELSE IF( ipack.EQ.6 ) THEN
1107*
1108 DO 460 j = 1, n
1109 DO 450 i = j - kuu, j
1110 a( i-j+kuu+1, j ) = slatm2( m, n, i, j, kl, ku, idist,
1111 $ iseed, d, igrade, dl, dr, ipvtng,
1112 $ iwork, sparse )
1113 450 CONTINUE
1114 460 CONTINUE
1115*
1116 ELSE IF( ipack.EQ.7 ) THEN
1117*
1118 IF( isym.EQ.0 ) THEN
1119 DO 480 j = 1, n
1120 DO 470 i = j - kuu, j
1121 a( i-j+kuu+1, j ) = slatm2( m, n, i, j, kl, ku,
1122 $ idist, iseed, d, igrade, dl,
1123 $ dr, ipvtng, iwork, sparse )
1124 IF( i.LT.1 )
1125 $ a( j-i+1+kuu, i+n ) = zero
1126 IF( i.GE.1 .AND. i.NE.j )
1127 $ a( j-i+1+kuu, i ) = a( i-j+kuu+1, j )
1128 470 CONTINUE
1129 480 CONTINUE
1130 ELSE IF( isym.EQ.1 ) THEN
1131 DO 500 j = 1, n
1132 DO 490 i = j - kuu, j + kll
1133 a( i-j+kuu+1, j ) = slatm2( m, n, i, j, kl, ku,
1134 $ idist, iseed, d, igrade, dl,
1135 $ dr, ipvtng, iwork, sparse )
1136 490 CONTINUE
1137 500 CONTINUE
1138 END IF
1139*
1140 END IF
1141*
1142 END IF
1143*
1144* 5) Scaling the norm
1145*
1146 IF( ipack.EQ.0 ) THEN
1147 onorm = slange( 'M', m, n, a, lda, tempa )
1148 ELSE IF( ipack.EQ.1 ) THEN
1149 onorm = slansy( 'M', 'U', n, a, lda, tempa )
1150 ELSE IF( ipack.EQ.2 ) THEN
1151 onorm = slansy( 'M', 'L', n, a, lda, tempa )
1152 ELSE IF( ipack.EQ.3 ) THEN
1153 onorm = slansp( 'M', 'U', n, a, tempa )
1154 ELSE IF( ipack.EQ.4 ) THEN
1155 onorm = slansp( 'M', 'L', n, a, tempa )
1156 ELSE IF( ipack.EQ.5 ) THEN
1157 onorm = slansb( 'M', 'L', n, kll, a, lda, tempa )
1158 ELSE IF( ipack.EQ.6 ) THEN
1159 onorm = slansb( 'M', 'U', n, kuu, a, lda, tempa )
1160 ELSE IF( ipack.EQ.7 ) THEN
1161 onorm = slangb( 'M', n, kll, kuu, a, lda, tempa )
1162 END IF
1163*
1164 IF( anorm.GE.zero ) THEN
1165*
1166 IF( anorm.GT.zero .AND. onorm.EQ.zero ) THEN
1167*
1168* Desired scaling impossible
1169*
1170 info = 5
1171 RETURN
1172*
1173 ELSE IF( ( anorm.GT.one .AND. onorm.LT.one ) .OR.
1174 $ ( anorm.LT.one .AND. onorm.GT.one ) ) THEN
1175*
1176* Scale carefully to avoid over / underflow
1177*
1178 IF( ipack.LE.2 ) THEN
1179 DO 510 j = 1, n
1180 CALL sscal( m, one / onorm, a( 1, j ), 1 )
1181 CALL sscal( m, anorm, a( 1, j ), 1 )
1182 510 CONTINUE
1183*
1184 ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1185*
1186 CALL sscal( n*( n+1 ) / 2, one / onorm, a, 1 )
1187 CALL sscal( n*( n+1 ) / 2, anorm, a, 1 )
1188*
1189 ELSE IF( ipack.GE.5 ) THEN
1190*
1191 DO 520 j = 1, n
1192 CALL sscal( kll+kuu+1, one / onorm, a( 1, j ), 1 )
1193 CALL sscal( kll+kuu+1, anorm, a( 1, j ), 1 )
1194 520 CONTINUE
1195*
1196 END IF
1197*
1198 ELSE
1199*
1200* Scale straightforwardly
1201*
1202 IF( ipack.LE.2 ) THEN
1203 DO 530 j = 1, n
1204 CALL sscal( m, anorm / onorm, a( 1, j ), 1 )
1205 530 CONTINUE
1206*
1207 ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1208*
1209 CALL sscal( n*( n+1 ) / 2, anorm / onorm, a, 1 )
1210*
1211 ELSE IF( ipack.GE.5 ) THEN
1212*
1213 DO 540 j = 1, n
1214 CALL sscal( kll+kuu+1, anorm / onorm, a( 1, j ), 1 )
1215 540 CONTINUE
1216 END IF
1217*
1218 END IF
1219*
1220 END IF
1221*
1222* End of SLATMR
1223*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
real function slangb(norm, n, kl, ku, ab, ldab, work)
SLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition slangb.f:124
real function slange(norm, m, n, a, lda, work)
SLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition slange.f:114
real function slansb(norm, uplo, n, k, ab, ldab, work)
SLANSB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansb.f:129
real function slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:122
real function slansp(norm, uplo, n, ap, work)
SLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansp.f:114
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine slatm1(mode, cond, irsign, idist, iseed, d, n, info)
SLATM1
Definition slatm1.f:135
real function slatm2(m, n, i, j, kl, ku, idist, iseed, d, igrade, dl, dr, ipvtng, iwork, sparse)
SLATM2
Definition slatm2.f:208
real function slatm3(m, n, i, j, isub, jsub, kl, ku, idist, iseed, d, igrade, dl, dr, ipvtng, iwork, sparse)
SLATM3
Definition slatm3.f:226
Here is the call graph for this function:
Here is the caller graph for this function: