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

◆ zlatme()

subroutine zlatme ( integer  n,
character  dist,
integer, dimension( 4 )  iseed,
complex*16, dimension( * )  d,
integer  mode,
double precision  cond,
complex*16  dmax,
character  rsign,
character  upper,
character  sim,
double precision, dimension( * )  ds,
integer  modes,
double precision  conds,
integer  kl,
integer  ku,
double precision  anorm,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( * )  work,
integer  info 
)

ZLATME

Purpose:
    ZLATME generates random non-symmetric square matrices with
    specified eigenvalues for testing LAPACK programs.

    ZLATME operates by applying the following sequence of
    operations:

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

    2. If UPPER='T', the upper triangle of A is set to random values
         out of distribution DIST.

    3. If SIM='T', A is multiplied on the left by a random matrix
         X, whose singular values are specified by DS, MODES, and
         CONDS, and on the right by X inverse.

    4. If KL < N-1, the lower bandwidth is reduced to KL using
         Householder transformations.  If KU < N-1, the upper
         bandwidth is reduced to KU.

    5. If ANORM is not negative, the matrix is scaled to have
         maximum-element-norm ANORM.

    (Note: since the matrix cannot be reduced beyond Hessenberg form,
     no packing options are available.)
Parameters
[in]N
          N is INTEGER
           The number of columns (or rows) of A. Not modified.
[in]DIST
          DIST is CHARACTER*1
           On entry, DIST specifies the type of distribution to be used
           to generate the random eigen-/singular values, and on the
           upper triangle (see UPPER).
           'U' => UNIFORM( 0, 1 )  ( 'U' for uniform )
           'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
           'N' => NORMAL( 0, 1 )   ( 'N' for normal )
           'D' => uniform on the complex disc |z| < 1.
           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 ZLATME
           to continue the same random number sequence.
           Changed on exit.
[in,out]D
          D is COMPLEX*16 array, dimension ( N )
           This array is used to specify the eigenvalues of A.  If
           MODE=0, then D is assumed to contain the eigenvalues
           otherwise they will be computed according to MODE, COND,
           DMAX, and RSIGN and placed in D.
           Modified if MODE is nonzero.
[in]MODE
          MODE is INTEGER
           On entry this describes how the eigenvalues are to
           be specified:
           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 between 1 and 4, D has entries ranging
              from 1 to 1/COND, if between -1 and -4, D has entries
              ranging from 1/COND to 1,
           Not modified.
[in]COND
          COND is DOUBLE PRECISION
           On entry, this is used as described under MODE above.
           If used, it must be >= 1. Not modified.
[in]DMAX
          DMAX is COMPLEX*16
           If MODE is neither -6, 0 nor 6, the contents of D, as
           computed according to MODE and COND, will be scaled by
           DMAX / max(abs(D(i))).  Note that DMAX need not be
           positive or real: if DMAX is negative or complex (or zero),
           D will be scaled by a negative or complex number (or zero).
           If RSIGN='F' then the largest (absolute) eigenvalue will be
           equal to DMAX.
           Not modified.
[in]RSIGN
          RSIGN is CHARACTER*1
           If MODE is not 0, 6, or -6, and RSIGN='T', then the
           elements of D, as computed according to MODE and COND, will
           be multiplied by a random complex number from the unit
           circle |z| = 1.  If RSIGN='F', they will not be.  RSIGN may
           only have the values 'T' or 'F'.
           Not modified.
[in]UPPER
          UPPER is CHARACTER*1
           If UPPER='T', then the elements of A above the diagonal
           will be set to random numbers out of DIST.  If UPPER='F',
           they will not.  UPPER may only have the values 'T' or 'F'.
           Not modified.
[in]SIM
          SIM is CHARACTER*1
           If SIM='T', then A will be operated on by a "similarity
           transform", i.e., multiplied on the left by a matrix X and
           on the right by X inverse.  X = U S V, where U and V are
           random unitary matrices and S is a (diagonal) matrix of
           singular values specified by DS, MODES, and CONDS.  If
           SIM='F', then A will not be transformed.
           Not modified.
[in,out]DS
          DS is DOUBLE PRECISION array, dimension ( N )
           This array is used to specify the singular values of X,
           in the same way that D specifies the eigenvalues of A.
           If MODE=0, the DS contains the singular values, which
           may not be zero.
           Modified if MODE is nonzero.
[in]MODES
          MODES is INTEGER
[in]CONDS
          CONDS is DOUBLE PRECISION
           Similar to MODE and COND, but for specifying the diagonal
           of S.  MODES=-6 and +6 are not allowed (since they would
           result in randomly ill-conditioned eigenvalues.)
[in]KL
          KL is INTEGER
           This specifies the lower bandwidth of the  matrix.  KL=1
           specifies upper Hessenberg form.  If KL is at least N-1,
           then A will have full lower bandwidth.
           Not modified.
[in]KU
          KU is INTEGER
           This specifies the upper bandwidth of the  matrix.  KU=1
           specifies lower Hessenberg form.  If KU is at least N-1,
           then A will have full upper bandwidth; if KU and KL
           are both at least N-1, then A will be dense.  Only one of
           KU and KL may be less than N-1.
           Not modified.
[in]ANORM
          ANORM is DOUBLE PRECISION
           If ANORM is not negative, then A will be scaled by a non-
           negative real number to make the maximum-element-norm of A
           to be ANORM.
           Not modified.
[out]A
          A is COMPLEX*16 array, dimension ( LDA, N )
           On exit A is the desired test matrix.
           Modified.
[in]LDA
          LDA is INTEGER
           LDA specifies the first dimension of A as declared in the
           calling program.  LDA must be at least M.
           Not modified.
[out]WORK
          WORK is COMPLEX*16 array, dimension ( 3*N )
           Workspace.
           Modified.
[out]INFO
          INFO is INTEGER
           Error code.  On exit, INFO will be set to one of the
           following values:
             0 => normal return
            -1 => N negative
            -2 => DIST illegal string
            -5 => MODE not in range -6 to 6
            -6 => COND less than 1.0, and MODE neither -6, 0 nor 6
            -9 => RSIGN is not 'T' or 'F'
           -10 => UPPER is not 'T' or 'F'
           -11 => SIM   is not 'T' or 'F'
           -12 => MODES=0 and DS has a zero singular value.
           -13 => MODES is not in the range -5 to 5.
           -14 => MODES is nonzero and CONDS is less than 1.
           -15 => KL is less than 1.
           -16 => KU is less than 1, or KL and KU are both less than
                  N-1.
           -19 => LDA is less than M.
            1  => Error return from ZLATM1 (computing D)
            2  => Cannot scale to DMAX (max. eigenvalue is 0)
            3  => Error return from DLATM1 (computing DS)
            4  => Error return from ZLARGE
            5  => Zero singular value from DLATM1.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 296 of file zlatme.f.

301*
302* -- LAPACK computational routine --
303* -- LAPACK is a software package provided by Univ. of Tennessee, --
304* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
305*
306* .. Scalar Arguments ..
307 CHARACTER DIST, RSIGN, SIM, UPPER
308 INTEGER INFO, KL, KU, LDA, MODE, MODES, N
309 DOUBLE PRECISION ANORM, COND, CONDS
310 COMPLEX*16 DMAX
311* ..
312* .. Array Arguments ..
313 INTEGER ISEED( 4 )
314 DOUBLE PRECISION DS( * )
315 COMPLEX*16 A( LDA, * ), D( * ), WORK( * )
316* ..
317*
318* =====================================================================
319*
320* .. Parameters ..
321 DOUBLE PRECISION ZERO
322 parameter( zero = 0.0d+0 )
323 DOUBLE PRECISION ONE
324 parameter( one = 1.0d+0 )
325 COMPLEX*16 CZERO
326 parameter( czero = ( 0.0d+0, 0.0d+0 ) )
327 COMPLEX*16 CONE
328 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
329* ..
330* .. Local Scalars ..
331 LOGICAL BADS
332 INTEGER I, IC, ICOLS, IDIST, IINFO, IR, IROWS, IRSIGN,
333 $ ISIM, IUPPER, J, JC, JCR
334 DOUBLE PRECISION RALPHA, TEMP
335 COMPLEX*16 ALPHA, TAU, XNORMS
336* ..
337* .. Local Arrays ..
338 DOUBLE PRECISION TEMPA( 1 )
339* ..
340* .. External Functions ..
341 LOGICAL LSAME
342 DOUBLE PRECISION ZLANGE
343 COMPLEX*16 ZLARND
344 EXTERNAL lsame, zlange, zlarnd
345* ..
346* .. External Subroutines ..
347 EXTERNAL dlatm1, xerbla, zcopy, zdscal, zgemv, zgerc,
349 $ zscal
350* ..
351* .. Intrinsic Functions ..
352 INTRINSIC abs, dconjg, max, mod
353* ..
354* .. Executable Statements ..
355*
356* 1) Decode and Test the input parameters.
357* Initialize flags & seed.
358*
359 info = 0
360*
361* Quick return if possible
362*
363 IF( n.EQ.0 )
364 $ RETURN
365*
366* Decode DIST
367*
368 IF( lsame( dist, 'U' ) ) THEN
369 idist = 1
370 ELSE IF( lsame( dist, 'S' ) ) THEN
371 idist = 2
372 ELSE IF( lsame( dist, 'N' ) ) THEN
373 idist = 3
374 ELSE IF( lsame( dist, 'D' ) ) THEN
375 idist = 4
376 ELSE
377 idist = -1
378 END IF
379*
380* Decode RSIGN
381*
382 IF( lsame( rsign, 'T' ) ) THEN
383 irsign = 1
384 ELSE IF( lsame( rsign, 'F' ) ) THEN
385 irsign = 0
386 ELSE
387 irsign = -1
388 END IF
389*
390* Decode UPPER
391*
392 IF( lsame( upper, 'T' ) ) THEN
393 iupper = 1
394 ELSE IF( lsame( upper, 'F' ) ) THEN
395 iupper = 0
396 ELSE
397 iupper = -1
398 END IF
399*
400* Decode SIM
401*
402 IF( lsame( sim, 'T' ) ) THEN
403 isim = 1
404 ELSE IF( lsame( sim, 'F' ) ) THEN
405 isim = 0
406 ELSE
407 isim = -1
408 END IF
409*
410* Check DS, if MODES=0 and ISIM=1
411*
412 bads = .false.
413 IF( modes.EQ.0 .AND. isim.EQ.1 ) THEN
414 DO 10 j = 1, n
415 IF( ds( j ).EQ.zero )
416 $ bads = .true.
417 10 CONTINUE
418 END IF
419*
420* Set INFO if an error
421*
422 IF( n.LT.0 ) THEN
423 info = -1
424 ELSE IF( idist.EQ.-1 ) THEN
425 info = -2
426 ELSE IF( abs( mode ).GT.6 ) THEN
427 info = -5
428 ELSE IF( ( mode.NE.0 .AND. abs( mode ).NE.6 ) .AND. cond.LT.one )
429 $ THEN
430 info = -6
431 ELSE IF( irsign.EQ.-1 ) THEN
432 info = -9
433 ELSE IF( iupper.EQ.-1 ) THEN
434 info = -10
435 ELSE IF( isim.EQ.-1 ) THEN
436 info = -11
437 ELSE IF( bads ) THEN
438 info = -12
439 ELSE IF( isim.EQ.1 .AND. abs( modes ).GT.5 ) THEN
440 info = -13
441 ELSE IF( isim.EQ.1 .AND. modes.NE.0 .AND. conds.LT.one ) THEN
442 info = -14
443 ELSE IF( kl.LT.1 ) THEN
444 info = -15
445 ELSE IF( ku.LT.1 .OR. ( ku.LT.n-1 .AND. kl.LT.n-1 ) ) THEN
446 info = -16
447 ELSE IF( lda.LT.max( 1, n ) ) THEN
448 info = -19
449 END IF
450*
451 IF( info.NE.0 ) THEN
452 CALL xerbla( 'ZLATME', -info )
453 RETURN
454 END IF
455*
456* Initialize random number generator
457*
458 DO 20 i = 1, 4
459 iseed( i ) = mod( abs( iseed( i ) ), 4096 )
460 20 CONTINUE
461*
462 IF( mod( iseed( 4 ), 2 ).NE.1 )
463 $ iseed( 4 ) = iseed( 4 ) + 1
464*
465* 2) Set up diagonal of A
466*
467* Compute D according to COND and MODE
468*
469 CALL zlatm1( mode, cond, irsign, idist, iseed, d, n, iinfo )
470 IF( iinfo.NE.0 ) THEN
471 info = 1
472 RETURN
473 END IF
474 IF( mode.NE.0 .AND. abs( mode ).NE.6 ) THEN
475*
476* Scale by DMAX
477*
478 temp = abs( d( 1 ) )
479 DO 30 i = 2, n
480 temp = max( temp, abs( d( i ) ) )
481 30 CONTINUE
482*
483 IF( temp.GT.zero ) THEN
484 alpha = dmax / temp
485 ELSE
486 info = 2
487 RETURN
488 END IF
489*
490 CALL zscal( n, alpha, d, 1 )
491*
492 END IF
493*
494 CALL zlaset( 'Full', n, n, czero, czero, a, lda )
495 CALL zcopy( n, d, 1, a, lda+1 )
496*
497* 3) If UPPER='T', set upper triangle of A to random numbers.
498*
499 IF( iupper.NE.0 ) THEN
500 DO 40 jc = 2, n
501 CALL zlarnv( idist, iseed, jc-1, a( 1, jc ) )
502 40 CONTINUE
503 END IF
504*
505* 4) If SIM='T', apply similarity transformation.
506*
507* -1
508* Transform is X A X , where X = U S V, thus
509*
510* it is U S V A V' (1/S) U'
511*
512 IF( isim.NE.0 ) THEN
513*
514* Compute S (singular values of the eigenvector matrix)
515* according to CONDS and MODES
516*
517 CALL dlatm1( modes, conds, 0, 0, iseed, ds, n, iinfo )
518 IF( iinfo.NE.0 ) THEN
519 info = 3
520 RETURN
521 END IF
522*
523* Multiply by V and V'
524*
525 CALL zlarge( n, a, lda, iseed, work, iinfo )
526 IF( iinfo.NE.0 ) THEN
527 info = 4
528 RETURN
529 END IF
530*
531* Multiply by S and (1/S)
532*
533 DO 50 j = 1, n
534 CALL zdscal( n, ds( j ), a( j, 1 ), lda )
535 IF( ds( j ).NE.zero ) THEN
536 CALL zdscal( n, one / ds( j ), a( 1, j ), 1 )
537 ELSE
538 info = 5
539 RETURN
540 END IF
541 50 CONTINUE
542*
543* Multiply by U and U'
544*
545 CALL zlarge( n, a, lda, iseed, work, iinfo )
546 IF( iinfo.NE.0 ) THEN
547 info = 4
548 RETURN
549 END IF
550 END IF
551*
552* 5) Reduce the bandwidth.
553*
554 IF( kl.LT.n-1 ) THEN
555*
556* Reduce bandwidth -- kill column
557*
558 DO 60 jcr = kl + 1, n - 1
559 ic = jcr - kl
560 irows = n + 1 - jcr
561 icols = n + kl - jcr
562*
563 CALL zcopy( irows, a( jcr, ic ), 1, work, 1 )
564 xnorms = work( 1 )
565 CALL zlarfg( irows, xnorms, work( 2 ), 1, tau )
566 tau = dconjg( tau )
567 work( 1 ) = cone
568 alpha = zlarnd( 5, iseed )
569*
570 CALL zgemv( 'C', irows, icols, cone, a( jcr, ic+1 ), lda,
571 $ work, 1, czero, work( irows+1 ), 1 )
572 CALL zgerc( irows, icols, -tau, work, 1, work( irows+1 ), 1,
573 $ a( jcr, ic+1 ), lda )
574*
575 CALL zgemv( 'N', n, irows, cone, a( 1, jcr ), lda, work, 1,
576 $ czero, work( irows+1 ), 1 )
577 CALL zgerc( n, irows, -dconjg( tau ), work( irows+1 ), 1,
578 $ work, 1, a( 1, jcr ), lda )
579*
580 a( jcr, ic ) = xnorms
581 CALL zlaset( 'Full', irows-1, 1, czero, czero,
582 $ a( jcr+1, ic ), lda )
583*
584 CALL zscal( icols+1, alpha, a( jcr, ic ), lda )
585 CALL zscal( n, dconjg( alpha ), a( 1, jcr ), 1 )
586 60 CONTINUE
587 ELSE IF( ku.LT.n-1 ) THEN
588*
589* Reduce upper bandwidth -- kill a row at a time.
590*
591 DO 70 jcr = ku + 1, n - 1
592 ir = jcr - ku
593 irows = n + ku - jcr
594 icols = n + 1 - jcr
595*
596 CALL zcopy( icols, a( ir, jcr ), lda, work, 1 )
597 xnorms = work( 1 )
598 CALL zlarfg( icols, xnorms, work( 2 ), 1, tau )
599 tau = dconjg( tau )
600 work( 1 ) = cone
601 CALL zlacgv( icols-1, work( 2 ), 1 )
602 alpha = zlarnd( 5, iseed )
603*
604 CALL zgemv( 'N', irows, icols, cone, a( ir+1, jcr ), lda,
605 $ work, 1, czero, work( icols+1 ), 1 )
606 CALL zgerc( irows, icols, -tau, work( icols+1 ), 1, work, 1,
607 $ a( ir+1, jcr ), lda )
608*
609 CALL zgemv( 'C', icols, n, cone, a( jcr, 1 ), lda, work, 1,
610 $ czero, work( icols+1 ), 1 )
611 CALL zgerc( icols, n, -dconjg( tau ), work, 1,
612 $ work( icols+1 ), 1, a( jcr, 1 ), lda )
613*
614 a( ir, jcr ) = xnorms
615 CALL zlaset( 'Full', 1, icols-1, czero, czero,
616 $ a( ir, jcr+1 ), lda )
617*
618 CALL zscal( irows+1, alpha, a( ir, jcr ), 1 )
619 CALL zscal( n, dconjg( alpha ), a( jcr, 1 ), lda )
620 70 CONTINUE
621 END IF
622*
623* Scale the matrix to have norm ANORM
624*
625 IF( anorm.GE.zero ) THEN
626 temp = zlange( 'M', n, n, a, lda, tempa )
627 IF( temp.GT.zero ) THEN
628 ralpha = anorm / temp
629 DO 80 j = 1, n
630 CALL zdscal( n, ralpha, a( 1, j ), 1 )
631 80 CONTINUE
632 END IF
633 END IF
634*
635 RETURN
636*
637* End of ZLATME
638*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlatm1(mode, cond, irsign, idist, iseed, d, n, info)
DLATM1
Definition dlatm1.f:135
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:160
subroutine zgerc(m, n, alpha, x, incx, y, incy, a, lda)
ZGERC
Definition zgerc.f:130
subroutine zlacgv(n, x, incx)
ZLACGV conjugates a complex vector.
Definition zlacgv.f:74
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:115
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
Definition zlarfg.f:106
subroutine zlarnv(idist, iseed, n, x)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition zlarnv.f:99
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zscal(n, za, zx, incx)
ZSCAL
Definition zscal.f:78
subroutine zlarge(n, a, lda, iseed, work, info)
ZLARGE
Definition zlarge.f:87
complex *16 function zlarnd(idist, iseed)
ZLARND
Definition zlarnd.f:75
subroutine zlatm1(mode, cond, irsign, idist, iseed, d, n, info)
ZLATM1
Definition zlatm1.f:137
Here is the call graph for this function:
Here is the caller graph for this function: