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