PLASMA  2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
zlatms.f
Go to the documentation of this file.
1  SUBROUTINE zlatms( M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX,
2  $ kl, ku, pack, a, lda, work, info )
3 *
4 * -- LAPACK test routine (version 3.1) --
5 * Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
6 * November 2006
7 *
8 * .. Scalar Arguments ..
9  CHARACTER dist, pack, sym
10  INTEGER info, kl, ku, lda, m, mode, n
11  DOUBLE PRECISION cond, dmax
12 * ..
13 * .. Array Arguments ..
14  INTEGER iseed( 4 )
15  DOUBLE PRECISION d( * )
16  COMPLEX*16 a( lda, * ), work( * )
17 * ..
18 *
19 * Purpose
20 * =======
21 *
22 * ZLATMS generates random matrices with specified singular values
23 * (or hermitian with specified eigenvalues)
24 * for testing LAPACK programs.
25 *
26 * ZLATMS operates by applying the following sequence of
27 * operations:
28 *
29 * Set the diagonal to D, where D may be input or
30 * computed according to MODE, COND, DMAX, and SYM
31 * as described below.
32 *
33 * Generate a matrix with the appropriate band structure, by one
34 * of two methods:
35 *
36 * Method A:
37 * Generate a dense M x N matrix by multiplying D on the left
38 * and the right by random unitary matrices, then:
39 *
40 * Reduce the bandwidth according to KL and KU, using
41 * Householder transformations.
42 *
43 * Method B:
44 * Convert the bandwidth-0 (i.e., diagonal) matrix to a
45 * bandwidth-1 matrix using Givens rotations, "chasing"
46 * out-of-band elements back, much as in QR; then convert
47 * the bandwidth-1 to a bandwidth-2 matrix, etc. Note
48 * that for reasonably small bandwidths (relative to M and
49 * N) this requires less storage, as a dense matrix is not
50 * generated. Also, for hermitian or symmetric matrices,
51 * only one triangle is generated.
52 *
53 * Method A is chosen if the bandwidth is a large fraction of the
54 * order of the matrix, and LDA is at least M (so a dense
55 * matrix can be stored.) Method B is chosen if the bandwidth
56 * is small (< 1/2 N for hermitian or symmetric, < .3 N+M for
57 * non-symmetric), or LDA is less than M and not less than the
58 * bandwidth.
59 *
60 * Pack the matrix if desired. Options specified by PACK are:
61 * no packing
62 * zero out upper half (if hermitian)
63 * zero out lower half (if hermitian)
64 * store the upper half columnwise (if hermitian or upper
65 * triangular)
66 * store the lower half columnwise (if hermitian or lower
67 * triangular)
68 * store the lower triangle in banded format (if hermitian or
69 * lower triangular)
70 * store the upper triangle in banded format (if hermitian or
71 * upper triangular)
72 * store the entire matrix in banded format
73 * If Method B is chosen, and band format is specified, then the
74 * matrix will be generated in the band format, so no repacking
75 * will be necessary.
76 *
77 * Arguments
78 * =========
79 *
80 * M - INTEGER
81 * The number of rows of A. Not modified.
82 *
83 * N - INTEGER
84 * The number of columns of A. N must equal M if the matrix
85 * is symmetric or hermitian (i.e., if SYM is not 'N')
86 * Not modified.
87 *
88 * DIST - CHARACTER*1
89 * On entry, DIST specifies the type of distribution to be used
90 * to generate the random eigen-/singular values.
91 * 'U' => UNIFORM( 0, 1 ) ( 'U' for uniform )
92 * 'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
93 * 'N' => NORMAL( 0, 1 ) ( 'N' for normal )
94 * Not modified.
95 *
96 * ISEED - INTEGER array, dimension ( 4 )
97 * On entry ISEED specifies the seed of the random number
98 * generator. They should lie between 0 and 4095 inclusive,
99 * and ISEED(4) should be odd. The random number generator
100 * uses a linear congruential sequence limited to small
101 * integers, and so should produce machine independent
102 * random numbers. The values of ISEED are changed on
103 * exit, and can be used in the next call to ZLATMS
104 * to continue the same random number sequence.
105 * Changed on exit.
106 *
107 * SYM - CHARACTER*1
108 * If SYM='H', the generated matrix is hermitian, with
109 * eigenvalues specified by D, COND, MODE, and DMAX; they
110 * may be positive, negative, or zero.
111 * If SYM='P', the generated matrix is hermitian, with
112 * eigenvalues (= singular values) specified by D, COND,
113 * MODE, and DMAX; they will not be negative.
114 * If SYM='N', the generated matrix is nonsymmetric, with
115 * singular values specified by D, COND, MODE, and DMAX;
116 * they will not be negative.
117 * If SYM='S', the generated matrix is (complex) symmetric,
118 * with singular values specified by D, COND, MODE, and
119 * DMAX; they will not be negative.
120 * Not modified.
121 *
122 * D - DOUBLE PRECISION array, dimension ( MIN( M, N ) )
123 * This array is used to specify the singular values or
124 * eigenvalues of A (see SYM, above.) If MODE=0, then D is
125 * assumed to contain the singular/eigenvalues, otherwise
126 * they will be computed according to MODE, COND, and DMAX,
127 * and placed in D.
128 * Modified if MODE is nonzero.
129 *
130 * MODE - INTEGER
131 * On entry this describes how the singular/eigenvalues are to
132 * be specified:
133 * MODE = 0 means use D as input
134 * MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
135 * MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
136 * MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
137 * MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
138 * MODE = 5 sets D to random numbers in the range
139 * ( 1/COND , 1 ) such that their logarithms
140 * are uniformly distributed.
141 * MODE = 6 set D to random numbers from same distribution
142 * as the rest of the matrix.
143 * MODE < 0 has the same meaning as ABS(MODE), except that
144 * the order of the elements of D is reversed.
145 * Thus if MODE is positive, D has entries ranging from
146 * 1 to 1/COND, if negative, from 1/COND to 1,
147 * If SYM='H', and MODE is neither 0, 6, nor -6, then
148 * the elements of D will also be multiplied by a random
149 * sign (i.e., +1 or -1.)
150 * Not modified.
151 *
152 * COND - DOUBLE PRECISION
153 * On entry, this is used as described under MODE above.
154 * If used, it must be >= 1. Not modified.
155 *
156 * DMAX - DOUBLE PRECISION
157 * If MODE is neither -6, 0 nor 6, the contents of D, as
158 * computed according to MODE and COND, will be scaled by
159 * DMAX / max(abs(D(i))); thus, the maximum absolute eigen- or
160 * singular value (which is to say the norm) will be abs(DMAX).
161 * Note that DMAX need not be positive: if DMAX is negative
162 * (or zero), D will be scaled by a negative number (or zero).
163 * Not modified.
164 *
165 * KL - INTEGER
166 * This specifies the lower bandwidth of the matrix. For
167 * example, KL=0 implies upper triangular, KL=1 implies upper
168 * Hessenberg, and KL being at least M-1 means that the matrix
169 * has full lower bandwidth. KL must equal KU if the matrix
170 * is symmetric or hermitian.
171 * Not modified.
172 *
173 * KU - INTEGER
174 * This specifies the upper bandwidth of the matrix. For
175 * example, KU=0 implies lower triangular, KU=1 implies lower
176 * Hessenberg, and KU being at least N-1 means that the matrix
177 * has full upper bandwidth. KL must equal KU if the matrix
178 * is symmetric or hermitian.
179 * Not modified.
180 *
181 * PACK - CHARACTER*1
182 * This specifies packing of matrix as follows:
183 * 'N' => no packing
184 * 'U' => zero out all subdiagonal entries (if symmetric
185 * or hermitian)
186 * 'L' => zero out all superdiagonal entries (if symmetric
187 * or hermitian)
188 * 'C' => store the upper triangle columnwise (only if the
189 * matrix is symmetric, hermitian, or upper triangular)
190 * 'R' => store the lower triangle columnwise (only if the
191 * matrix is symmetric, hermitian, or lower triangular)
192 * 'B' => store the lower triangle in band storage scheme
193 * (only if the matrix is symmetric, hermitian, or
194 * lower triangular)
195 * 'Q' => store the upper triangle in band storage scheme
196 * (only if the matrix is symmetric, hermitian, or
197 * upper triangular)
198 * 'Z' => store the entire matrix in band storage scheme
199 * (pivoting can be provided for by using this
200 * option to store A in the trailing rows of
201 * the allocated storage)
202 *
203 * Using these options, the various LAPACK packed and banded
204 * storage schemes can be obtained:
205 * GB - use 'Z'
206 * PB, SB, HB, or TB - use 'B' or 'Q'
207 * PP, SP, HB, or TP - use 'C' or 'R'
208 *
209 * If two calls to ZLATMS differ only in the PACK parameter,
210 * they will generate mathematically equivalent matrices.
211 * Not modified.
212 *
213 * A - COMPLEX*16 array, dimension ( LDA, N )
214 * On exit A is the desired test matrix. A is first generated
215 * in full (unpacked) form, and then packed, if so specified
216 * by PACK. Thus, the first M elements of the first N
217 * columns will always be modified. If PACK specifies a
218 * packed or banded storage scheme, all LDA elements of the
219 * first N columns will be modified; the elements of the
220 * array which do not correspond to elements of the generated
221 * matrix are set to zero.
222 * Modified.
223 *
224 * LDA - INTEGER
225 * LDA specifies the first dimension of A as declared in the
226 * calling program. If PACK='N', 'U', 'L', 'C', or 'R', then
227 * LDA must be at least M. If PACK='B' or 'Q', then LDA must
228 * be at least MIN( KL, M-1) (which is equal to MIN(KU,N-1)).
229 * If PACK='Z', LDA must be large enough to hold the packed
230 * array: MIN( KU, N-1) + MIN( KL, M-1) + 1.
231 * Not modified.
232 *
233 * WORK - COMPLEX*16 array, dimension ( 3*MAX( N, M ) )
234 * Workspace.
235 * Modified.
236 *
237 * INFO - INTEGER
238 * Error code. On exit, INFO will be set to one of the
239 * following values:
240 * 0 => normal return
241 * -1 => M negative or unequal to N and SYM='S', 'H', or 'P'
242 * -2 => N negative
243 * -3 => DIST illegal string
244 * -5 => SYM illegal string
245 * -7 => MODE not in range -6 to 6
246 * -8 => COND less than 1.0, and MODE neither -6, 0 nor 6
247 * -10 => KL negative
248 * -11 => KU negative, or SYM is not 'N' and KU is not equal to
249 * KL
250 * -12 => PACK illegal string, or PACK='U' or 'L', and SYM='N';
251 * or PACK='C' or 'Q' and SYM='N' and KL is not zero;
252 * or PACK='R' or 'B' and SYM='N' and KU is not zero;
253 * or PACK='U', 'L', 'C', 'R', 'B', or 'Q', and M is not
254 * N.
255 * -14 => LDA is less than M, or PACK='Z' and LDA is less than
256 * MIN(KU,N-1) + MIN(KL,M-1) + 1.
257 * 1 => Error return from DLATM1
258 * 2 => Cannot scale to DMAX (max. sing. value is 0)
259 * 3 => Error return from ZLAGGE, CLAGHE or CLAGSY
260 *
261 * =====================================================================
262 *
263 * .. Parameters ..
264  DOUBLE PRECISION zero
265  parameter( zero = 0.0d+0 )
266  DOUBLE PRECISION one
267  parameter( one = 1.0d+0 )
268  COMPLEX*16 czero
269  parameter( czero = ( 0.0d+0, 0.0d+0 ) )
270  DOUBLE PRECISION twopi
271  parameter( twopi = 6.2831853071795864769252867663d+0 )
272 * ..
273 * .. Local Scalars ..
274  LOGICAL givens, ilextr, iltemp, topdwn, zsym
275  INTEGER i, ic, icol, idist, iendch, iinfo, il, ilda,
276  $ ioffg, ioffst, ipack, ipackg, ir, ir1, ir2,
277  $ irow, irsign, iskew, isym, isympk, j, jc, jch,
278  $ jkl, jku, jr, k, llb, minlda, mnmin, mr, nc,
279  $ uub
280  DOUBLE PRECISION alpha, angle, realc, temp
281  COMPLEX*16 c, ct, ctemp, dummy, extra, s, st
282 * ..
283 * .. External Functions ..
284  LOGICAL lsame
285  DOUBLE PRECISION dlarnd
286  COMPLEX*16 zlarnd
287  EXTERNAL lsame, dlarnd, zlarnd
288 * ..
289 * .. External Subroutines ..
290  EXTERNAL dlatm1, dscal, xerbla, zlagge, zlaghe, zlagsy,
291  $ zlarot, zlartg, zlaset
292 * ..
293 * .. Intrinsic Functions ..
294  INTRINSIC abs, cos, dble, dcmplx, dconjg, max, min, mod,
295  $ sin
296 * ..
297 * .. Executable Statements ..
298 *
299 * 1) Decode and Test the input parameters.
300 * Initialize flags & seed.
301 *
302  info = 0
303 *
304 * Quick return if possible
305 *
306  IF( m.EQ.0 .OR. n.EQ.0 )
307  $ return
308 *
309 * Decode DIST
310 *
311  IF( lsame( dist, 'U' ) ) THEN
312  idist = 1
313  ELSE IF( lsame( dist, 'S' ) ) THEN
314  idist = 2
315  ELSE IF( lsame( dist, 'N' ) ) THEN
316  idist = 3
317  ELSE
318  idist = -1
319  END IF
320 *
321 * Decode SYM
322 *
323  IF( lsame( sym, 'N' ) ) THEN
324  isym = 1
325  irsign = 0
326  zsym = .false.
327  ELSE IF( lsame( sym, 'P' ) ) THEN
328  isym = 2
329  irsign = 0
330  zsym = .false.
331  ELSE IF( lsame( sym, 'S' ) ) THEN
332  isym = 2
333  irsign = 0
334  zsym = .true.
335  ELSE IF( lsame( sym, 'H' ) ) THEN
336  isym = 2
337  irsign = 1
338  zsym = .false.
339  ELSE
340  isym = -1
341  END IF
342 *
343 * Decode PACK
344 *
345  isympk = 0
346  IF( lsame( pack, 'N' ) ) THEN
347  ipack = 0
348  ELSE IF( lsame( pack, 'U' ) ) THEN
349  ipack = 1
350  isympk = 1
351  ELSE IF( lsame( pack, 'L' ) ) THEN
352  ipack = 2
353  isympk = 1
354  ELSE IF( lsame( pack, 'C' ) ) THEN
355  ipack = 3
356  isympk = 2
357  ELSE IF( lsame( pack, 'R' ) ) THEN
358  ipack = 4
359  isympk = 3
360  ELSE IF( lsame( pack, 'B' ) ) THEN
361  ipack = 5
362  isympk = 3
363  ELSE IF( lsame( pack, 'Q' ) ) THEN
364  ipack = 6
365  isympk = 2
366  ELSE IF( lsame( pack, 'Z' ) ) THEN
367  ipack = 7
368  ELSE
369  ipack = -1
370  END IF
371 *
372 * Set certain internal parameters
373 *
374  mnmin = min( m, n )
375  llb = min( kl, m-1 )
376  uub = min( ku, n-1 )
377  mr = min( m, n+llb )
378  nc = min( n, m+uub )
379 *
380  IF( ipack.EQ.5 .OR. ipack.EQ.6 ) THEN
381  minlda = uub + 1
382  ELSE IF( ipack.EQ.7 ) THEN
383  minlda = llb + uub + 1
384  ELSE
385  minlda = m
386  END IF
387 *
388 * Use Givens rotation method if bandwidth small enough,
389 * or if LDA is too small to store the matrix unpacked.
390 *
391  givens = .false.
392  IF( isym.EQ.1 ) THEN
393  IF( dble( llb+uub ).LT.0.3d0*dble( max( 1, mr+nc ) ) )
394  $ givens = .true.
395  ELSE
396  IF( 2*llb.LT.m )
397  $ givens = .true.
398  END IF
399  IF( lda.LT.m .AND. lda.GE.minlda )
400  $ givens = .true.
401 *
402 * Set INFO if an error
403 *
404  IF( m.LT.0 ) THEN
405  info = -1
406  ELSE IF( m.NE.n .AND. isym.NE.1 ) THEN
407  info = -1
408  ELSE IF( n.LT.0 ) THEN
409  info = -2
410  ELSE IF( idist.EQ.-1 ) THEN
411  info = -3
412  ELSE IF( isym.EQ.-1 ) THEN
413  info = -5
414  ELSE IF( abs( mode ).GT.6 ) THEN
415  info = -7
416  ELSE IF( ( mode.NE.0 .AND. abs( mode ).NE.6 ) .AND. cond.LT.one )
417  $ THEN
418  info = -8
419  ELSE IF( kl.LT.0 ) THEN
420  info = -10
421  ELSE IF( ku.LT.0 .OR. ( isym.NE.1 .AND. kl.NE.ku ) ) THEN
422  info = -11
423  ELSE IF( ipack.EQ.-1 .OR. ( isympk.EQ.1 .AND. isym.EQ.1 ) .OR.
424  $ ( isympk.EQ.2 .AND. isym.EQ.1 .AND. kl.GT.0 ) .OR.
425  $ ( isympk.EQ.3 .AND. isym.EQ.1 .AND. ku.GT.0 ) .OR.
426  $ ( isympk.NE.0 .AND. m.NE.n ) ) THEN
427  info = -12
428  ELSE IF( lda.LT.max( 1, minlda ) ) THEN
429  info = -14
430  END IF
431 *
432  IF( info.NE.0 ) THEN
433  CALL xerbla( 'ZLATMS', -info )
434  return
435  END IF
436 *
437 * Initialize random number generator
438 *
439  DO 10 i = 1, 4
440  iseed( i ) = mod( abs( iseed( i ) ), 4096 )
441  10 continue
442 *
443  IF( mod( iseed( 4 ), 2 ).NE.1 )
444  $ iseed( 4 ) = iseed( 4 ) + 1
445 *
446 * 2) Set up D if indicated.
447 *
448 * Compute D according to COND and MODE
449 *
450  CALL dlatm1( mode, cond, irsign, idist, iseed, d, mnmin, iinfo )
451  IF( iinfo.NE.0 ) THEN
452  info = 1
453  return
454  END IF
455 *
456 * Choose Top-Down if D is (apparently) increasing,
457 * Bottom-Up if D is (apparently) decreasing.
458 *
459  IF( abs( d( 1 ) ).LE.abs( d( mnmin ) ) ) THEN
460  topdwn = .true.
461  ELSE
462  topdwn = .false.
463  END IF
464 *
465  IF( mode.NE.0 .AND. abs( mode ).NE.6 ) THEN
466 *
467 * Scale by DMAX
468 *
469  temp = abs( d( 1 ) )
470  DO 20 i = 2, mnmin
471  temp = max( temp, abs( d( i ) ) )
472  20 continue
473 *
474  IF( temp.GT.zero ) THEN
475  alpha = dmax / temp
476  ELSE
477  info = 2
478  return
479  END IF
480 *
481  CALL dscal( mnmin, alpha, d, 1 )
482 *
483  END IF
484 *
485  CALL zlaset( 'Full', lda, n, czero, czero, a, lda )
486 *
487 * 3) Generate Banded Matrix using Givens rotations.
488 * Also the special case of UUB=LLB=0
489 *
490 * Compute Addressing constants to cover all
491 * storage formats. Whether GE, HE, SY, GB, HB, or SB,
492 * upper or lower triangle or both,
493 * the (i,j)-th element is in
494 * A( i - ISKEW*j + IOFFST, j )
495 *
496  IF( ipack.GT.4 ) THEN
497  ilda = lda - 1
498  iskew = 1
499  IF( ipack.GT.5 ) THEN
500  ioffst = uub + 1
501  ELSE
502  ioffst = 1
503  END IF
504  ELSE
505  ilda = lda
506  iskew = 0
507  ioffst = 0
508  END IF
509 *
510 * IPACKG is the format that the matrix is generated in. If this is
511 * different from IPACK, then the matrix must be repacked at the
512 * end. It also signals how to compute the norm, for scaling.
513 *
514  ipackg = 0
515 *
516 * Diagonal Matrix -- We are done, unless it
517 * is to be stored HP/SP/PP/TP (PACK='R' or 'C')
518 *
519  IF( llb.EQ.0 .AND. uub.EQ.0 ) THEN
520  DO 30 j = 1, mnmin
521  a( ( 1-iskew )*j+ioffst, j ) = dcmplx( d( j ) )
522  30 continue
523 *
524  IF( ipack.LE.2 .OR. ipack.GE.5 )
525  $ ipackg = ipack
526 *
527  ELSE IF( givens ) THEN
528 *
529 * Check whether to use Givens rotations,
530 * Householder transformations, or nothing.
531 *
532  IF( isym.EQ.1 ) THEN
533 *
534 * Non-symmetric -- A = U D V
535 *
536  IF( ipack.GT.4 ) THEN
537  ipackg = ipack
538  ELSE
539  ipackg = 0
540  END IF
541 *
542  DO 40 j = 1, mnmin
543  a( ( 1-iskew )*j+ioffst, j ) = dcmplx( d( j ) )
544  40 continue
545 *
546  IF( topdwn ) THEN
547  jkl = 0
548  DO 70 jku = 1, uub
549 *
550 * Transform from bandwidth JKL, JKU-1 to JKL, JKU
551 *
552 * Last row actually rotated is M
553 * Last column actually rotated is MIN( M+JKU, N )
554 *
555  DO 60 jr = 1, min( m+jku, n ) + jkl - 1
556  extra = czero
557  angle = twopi*dlarnd( 1, iseed )
558  c = cos( angle )*zlarnd( 5, iseed )
559  s = sin( angle )*zlarnd( 5, iseed )
560  icol = max( 1, jr-jkl )
561  IF( jr.LT.m ) THEN
562  il = min( n, jr+jku ) + 1 - icol
563  CALL zlarot( .true., jr.GT.jkl, .false., il, c,
564  $ s, a( jr-iskew*icol+ioffst, icol ),
565  $ ilda, extra, dummy )
566  END IF
567 *
568 * Chase "EXTRA" back up
569 *
570  ir = jr
571  ic = icol
572  DO 50 jch = jr - jkl, 1, -jkl - jku
573  IF( ir.LT.m ) THEN
574  CALL zlartg( a( ir+1-iskew*( ic+1 )+ioffst,
575  $ ic+1 ), extra, realc, s, dummy )
576  dummy = zlarnd( 5, iseed )
577  c = dconjg( realc*dummy )
578  s = dconjg( -s*dummy )
579  END IF
580  irow = max( 1, jch-jku )
581  il = ir + 2 - irow
582  ctemp = czero
583  iltemp = jch.GT.jku
584  CALL zlarot( .false., iltemp, .true., il, c, s,
585  $ a( irow-iskew*ic+ioffst, ic ),
586  $ ilda, ctemp, extra )
587  IF( iltemp ) THEN
588  CALL zlartg( a( irow+1-iskew*( ic+1 )+ioffst,
589  $ ic+1 ), ctemp, realc, s, dummy )
590  dummy = zlarnd( 5, iseed )
591  c = dconjg( realc*dummy )
592  s = dconjg( -s*dummy )
593 *
594  icol = max( 1, jch-jku-jkl )
595  il = ic + 2 - icol
596  extra = czero
597  CALL zlarot( .true., jch.GT.jku+jkl, .true.,
598  $ il, c, s, a( irow-iskew*icol+
599  $ ioffst, icol ), ilda, extra,
600  $ ctemp )
601  ic = icol
602  ir = irow
603  END IF
604  50 continue
605  60 continue
606  70 continue
607 *
608  jku = uub
609  DO 100 jkl = 1, llb
610 *
611 * Transform from bandwidth JKL-1, JKU to JKL, JKU
612 *
613  DO 90 jc = 1, min( n+jkl, m ) + jku - 1
614  extra = czero
615  angle = twopi*dlarnd( 1, iseed )
616  c = cos( angle )*zlarnd( 5, iseed )
617  s = sin( angle )*zlarnd( 5, iseed )
618  irow = max( 1, jc-jku )
619  IF( jc.LT.n ) THEN
620  il = min( m, jc+jkl ) + 1 - irow
621  CALL zlarot( .false., jc.GT.jku, .false., il, c,
622  $ s, a( irow-iskew*jc+ioffst, jc ),
623  $ ilda, extra, dummy )
624  END IF
625 *
626 * Chase "EXTRA" back up
627 *
628  ic = jc
629  ir = irow
630  DO 80 jch = jc - jku, 1, -jkl - jku
631  IF( ic.LT.n ) THEN
632  CALL zlartg( a( ir+1-iskew*( ic+1 )+ioffst,
633  $ ic+1 ), extra, realc, s, dummy )
634  dummy = zlarnd( 5, iseed )
635  c = dconjg( realc*dummy )
636  s = dconjg( -s*dummy )
637  END IF
638  icol = max( 1, jch-jkl )
639  il = ic + 2 - icol
640  ctemp = czero
641  iltemp = jch.GT.jkl
642  CALL zlarot( .true., iltemp, .true., il, c, s,
643  $ a( ir-iskew*icol+ioffst, icol ),
644  $ ilda, ctemp, extra )
645  IF( iltemp ) THEN
646  CALL zlartg( a( ir+1-iskew*( icol+1 )+ioffst,
647  $ icol+1 ), ctemp, realc, s,
648  $ dummy )
649  dummy = zlarnd( 5, iseed )
650  c = dconjg( realc*dummy )
651  s = dconjg( -s*dummy )
652  irow = max( 1, jch-jkl-jku )
653  il = ir + 2 - irow
654  extra = czero
655  CALL zlarot( .false., jch.GT.jkl+jku, .true.,
656  $ il, c, s, a( irow-iskew*icol+
657  $ ioffst, icol ), ilda, extra,
658  $ ctemp )
659  ic = icol
660  ir = irow
661  END IF
662  80 continue
663  90 continue
664  100 continue
665 *
666  ELSE
667 *
668 * Bottom-Up -- Start at the bottom right.
669 *
670  jkl = 0
671  DO 130 jku = 1, uub
672 *
673 * Transform from bandwidth JKL, JKU-1 to JKL, JKU
674 *
675 * First row actually rotated is M
676 * First column actually rotated is MIN( M+JKU, N )
677 *
678  iendch = min( m, n+jkl ) - 1
679  DO 120 jc = min( m+jku, n ) - 1, 1 - jkl, -1
680  extra = czero
681  angle = twopi*dlarnd( 1, iseed )
682  c = cos( angle )*zlarnd( 5, iseed )
683  s = sin( angle )*zlarnd( 5, iseed )
684  irow = max( 1, jc-jku+1 )
685  IF( jc.GT.0 ) THEN
686  il = min( m, jc+jkl+1 ) + 1 - irow
687  CALL zlarot( .false., .false., jc+jkl.LT.m, il,
688  $ c, s, a( irow-iskew*jc+ioffst,
689  $ jc ), ilda, dummy, extra )
690  END IF
691 *
692 * Chase "EXTRA" back down
693 *
694  ic = jc
695  DO 110 jch = jc + jkl, iendch, jkl + jku
696  ilextr = ic.GT.0
697  IF( ilextr ) THEN
698  CALL zlartg( a( jch-iskew*ic+ioffst, ic ),
699  $ extra, realc, s, dummy )
700  dummy = zlarnd( 5, iseed )
701  c = realc*dummy
702  s = s*dummy
703  END IF
704  ic = max( 1, ic )
705  icol = min( n-1, jch+jku )
706  iltemp = jch + jku.LT.n
707  ctemp = czero
708  CALL zlarot( .true., ilextr, iltemp, icol+2-ic,
709  $ c, s, a( jch-iskew*ic+ioffst, ic ),
710  $ ilda, extra, ctemp )
711  IF( iltemp ) THEN
712  CALL zlartg( a( jch-iskew*icol+ioffst,
713  $ icol ), ctemp, realc, s, dummy )
714  dummy = zlarnd( 5, iseed )
715  c = realc*dummy
716  s = s*dummy
717  il = min( iendch, jch+jkl+jku ) + 2 - jch
718  extra = czero
719  CALL zlarot( .false., .true.,
720  $ jch+jkl+jku.LE.iendch, il, c, s,
721  $ a( jch-iskew*icol+ioffst,
722  $ icol ), ilda, ctemp, extra )
723  ic = icol
724  END IF
725  110 continue
726  120 continue
727  130 continue
728 *
729  jku = uub
730  DO 160 jkl = 1, llb
731 *
732 * Transform from bandwidth JKL-1, JKU to JKL, JKU
733 *
734 * First row actually rotated is MIN( N+JKL, M )
735 * First column actually rotated is N
736 *
737  iendch = min( n, m+jku ) - 1
738  DO 150 jr = min( n+jkl, m ) - 1, 1 - jku, -1
739  extra = czero
740  angle = twopi*dlarnd( 1, iseed )
741  c = cos( angle )*zlarnd( 5, iseed )
742  s = sin( angle )*zlarnd( 5, iseed )
743  icol = max( 1, jr-jkl+1 )
744  IF( jr.GT.0 ) THEN
745  il = min( n, jr+jku+1 ) + 1 - icol
746  CALL zlarot( .true., .false., jr+jku.LT.n, il,
747  $ c, s, a( jr-iskew*icol+ioffst,
748  $ icol ), ilda, dummy, extra )
749  END IF
750 *
751 * Chase "EXTRA" back down
752 *
753  ir = jr
754  DO 140 jch = jr + jku, iendch, jkl + jku
755  ilextr = ir.GT.0
756  IF( ilextr ) THEN
757  CALL zlartg( a( ir-iskew*jch+ioffst, jch ),
758  $ extra, realc, s, dummy )
759  dummy = zlarnd( 5, iseed )
760  c = realc*dummy
761  s = s*dummy
762  END IF
763  ir = max( 1, ir )
764  irow = min( m-1, jch+jkl )
765  iltemp = jch + jkl.LT.m
766  ctemp = czero
767  CALL zlarot( .false., ilextr, iltemp, irow+2-ir,
768  $ c, s, a( ir-iskew*jch+ioffst,
769  $ jch ), ilda, extra, ctemp )
770  IF( iltemp ) THEN
771  CALL zlartg( a( irow-iskew*jch+ioffst, jch ),
772  $ ctemp, realc, s, dummy )
773  dummy = zlarnd( 5, iseed )
774  c = realc*dummy
775  s = s*dummy
776  il = min( iendch, jch+jkl+jku ) + 2 - jch
777  extra = czero
778  CALL zlarot( .true., .true.,
779  $ jch+jkl+jku.LE.iendch, il, c, s,
780  $ a( irow-iskew*jch+ioffst, jch ),
781  $ ilda, ctemp, extra )
782  ir = irow
783  END IF
784  140 continue
785  150 continue
786  160 continue
787 *
788  END IF
789 *
790  ELSE
791 *
792 * Symmetric -- A = U D U'
793 * Hermitian -- A = U D U*
794 *
795  ipackg = ipack
796  ioffg = ioffst
797 *
798  IF( topdwn ) THEN
799 *
800 * Top-Down -- Generate Upper triangle only
801 *
802  IF( ipack.GE.5 ) THEN
803  ipackg = 6
804  ioffg = uub + 1
805  ELSE
806  ipackg = 1
807  END IF
808 *
809  DO 170 j = 1, mnmin
810  a( ( 1-iskew )*j+ioffg, j ) = dcmplx( d( j ) )
811  170 continue
812 *
813  DO 200 k = 1, uub
814  DO 190 jc = 1, n - 1
815  irow = max( 1, jc-k )
816  il = min( jc+1, k+2 )
817  extra = czero
818  ctemp = a( jc-iskew*( jc+1 )+ioffg, jc+1 )
819  angle = twopi*dlarnd( 1, iseed )
820  c = cos( angle )*zlarnd( 5, iseed )
821  s = sin( angle )*zlarnd( 5, iseed )
822  IF( zsym ) THEN
823  ct = c
824  st = s
825  ELSE
826  ctemp = dconjg( ctemp )
827  ct = dconjg( c )
828  st = dconjg( s )
829  END IF
830  CALL zlarot( .false., jc.GT.k, .true., il, c, s,
831  $ a( irow-iskew*jc+ioffg, jc ), ilda,
832  $ extra, ctemp )
833  CALL zlarot( .true., .true., .false.,
834  $ min( k, n-jc )+1, ct, st,
835  $ a( ( 1-iskew )*jc+ioffg, jc ), ilda,
836  $ ctemp, dummy )
837 *
838 * Chase EXTRA back up the matrix
839 *
840  icol = jc
841  DO 180 jch = jc - k, 1, -k
842  CALL zlartg( a( jch+1-iskew*( icol+1 )+ioffg,
843  $ icol+1 ), extra, realc, s, dummy )
844  dummy = zlarnd( 5, iseed )
845  c = dconjg( realc*dummy )
846  s = dconjg( -s*dummy )
847  ctemp = a( jch-iskew*( jch+1 )+ioffg, jch+1 )
848  IF( zsym ) THEN
849  ct = c
850  st = s
851  ELSE
852  ctemp = dconjg( ctemp )
853  ct = dconjg( c )
854  st = dconjg( s )
855  END IF
856  CALL zlarot( .true., .true., .true., k+2, c, s,
857  $ a( ( 1-iskew )*jch+ioffg, jch ),
858  $ ilda, ctemp, extra )
859  irow = max( 1, jch-k )
860  il = min( jch+1, k+2 )
861  extra = czero
862  CALL zlarot( .false., jch.GT.k, .true., il, ct,
863  $ st, a( irow-iskew*jch+ioffg, jch ),
864  $ ilda, extra, ctemp )
865  icol = jch
866  180 continue
867  190 continue
868  200 continue
869 *
870 * If we need lower triangle, copy from upper. Note that
871 * the order of copying is chosen to work for 'q' -> 'b'
872 *
873  IF( ipack.NE.ipackg .AND. ipack.NE.3 ) THEN
874  DO 230 jc = 1, n
875  irow = ioffst - iskew*jc
876  IF( zsym ) THEN
877  DO 210 jr = jc, min( n, jc+uub )
878  a( jr+irow, jc ) = a( jc-iskew*jr+ioffg, jr )
879  210 continue
880  ELSE
881  DO 220 jr = jc, min( n, jc+uub )
882  a( jr+irow, jc ) = dconjg( a( jc-iskew*jr+
883  $ ioffg, jr ) )
884  220 continue
885  END IF
886  230 continue
887  IF( ipack.EQ.5 ) THEN
888  DO 250 jc = n - uub + 1, n
889  DO 240 jr = n + 2 - jc, uub + 1
890  a( jr, jc ) = czero
891  240 continue
892  250 continue
893  END IF
894  IF( ipackg.EQ.6 ) THEN
895  ipackg = ipack
896  ELSE
897  ipackg = 0
898  END IF
899  END IF
900  ELSE
901 *
902 * Bottom-Up -- Generate Lower triangle only
903 *
904  IF( ipack.GE.5 ) THEN
905  ipackg = 5
906  IF( ipack.EQ.6 )
907  $ ioffg = 1
908  ELSE
909  ipackg = 2
910  END IF
911 *
912  DO 260 j = 1, mnmin
913  a( ( 1-iskew )*j+ioffg, j ) = dcmplx( d( j ) )
914  260 continue
915 *
916  DO 290 k = 1, uub
917  DO 280 jc = n - 1, 1, -1
918  il = min( n+1-jc, k+2 )
919  extra = czero
920  ctemp = a( 1+( 1-iskew )*jc+ioffg, jc )
921  angle = twopi*dlarnd( 1, iseed )
922  c = cos( angle )*zlarnd( 5, iseed )
923  s = sin( angle )*zlarnd( 5, iseed )
924  IF( zsym ) THEN
925  ct = c
926  st = s
927  ELSE
928  ctemp = dconjg( ctemp )
929  ct = dconjg( c )
930  st = dconjg( s )
931  END IF
932  CALL zlarot( .false., .true., n-jc.GT.k, il, c, s,
933  $ a( ( 1-iskew )*jc+ioffg, jc ), ilda,
934  $ ctemp, extra )
935  icol = max( 1, jc-k+1 )
936  CALL zlarot( .true., .false., .true., jc+2-icol,
937  $ ct, st, a( jc-iskew*icol+ioffg,
938  $ icol ), ilda, dummy, ctemp )
939 *
940 * Chase EXTRA back down the matrix
941 *
942  icol = jc
943  DO 270 jch = jc + k, n - 1, k
944  CALL zlartg( a( jch-iskew*icol+ioffg, icol ),
945  $ extra, realc, s, dummy )
946  dummy = zlarnd( 5, iseed )
947  c = realc*dummy
948  s = s*dummy
949  ctemp = a( 1+( 1-iskew )*jch+ioffg, jch )
950  IF( zsym ) THEN
951  ct = c
952  st = s
953  ELSE
954  ctemp = dconjg( ctemp )
955  ct = dconjg( c )
956  st = dconjg( s )
957  END IF
958  CALL zlarot( .true., .true., .true., k+2, c, s,
959  $ a( jch-iskew*icol+ioffg, icol ),
960  $ ilda, extra, ctemp )
961  il = min( n+1-jch, k+2 )
962  extra = czero
963  CALL zlarot( .false., .true., n-jch.GT.k, il,
964  $ ct, st, a( ( 1-iskew )*jch+ioffg,
965  $ jch ), ilda, ctemp, extra )
966  icol = jch
967  270 continue
968  280 continue
969  290 continue
970 *
971 * If we need upper triangle, copy from lower. Note that
972 * the order of copying is chosen to work for 'b' -> 'q'
973 *
974  IF( ipack.NE.ipackg .AND. ipack.NE.4 ) THEN
975  DO 320 jc = n, 1, -1
976  irow = ioffst - iskew*jc
977  IF( zsym ) THEN
978  DO 300 jr = jc, max( 1, jc-uub ), -1
979  a( jr+irow, jc ) = a( jc-iskew*jr+ioffg, jr )
980  300 continue
981  ELSE
982  DO 310 jr = jc, max( 1, jc-uub ), -1
983  a( jr+irow, jc ) = dconjg( a( jc-iskew*jr+
984  $ ioffg, jr ) )
985  310 continue
986  END IF
987  320 continue
988  IF( ipack.EQ.6 ) THEN
989  DO 340 jc = 1, uub
990  DO 330 jr = 1, uub + 1 - jc
991  a( jr, jc ) = czero
992  330 continue
993  340 continue
994  END IF
995  IF( ipackg.EQ.5 ) THEN
996  ipackg = ipack
997  ELSE
998  ipackg = 0
999  END IF
1000  END IF
1001  END IF
1002 *
1003 * Ensure that the diagonal is real if Hermitian
1004 *
1005  IF( .NOT.zsym ) THEN
1006  DO 350 jc = 1, n
1007  irow = ioffst + ( 1-iskew )*jc
1008  a( irow, jc ) = dcmplx( dble( a( irow, jc ) ) )
1009  350 continue
1010  END IF
1011 *
1012  END IF
1013 *
1014  ELSE
1015 *
1016 * 4) Generate Banded Matrix by first
1017 * Rotating by random Unitary matrices,
1018 * then reducing the bandwidth using Householder
1019 * transformations.
1020 *
1021 * Note: we should get here only if LDA .ge. N
1022 *
1023  IF( isym.EQ.1 ) THEN
1024 *
1025 * Non-symmetric -- A = U D V
1026 *
1027  CALL zlagge( mr, nc, llb, uub, d, a, lda, iseed, work,
1028  $ iinfo )
1029  ELSE
1030 *
1031 * Symmetric -- A = U D U' or
1032 * Hermitian -- A = U D U*
1033 *
1034  IF( zsym ) THEN
1035  CALL zlagsy( m, llb, d, a, lda, iseed, work, iinfo )
1036  ELSE
1037  CALL zlaghe( m, llb, d, a, lda, iseed, work, iinfo )
1038  END IF
1039  END IF
1040 *
1041  IF( iinfo.NE.0 ) THEN
1042  info = 3
1043  return
1044  END IF
1045  END IF
1046 *
1047 * 5) Pack the matrix
1048 *
1049  IF( ipack.NE.ipackg ) THEN
1050  IF( ipack.EQ.1 ) THEN
1051 *
1052 * 'U' -- Upper triangular, not packed
1053 *
1054  DO 370 j = 1, m
1055  DO 360 i = j + 1, m
1056  a( i, j ) = czero
1057  360 continue
1058  370 continue
1059 *
1060  ELSE IF( ipack.EQ.2 ) THEN
1061 *
1062 * 'L' -- Lower triangular, not packed
1063 *
1064  DO 390 j = 2, m
1065  DO 380 i = 1, j - 1
1066  a( i, j ) = czero
1067  380 continue
1068  390 continue
1069 *
1070  ELSE IF( ipack.EQ.3 ) THEN
1071 *
1072 * 'C' -- Upper triangle packed Columnwise.
1073 *
1074  icol = 1
1075  irow = 0
1076  DO 410 j = 1, m
1077  DO 400 i = 1, j
1078  irow = irow + 1
1079  IF( irow.GT.lda ) THEN
1080  irow = 1
1081  icol = icol + 1
1082  END IF
1083  a( irow, icol ) = a( i, j )
1084  400 continue
1085  410 continue
1086 *
1087  ELSE IF( ipack.EQ.4 ) THEN
1088 *
1089 * 'R' -- Lower triangle packed Columnwise.
1090 *
1091  icol = 1
1092  irow = 0
1093  DO 430 j = 1, m
1094  DO 420 i = j, m
1095  irow = irow + 1
1096  IF( irow.GT.lda ) THEN
1097  irow = 1
1098  icol = icol + 1
1099  END IF
1100  a( irow, icol ) = a( i, j )
1101  420 continue
1102  430 continue
1103 *
1104  ELSE IF( ipack.GE.5 ) THEN
1105 *
1106 * 'B' -- The lower triangle is packed as a band matrix.
1107 * 'Q' -- The upper triangle is packed as a band matrix.
1108 * 'Z' -- The whole matrix is packed as a band matrix.
1109 *
1110  IF( ipack.EQ.5 )
1111  $ uub = 0
1112  IF( ipack.EQ.6 )
1113  $ llb = 0
1114 *
1115  DO 450 j = 1, uub
1116  DO 440 i = min( j+llb, m ), 1, -1
1117  a( i-j+uub+1, j ) = a( i, j )
1118  440 continue
1119  450 continue
1120 *
1121  DO 470 j = uub + 2, n
1122  DO 460 i = j - uub, min( j+llb, m )
1123  a( i-j+uub+1, j ) = a( i, j )
1124  460 continue
1125  470 continue
1126  END IF
1127 *
1128 * If packed, zero out extraneous elements.
1129 *
1130 * Symmetric/Triangular Packed --
1131 * zero out everything after A(IROW,ICOL)
1132 *
1133  IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1134  DO 490 jc = icol, m
1135  DO 480 jr = irow + 1, lda
1136  a( jr, jc ) = czero
1137  480 continue
1138  irow = 0
1139  490 continue
1140 *
1141  ELSE IF( ipack.GE.5 ) THEN
1142 *
1143 * Packed Band --
1144 * 1st row is now in A( UUB+2-j, j), zero above it
1145 * m-th row is now in A( M+UUB-j,j), zero below it
1146 * last non-zero diagonal is now in A( UUB+LLB+1,j ),
1147 * zero below it, too.
1148 *
1149  ir1 = uub + llb + 2
1150  ir2 = uub + m + 2
1151  DO 520 jc = 1, n
1152  DO 500 jr = 1, uub + 1 - jc
1153  a( jr, jc ) = czero
1154  500 continue
1155  DO 510 jr = max( 1, min( ir1, ir2-jc ) ), lda
1156  a( jr, jc ) = czero
1157  510 continue
1158  520 continue
1159  END IF
1160  END IF
1161 *
1162  return
1163 *
1164 * End of ZLATMS
1165 *
1166  END