LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
slatme.f
Go to the documentation of this file.
1*> \brief \b SLATME
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE SLATME( N, DIST, ISEED, D, MODE, COND, DMAX, EI,
12* RSIGN,
13* UPPER, SIM, DS, MODES, CONDS, KL, KU, ANORM,
14* A,
15* LDA, WORK, INFO )
16*
17* .. Scalar Arguments ..
18* CHARACTER DIST, RSIGN, SIM, UPPER
19* INTEGER INFO, KL, KU, LDA, MODE, MODES, N
20* REAL ANORM, COND, CONDS, DMAX
21* ..
22* .. Array Arguments ..
23* CHARACTER EI( * )
24* INTEGER ISEED( 4 )
25* REAL A( LDA, * ), D( * ), DS( * ), WORK( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> SLATME generates random non-symmetric square matrices with
35*> specified eigenvalues for testing LAPACK programs.
36*>
37*> SLATME operates by applying the following sequence of
38*> operations:
39*>
40*> 1. Set the diagonal to D, where D may be input or
41*> computed according to MODE, COND, DMAX, and RSIGN
42*> as described below.
43*>
44*> 2. If complex conjugate pairs are desired (MODE=0 and EI(1)='R',
45*> or MODE=5), certain pairs of adjacent elements of D are
46*> interpreted as the real and complex parts of a complex
47*> conjugate pair; A thus becomes block diagonal, with 1x1
48*> and 2x2 blocks.
49*>
50*> 3. If UPPER='T', the upper triangle of A is set to random values
51*> out of distribution DIST.
52*>
53*> 4. If SIM='T', A is multiplied on the left by a random matrix
54*> X, whose singular values are specified by DS, MODES, and
55*> CONDS, and on the right by X inverse.
56*>
57*> 5. If KL < N-1, the lower bandwidth is reduced to KL using
58*> Householder transformations. If KU < N-1, the upper
59*> bandwidth is reduced to KU.
60*>
61*> 6. If ANORM is not negative, the matrix is scaled to have
62*> maximum-element-norm ANORM.
63*>
64*> (Note: since the matrix cannot be reduced beyond Hessenberg form,
65*> no packing options are available.)
66*> \endverbatim
67*
68* Arguments:
69* ==========
70*
71*> \param[in] N
72*> \verbatim
73*> N is INTEGER
74*> The number of columns (or rows) of A. Not modified.
75*> \endverbatim
76*>
77*> \param[in] DIST
78*> \verbatim
79*> DIST is CHARACTER*1
80*> On entry, DIST specifies the type of distribution to be used
81*> to generate the random eigen-/singular values, and for the
82*> upper triangle (see UPPER).
83*> 'U' => UNIFORM( 0, 1 ) ( 'U' for uniform )
84*> 'S' => UNIFORM( -1, 1 ) ( 'S' for symmetric )
85*> 'N' => NORMAL( 0, 1 ) ( 'N' for normal )
86*> Not modified.
87*> \endverbatim
88*>
89*> \param[in,out] ISEED
90*> \verbatim
91*> ISEED is INTEGER array, dimension ( 4 )
92*> On entry ISEED specifies the seed of the random number
93*> generator. They should lie between 0 and 4095 inclusive,
94*> and ISEED(4) should be odd. The random number generator
95*> uses a linear congruential sequence limited to small
96*> integers, and so should produce machine independent
97*> random numbers. The values of ISEED are changed on
98*> exit, and can be used in the next call to SLATME
99*> to continue the same random number sequence.
100*> Changed on exit.
101*> \endverbatim
102*>
103*> \param[in,out] D
104*> \verbatim
105*> D is REAL array, dimension ( N )
106*> This array is used to specify the eigenvalues of A. If
107*> MODE=0, then D is assumed to contain the eigenvalues (but
108*> see the description of EI), otherwise they will be
109*> computed according to MODE, COND, DMAX, and RSIGN and
110*> placed in D.
111*> Modified if MODE is nonzero.
112*> \endverbatim
113*>
114*> \param[in] MODE
115*> \verbatim
116*> MODE is INTEGER
117*> On entry this describes how the eigenvalues are to
118*> be specified:
119*> MODE = 0 means use D (with EI) as input
120*> MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
121*> MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
122*> MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
123*> MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
124*> MODE = 5 sets D to random numbers in the range
125*> ( 1/COND , 1 ) such that their logarithms
126*> are uniformly distributed. Each odd-even pair
127*> of elements will be either used as two real
128*> eigenvalues or as the real and imaginary part
129*> of a complex conjugate pair of eigenvalues;
130*> the choice of which is done is random, with
131*> 50-50 probability, for each pair.
132*> MODE = 6 set D to random numbers from same distribution
133*> as the rest of the matrix.
134*> MODE < 0 has the same meaning as ABS(MODE), except that
135*> the order of the elements of D is reversed.
136*> Thus if MODE is between 1 and 4, D has entries ranging
137*> from 1 to 1/COND, if between -1 and -4, D has entries
138*> ranging from 1/COND to 1,
139*> Not modified.
140*> \endverbatim
141*>
142*> \param[in] COND
143*> \verbatim
144*> COND is REAL
145*> On entry, this is used as described under MODE above.
146*> If used, it must be >= 1. Not modified.
147*> \endverbatim
148*>
149*> \param[in] DMAX
150*> \verbatim
151*> DMAX is REAL
152*> If MODE is neither -6, 0 nor 6, the contents of D, as
153*> computed according to MODE and COND, will be scaled by
154*> DMAX / max(abs(D(i))). Note that DMAX need not be
155*> positive: if DMAX is negative (or zero), D will be
156*> scaled by a negative number (or zero).
157*> Not modified.
158*> \endverbatim
159*>
160*> \param[in] EI
161*> \verbatim
162*> EI is CHARACTER*1 array, dimension ( N )
163*> If MODE is 0, and EI(1) is not ' ' (space character),
164*> this array specifies which elements of D (on input) are
165*> real eigenvalues and which are the real and imaginary parts
166*> of a complex conjugate pair of eigenvalues. The elements
167*> of EI may then only have the values 'R' and 'I'. If
168*> EI(j)='R' and EI(j+1)='I', then the j-th eigenvalue is
169*> CMPLX( D(j) , D(j+1) ), and the (j+1)-th is the complex
170*> conjugate thereof. If EI(j)=EI(j+1)='R', then the j-th
171*> eigenvalue is D(j) (i.e., real). EI(1) may not be 'I',
172*> nor may two adjacent elements of EI both have the value 'I'.
173*> If MODE is not 0, then EI is ignored. If MODE is 0 and
174*> EI(1)=' ', then the eigenvalues will all be real.
175*> Not modified.
176*> \endverbatim
177*>
178*> \param[in] RSIGN
179*> \verbatim
180*> RSIGN is CHARACTER*1
181*> If MODE is not 0, 6, or -6, and RSIGN='T', then the
182*> elements of D, as computed according to MODE and COND, will
183*> be multiplied by a random sign (+1 or -1). If RSIGN='F',
184*> they will not be. RSIGN may only have the values 'T' or
185*> 'F'.
186*> Not modified.
187*> \endverbatim
188*>
189*> \param[in] UPPER
190*> \verbatim
191*> UPPER is CHARACTER*1
192*> If UPPER='T', then the elements of A above the diagonal
193*> (and above the 2x2 diagonal blocks, if A has complex
194*> eigenvalues) will be set to random numbers out of DIST.
195*> If UPPER='F', they will not. UPPER may only have the
196*> values 'T' or 'F'.
197*> Not modified.
198*> \endverbatim
199*>
200*> \param[in] SIM
201*> \verbatim
202*> SIM is CHARACTER*1
203*> If SIM='T', then A will be operated on by a "similarity
204*> transform", i.e., multiplied on the left by a matrix X and
205*> on the right by X inverse. X = U S V, where U and V are
206*> random unitary matrices and S is a (diagonal) matrix of
207*> singular values specified by DS, MODES, and CONDS. If
208*> SIM='F', then A will not be transformed.
209*> Not modified.
210*> \endverbatim
211*>
212*> \param[in,out] DS
213*> \verbatim
214*> DS is REAL array, dimension ( N )
215*> This array is used to specify the singular values of X,
216*> in the same way that D specifies the eigenvalues of A.
217*> If MODE=0, the DS contains the singular values, which
218*> may not be zero.
219*> Modified if MODE is nonzero.
220*> \endverbatim
221*>
222*> \param[in] MODES
223*> \verbatim
224*> MODES is INTEGER
225*> \endverbatim
226*>
227*> \param[in] CONDS
228*> \verbatim
229*> CONDS is REAL
230*> Same as MODE and COND, but for specifying the diagonal
231*> of S. MODES=-6 and +6 are not allowed (since they would
232*> result in randomly ill-conditioned eigenvalues.)
233*> \endverbatim
234*>
235*> \param[in] KL
236*> \verbatim
237*> KL is INTEGER
238*> This specifies the lower bandwidth of the matrix. KL=1
239*> specifies upper Hessenberg form. If KL is at least N-1,
240*> then A will have full lower bandwidth. KL must be at
241*> least 1.
242*> Not modified.
243*> \endverbatim
244*>
245*> \param[in] KU
246*> \verbatim
247*> KU is INTEGER
248*> This specifies the upper bandwidth of the matrix. KU=1
249*> specifies lower Hessenberg form. If KU is at least N-1,
250*> then A will have full upper bandwidth; if KU and KL
251*> are both at least N-1, then A will be dense. Only one of
252*> KU and KL may be less than N-1. KU must be at least 1.
253*> Not modified.
254*> \endverbatim
255*>
256*> \param[in] ANORM
257*> \verbatim
258*> ANORM is REAL
259*> If ANORM is not negative, then A will be scaled by a non-
260*> negative real number to make the maximum-element-norm of A
261*> to be ANORM.
262*> Not modified.
263*> \endverbatim
264*>
265*> \param[out] A
266*> \verbatim
267*> A is REAL array, dimension ( LDA, N )
268*> On exit A is the desired test matrix.
269*> Modified.
270*> \endverbatim
271*>
272*> \param[in] LDA
273*> \verbatim
274*> LDA is INTEGER
275*> LDA specifies the first dimension of A as declared in the
276*> calling program. LDA must be at least N.
277*> Not modified.
278*> \endverbatim
279*>
280*> \param[out] WORK
281*> \verbatim
282*> WORK is REAL array, dimension ( 3*N )
283*> Workspace.
284*> Modified.
285*> \endverbatim
286*>
287*> \param[out] INFO
288*> \verbatim
289*> INFO is INTEGER
290*> Error code. On exit, INFO will be set to one of the
291*> following values:
292*> 0 => normal return
293*> -1 => N negative
294*> -2 => DIST illegal string
295*> -5 => MODE not in range -6 to 6
296*> -6 => COND less than 1.0, and MODE neither -6, 0 nor 6
297*> -8 => EI(1) is not ' ' or 'R', EI(j) is not 'R' or 'I', or
298*> two adjacent elements of EI are 'I'.
299*> -9 => RSIGN is not 'T' or 'F'
300*> -10 => UPPER is not 'T' or 'F'
301*> -11 => SIM is not 'T' or 'F'
302*> -12 => MODES=0 and DS has a zero singular value.
303*> -13 => MODES is not in the range -5 to 5.
304*> -14 => MODES is nonzero and CONDS is less than 1.
305*> -15 => KL is less than 1.
306*> -16 => KU is less than 1, or KL and KU are both less than
307*> N-1.
308*> -19 => LDA is less than N.
309*> 1 => Error return from SLATM1 (computing D)
310*> 2 => Cannot scale to DMAX (max. eigenvalue is 0)
311*> 3 => Error return from SLATM1 (computing DS)
312*> 4 => Error return from SLARGE
313*> 5 => Zero singular value from SLATM1.
314*> \endverbatim
315*
316* Authors:
317* ========
318*
319*> \author Univ. of Tennessee
320*> \author Univ. of California Berkeley
321*> \author Univ. of Colorado Denver
322*> \author NAG Ltd.
323*
324*> \ingroup real_matgen
325*
326* =====================================================================
327 SUBROUTINE slatme( N, DIST, ISEED, D, MODE, COND, DMAX, EI,
328 $ RSIGN,
329 $ UPPER, SIM, DS, MODES, CONDS, KL, KU, ANORM,
330 $ A,
331 $ LDA, WORK, INFO )
332*
333* -- LAPACK computational routine --
334* -- LAPACK is a software package provided by Univ. of Tennessee, --
335* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
336*
337* .. Scalar Arguments ..
338 CHARACTER DIST, RSIGN, SIM, UPPER
339 INTEGER INFO, KL, KU, LDA, MODE, MODES, N
340 REAL ANORM, COND, CONDS, DMAX
341* ..
342* .. Array Arguments ..
343 CHARACTER EI( * )
344 INTEGER ISEED( 4 )
345 REAL A( LDA, * ), D( * ), DS( * ), WORK( * )
346* ..
347*
348* =====================================================================
349*
350* .. Parameters ..
351 REAL ZERO
352 PARAMETER ( ZERO = 0.0e0 )
353 REAL ONE
354 PARAMETER ( ONE = 1.0e0 )
355 REAL HALF
356 parameter( half = 1.0e0 / 2.0e0 )
357* ..
358* .. Local Scalars ..
359 LOGICAL BADEI, BADS, USEEI
360 INTEGER I, IC, ICOLS, IDIST, IINFO, IR, IROWS, IRSIGN,
361 $ ISIM, IUPPER, J, JC, JCR, JR
362 REAL ALPHA, TAU, TEMP, XNORMS
363* ..
364* .. Local Arrays ..
365 REAL TEMPA( 1 )
366* ..
367* .. External Functions ..
368 LOGICAL LSAME
369 REAL SLANGE, SLARAN
370 EXTERNAL LSAME, SLANGE, SLARAN
371* ..
372* .. External Subroutines ..
373 EXTERNAL scopy, sgemv, sger, slarfg, slarge, slarnv,
375* ..
376* .. Intrinsic Functions ..
377 INTRINSIC abs, max, mod
378* ..
379* .. Executable Statements ..
380*
381* 1) Decode and Test the input parameters.
382* Initialize flags & seed.
383*
384 info = 0
385*
386* Quick return if possible
387*
388 IF( n.EQ.0 )
389 $ RETURN
390*
391* Decode DIST
392*
393 IF( lsame( dist, 'U' ) ) THEN
394 idist = 1
395 ELSE IF( lsame( dist, 'S' ) ) THEN
396 idist = 2
397 ELSE IF( lsame( dist, 'N' ) ) THEN
398 idist = 3
399 ELSE
400 idist = -1
401 END IF
402*
403* Check EI
404*
405 useei = .true.
406 badei = .false.
407 IF( lsame( ei( 1 ), ' ' ) .OR. mode.NE.0 ) THEN
408 useei = .false.
409 ELSE
410 IF( lsame( ei( 1 ), 'R' ) ) THEN
411 DO 10 j = 2, n
412 IF( lsame( ei( j ), 'I' ) ) THEN
413 IF( lsame( ei( j-1 ), 'I' ) )
414 $ badei = .true.
415 ELSE
416 IF( .NOT.lsame( ei( j ), 'R' ) )
417 $ badei = .true.
418 END IF
419 10 CONTINUE
420 ELSE
421 badei = .true.
422 END IF
423 END IF
424*
425* Decode RSIGN
426*
427 IF( lsame( rsign, 'T' ) ) THEN
428 irsign = 1
429 ELSE IF( lsame( rsign, 'F' ) ) THEN
430 irsign = 0
431 ELSE
432 irsign = -1
433 END IF
434*
435* Decode UPPER
436*
437 IF( lsame( upper, 'T' ) ) THEN
438 iupper = 1
439 ELSE IF( lsame( upper, 'F' ) ) THEN
440 iupper = 0
441 ELSE
442 iupper = -1
443 END IF
444*
445* Decode SIM
446*
447 IF( lsame( sim, 'T' ) ) THEN
448 isim = 1
449 ELSE IF( lsame( sim, 'F' ) ) THEN
450 isim = 0
451 ELSE
452 isim = -1
453 END IF
454*
455* Check DS, if MODES=0 and ISIM=1
456*
457 bads = .false.
458 IF( modes.EQ.0 .AND. isim.EQ.1 ) THEN
459 DO 20 j = 1, n
460 IF( ds( j ).EQ.zero )
461 $ bads = .true.
462 20 CONTINUE
463 END IF
464*
465* Set INFO if an error
466*
467 IF( n.LT.0 ) THEN
468 info = -1
469 ELSE IF( idist.EQ.-1 ) THEN
470 info = -2
471 ELSE IF( abs( mode ).GT.6 ) THEN
472 info = -5
473 ELSE IF( ( mode.NE.0 .AND. abs( mode ).NE.6 ) .AND. cond.LT.one )
474 $ THEN
475 info = -6
476 ELSE IF( badei ) THEN
477 info = -8
478 ELSE IF( irsign.EQ.-1 ) THEN
479 info = -9
480 ELSE IF( iupper.EQ.-1 ) THEN
481 info = -10
482 ELSE IF( isim.EQ.-1 ) THEN
483 info = -11
484 ELSE IF( bads ) THEN
485 info = -12
486 ELSE IF( isim.EQ.1 .AND. abs( modes ).GT.5 ) THEN
487 info = -13
488 ELSE IF( isim.EQ.1 .AND. modes.NE.0 .AND. conds.LT.one ) THEN
489 info = -14
490 ELSE IF( kl.LT.1 ) THEN
491 info = -15
492 ELSE IF( ku.LT.1 .OR. ( ku.LT.n-1 .AND. kl.LT.n-1 ) ) THEN
493 info = -16
494 ELSE IF( lda.LT.max( 1, n ) ) THEN
495 info = -19
496 END IF
497*
498 IF( info.NE.0 ) THEN
499 CALL xerbla( 'SLATME', -info )
500 RETURN
501 END IF
502*
503* Initialize random number generator
504*
505 DO 30 i = 1, 4
506 iseed( i ) = mod( abs( iseed( i ) ), 4096 )
507 30 CONTINUE
508*
509 IF( mod( iseed( 4 ), 2 ).NE.1 )
510 $ iseed( 4 ) = iseed( 4 ) + 1
511*
512* 2) Set up diagonal of A
513*
514* Compute D according to COND and MODE
515*
516 CALL slatm1( mode, cond, irsign, idist, iseed, d, n, iinfo )
517 IF( iinfo.NE.0 ) THEN
518 info = 1
519 RETURN
520 END IF
521 IF( mode.NE.0 .AND. abs( mode ).NE.6 ) THEN
522*
523* Scale by DMAX
524*
525 temp = abs( d( 1 ) )
526 DO 40 i = 2, n
527 temp = max( temp, abs( d( i ) ) )
528 40 CONTINUE
529*
530 IF( temp.GT.zero ) THEN
531 alpha = dmax / temp
532 ELSE IF( dmax.NE.zero ) THEN
533 info = 2
534 RETURN
535 ELSE
536 alpha = zero
537 END IF
538*
539 CALL sscal( n, alpha, d, 1 )
540*
541 END IF
542*
543 CALL slaset( 'Full', n, n, zero, zero, a, lda )
544 CALL scopy( n, d, 1, a, lda+1 )
545*
546* Set up complex conjugate pairs
547*
548 IF( mode.EQ.0 ) THEN
549 IF( useei ) THEN
550 DO 50 j = 2, n
551 IF( lsame( ei( j ), 'I' ) ) THEN
552 a( j-1, j ) = a( j, j )
553 a( j, j-1 ) = -a( j, j )
554 a( j, j ) = a( j-1, j-1 )
555 END IF
556 50 CONTINUE
557 END IF
558*
559 ELSE IF( abs( mode ).EQ.5 ) THEN
560*
561 DO 60 j = 2, n, 2
562 IF( slaran( iseed ).GT.half ) THEN
563 a( j-1, j ) = a( j, j )
564 a( j, j-1 ) = -a( j, j )
565 a( j, j ) = a( j-1, j-1 )
566 END IF
567 60 CONTINUE
568 END IF
569*
570* 3) If UPPER='T', set upper triangle of A to random numbers.
571* (but don't modify the corners of 2x2 blocks.)
572*
573 IF( iupper.NE.0 ) THEN
574 DO 70 jc = 2, n
575 IF( a( jc-1, jc ).NE.zero ) THEN
576 jr = jc - 2
577 ELSE
578 jr = jc - 1
579 END IF
580 CALL slarnv( idist, iseed, jr, a( 1, jc ) )
581 70 CONTINUE
582 END IF
583*
584* 4) If SIM='T', apply similarity transformation.
585*
586* -1
587* Transform is X A X , where X = U S V, thus
588*
589* it is U S V A V' (1/S) U'
590*
591 IF( isim.NE.0 ) THEN
592*
593* Compute S (singular values of the eigenvector matrix)
594* according to CONDS and MODES
595*
596 CALL slatm1( modes, conds, 0, 0, iseed, ds, n, iinfo )
597 IF( iinfo.NE.0 ) THEN
598 info = 3
599 RETURN
600 END IF
601*
602* Multiply by V and V'
603*
604 CALL slarge( n, a, lda, iseed, work, iinfo )
605 IF( iinfo.NE.0 ) THEN
606 info = 4
607 RETURN
608 END IF
609*
610* Multiply by S and (1/S)
611*
612 DO 80 j = 1, n
613 CALL sscal( n, ds( j ), a( j, 1 ), lda )
614 IF( ds( j ).NE.zero ) THEN
615 CALL sscal( n, one / ds( j ), a( 1, j ), 1 )
616 ELSE
617 info = 5
618 RETURN
619 END IF
620 80 CONTINUE
621*
622* Multiply by U and U'
623*
624 CALL slarge( n, a, lda, iseed, work, iinfo )
625 IF( iinfo.NE.0 ) THEN
626 info = 4
627 RETURN
628 END IF
629 END IF
630*
631* 5) Reduce the bandwidth.
632*
633 IF( kl.LT.n-1 ) THEN
634*
635* Reduce bandwidth -- kill column
636*
637 DO 90 jcr = kl + 1, n - 1
638 ic = jcr - kl
639 irows = n + 1 - jcr
640 icols = n + kl - jcr
641*
642 CALL scopy( irows, a( jcr, ic ), 1, work, 1 )
643 xnorms = work( 1 )
644 CALL slarfg( irows, xnorms, work( 2 ), 1, tau )
645 work( 1 ) = one
646*
647 CALL sgemv( 'T', irows, icols, one, a( jcr, ic+1 ), lda,
648 $ work, 1, zero, work( irows+1 ), 1 )
649 CALL sger( irows, icols, -tau, work, 1, work( irows+1 ), 1,
650 $ a( jcr, ic+1 ), lda )
651*
652 CALL sgemv( 'N', n, irows, one, a( 1, jcr ), lda, work, 1,
653 $ zero, work( irows+1 ), 1 )
654 CALL sger( n, irows, -tau, work( irows+1 ), 1, work, 1,
655 $ a( 1, jcr ), lda )
656*
657 a( jcr, ic ) = xnorms
658 CALL slaset( 'Full', irows-1, 1, zero, zero, a( jcr+1, ic ),
659 $ lda )
660 90 CONTINUE
661 ELSE IF( ku.LT.n-1 ) THEN
662*
663* Reduce upper bandwidth -- kill a row at a time.
664*
665 DO 100 jcr = ku + 1, n - 1
666 ir = jcr - ku
667 irows = n + ku - jcr
668 icols = n + 1 - jcr
669*
670 CALL scopy( icols, a( ir, jcr ), lda, work, 1 )
671 xnorms = work( 1 )
672 CALL slarfg( icols, xnorms, work( 2 ), 1, tau )
673 work( 1 ) = one
674*
675 CALL sgemv( 'N', irows, icols, one, a( ir+1, jcr ), lda,
676 $ work, 1, zero, work( icols+1 ), 1 )
677 CALL sger( irows, icols, -tau, work( icols+1 ), 1, work, 1,
678 $ a( ir+1, jcr ), lda )
679*
680 CALL sgemv( 'C', icols, n, one, a( jcr, 1 ), lda, work, 1,
681 $ zero, work( icols+1 ), 1 )
682 CALL sger( icols, n, -tau, work, 1, work( icols+1 ), 1,
683 $ a( jcr, 1 ), lda )
684*
685 a( ir, jcr ) = xnorms
686 CALL slaset( 'Full', 1, icols-1, zero, zero, a( ir, jcr+1 ),
687 $ lda )
688 100 CONTINUE
689 END IF
690*
691* Scale the matrix to have norm ANORM
692*
693 IF( anorm.GE.zero ) THEN
694 temp = slange( 'M', n, n, a, lda, tempa )
695 IF( temp.GT.zero ) THEN
696 alpha = anorm / temp
697 DO 110 j = 1, n
698 CALL sscal( n, alpha, a( 1, j ), 1 )
699 110 CONTINUE
700 END IF
701 END IF
702*
703 RETURN
704*
705* End of SLATME
706*
707 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
SGEMV
Definition sgemv.f:158
subroutine sger(m, n, alpha, x, incx, y, incy, a, lda)
SGER
Definition sger.f:130
subroutine slarfg(n, alpha, x, incx, tau)
SLARFG generates an elementary reflector (Householder matrix).
Definition slarfg.f:106
subroutine slarnv(idist, iseed, n, x)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition slarnv.f:97
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine slarge(n, a, lda, iseed, work, info)
SLARGE
Definition slarge.f:87
subroutine slatm1(mode, cond, irsign, idist, iseed, d, n, info)
SLATM1
Definition slatm1.f:135
subroutine slatme(n, dist, iseed, d, mode, cond, dmax, ei, rsign, upper, sim, ds, modes, conds, kl, ku, anorm, a, lda, work, info)
SLATME
Definition slatme.f:332