LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
clatmr.f
Go to the documentation of this file.
1*> \brief \b CLATMR
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 CLATMR( M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX,
12* RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER,
13* CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM,
14* PACK, A, LDA, IWORK, INFO )
15*
16* .. Scalar Arguments ..
17* CHARACTER DIST, GRADE, PACK, PIVTNG, RSIGN, SYM
18* INTEGER INFO, KL, KU, LDA, M, MODE, MODEL, MODER, N
19* REAL ANORM, COND, CONDL, CONDR, SPARSE
20* COMPLEX DMAX
21* ..
22* .. Array Arguments ..
23* INTEGER IPIVOT( * ), ISEED( 4 ), IWORK( * )
24* COMPLEX A( LDA, * ), D( * ), DL( * ), DR( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> CLATMR generates random matrices of various types for testing
34*> LAPACK programs.
35*>
36*> CLATMR operates by applying the following sequence of
37*> operations:
38*>
39*> Generate a matrix A with random entries of distribution DIST
40*> which is symmetric if SYM='S', Hermitian if SYM='H', and
41*> nonsymmetric if SYM='N'.
42*>
43*> Set the diagonal to D, where D may be input or
44*> computed according to MODE, COND, DMAX and RSIGN
45*> as described below.
46*>
47*> Grade the matrix, if desired, from the left and/or right
48*> as specified by GRADE. The inputs DL, MODEL, CONDL, DR,
49*> MODER and CONDR also determine the grading as described
50*> below.
51*>
52*> Permute, if desired, the rows and/or columns as specified by
53*> PIVTNG and IPIVOT.
54*>
55*> Set random entries to zero, if desired, to get a random sparse
56*> matrix as specified by SPARSE.
57*>
58*> Make A a band matrix, if desired, by zeroing out the matrix
59*> outside a band of lower bandwidth KL and upper bandwidth KU.
60*>
61*> Scale A, if desired, to have maximum entry ANORM.
62*>
63*> Pack the matrix if desired. Options specified by PACK are:
64*> no packing
65*> zero out upper half (if symmetric or Hermitian)
66*> zero out lower half (if symmetric or Hermitian)
67*> store the upper half columnwise (if symmetric or Hermitian
68*> or square upper triangular)
69*> store the lower half columnwise (if symmetric or Hermitian
70*> or square lower triangular)
71*> same as upper half rowwise if symmetric
72*> same as conjugate upper half rowwise if Hermitian
73*> store the lower triangle in banded format
74*> (if symmetric or Hermitian)
75*> store the upper triangle in banded format
76*> (if symmetric or Hermitian)
77*> store the entire matrix in banded format
78*>
79*> Note: If two calls to CLATMR differ only in the PACK parameter,
80*> they will generate mathematically equivalent matrices.
81*>
82*> If two calls to CLATMR both have full bandwidth (KL = M-1
83*> and KU = N-1), and differ only in the PIVTNG and PACK
84*> parameters, then the matrices generated will differ only
85*> in the order of the rows and/or columns, and otherwise
86*> contain the same data. This consistency cannot be and
87*> is not maintained with less than full bandwidth.
88*> \endverbatim
89*
90* Arguments:
91* ==========
92*
93*> \param[in] M
94*> \verbatim
95*> M is INTEGER
96*> Number of rows of A. Not modified.
97*> \endverbatim
98*>
99*> \param[in] N
100*> \verbatim
101*> N is INTEGER
102*> Number of columns of A. Not modified.
103*> \endverbatim
104*>
105*> \param[in] DIST
106*> \verbatim
107*> DIST is CHARACTER*1
108*> On entry, DIST specifies the type of distribution to be used
109*> to generate a random matrix .
110*> 'U' => real and imaginary parts are independent
111*> UNIFORM( 0, 1 ) ( 'U' for uniform )
112*> 'S' => real and imaginary parts are independent
113*> UNIFORM( -1, 1 ) ( 'S' for symmetric )
114*> 'N' => real and imaginary parts are independent
115*> NORMAL( 0, 1 ) ( 'N' for normal )
116*> 'D' => uniform on interior of unit disk ( 'D' for disk )
117*> Not modified.
118*> \endverbatim
119*>
120*> \param[in,out] ISEED
121*> \verbatim
122*> ISEED is INTEGER array, dimension (4)
123*> On entry ISEED specifies the seed of the random number
124*> generator. They should lie between 0 and 4095 inclusive,
125*> and ISEED(4) should be odd. The random number generator
126*> uses a linear congruential sequence limited to small
127*> integers, and so should produce machine independent
128*> random numbers. The values of ISEED are changed on
129*> exit, and can be used in the next call to CLATMR
130*> to continue the same random number sequence.
131*> Changed on exit.
132*> \endverbatim
133*>
134*> \param[in] SYM
135*> \verbatim
136*> SYM is CHARACTER*1
137*> If SYM='S', generated matrix is symmetric.
138*> If SYM='H', generated matrix is Hermitian.
139*> If SYM='N', generated matrix is nonsymmetric.
140*> Not modified.
141*> \endverbatim
142*>
143*> \param[in,out] D
144*> \verbatim
145*> D is COMPLEX array, dimension (min(M,N))
146*> On entry this array specifies the diagonal entries
147*> of the diagonal of A. D may either be specified
148*> on entry, or set according to MODE and COND as described
149*> below. If the matrix is Hermitian, the real part of D
150*> will be taken. May be changed on exit if MODE is nonzero.
151*> \endverbatim
152*>
153*> \param[in] MODE
154*> \verbatim
155*> MODE is INTEGER
156*> On entry describes how D is to be used:
157*> MODE = 0 means use D as input
158*> MODE = 1 sets D(1)=1 and D(2:N)=1.0/COND
159*> MODE = 2 sets D(1:N-1)=1 and D(N)=1.0/COND
160*> MODE = 3 sets D(I)=COND**(-(I-1)/(N-1))
161*> MODE = 4 sets D(i)=1 - (i-1)/(N-1)*(1 - 1/COND)
162*> MODE = 5 sets D to random numbers in the range
163*> ( 1/COND , 1 ) such that their logarithms
164*> are uniformly distributed.
165*> MODE = 6 set D to random numbers from same distribution
166*> as the rest of the matrix.
167*> MODE < 0 has the same meaning as ABS(MODE), except that
168*> the order of the elements of D is reversed.
169*> Thus if MODE is positive, D has entries ranging from
170*> 1 to 1/COND, if negative, from 1/COND to 1,
171*> Not modified.
172*> \endverbatim
173*>
174*> \param[in] COND
175*> \verbatim
176*> COND is REAL
177*> On entry, used as described under MODE above.
178*> If used, it must be >= 1. Not modified.
179*> \endverbatim
180*>
181*> \param[in] DMAX
182*> \verbatim
183*> DMAX is COMPLEX
184*> If MODE neither -6, 0 nor 6, the diagonal is scaled by
185*> DMAX / max(abs(D(i))), so that maximum absolute entry
186*> of diagonal is abs(DMAX). If DMAX is complex (or zero),
187*> diagonal will be scaled by a complex number (or zero).
188*> \endverbatim
189*>
190*> \param[in] RSIGN
191*> \verbatim
192*> RSIGN is CHARACTER*1
193*> If MODE neither -6, 0 nor 6, specifies sign of diagonal
194*> as follows:
195*> 'T' => diagonal entries are multiplied by a random complex
196*> number uniformly distributed with absolute value 1
197*> 'F' => diagonal unchanged
198*> Not modified.
199*> \endverbatim
200*>
201*> \param[in] GRADE
202*> \verbatim
203*> GRADE is CHARACTER*1
204*> Specifies grading of matrix as follows:
205*> 'N' => no grading
206*> 'L' => matrix premultiplied by diag( DL )
207*> (only if matrix nonsymmetric)
208*> 'R' => matrix postmultiplied by diag( DR )
209*> (only if matrix nonsymmetric)
210*> 'B' => matrix premultiplied by diag( DL ) and
211*> postmultiplied by diag( DR )
212*> (only if matrix nonsymmetric)
213*> 'H' => matrix premultiplied by diag( DL ) and
214*> postmultiplied by diag( CONJG(DL) )
215*> (only if matrix Hermitian or nonsymmetric)
216*> 'S' => matrix premultiplied by diag( DL ) and
217*> postmultiplied by diag( DL )
218*> (only if matrix symmetric or nonsymmetric)
219*> 'E' => matrix premultiplied by diag( DL ) and
220*> postmultiplied by inv( diag( DL ) )
221*> ( 'S' for similarity )
222*> (only if matrix nonsymmetric)
223*> Note: if GRADE='S', then M must equal N.
224*> Not modified.
225*> \endverbatim
226*>
227*> \param[in,out] DL
228*> \verbatim
229*> DL is COMPLEX array, dimension (M)
230*> If MODEL=0, then on entry this array specifies the diagonal
231*> entries of a diagonal matrix used as described under GRADE
232*> above. If MODEL is not zero, then DL will be set according
233*> to MODEL and CONDL, analogous to the way D is set according
234*> to MODE and COND (except there is no DMAX parameter for DL).
235*> If GRADE='E', then DL cannot have zero entries.
236*> Not referenced if GRADE = 'N' or 'R'. Changed on exit.
237*> \endverbatim
238*>
239*> \param[in] MODEL
240*> \verbatim
241*> MODEL is INTEGER
242*> This specifies how the diagonal array DL is to be computed,
243*> just as MODE specifies how D is to be computed.
244*> Not modified.
245*> \endverbatim
246*>
247*> \param[in] CONDL
248*> \verbatim
249*> CONDL is REAL
250*> When MODEL is not zero, this specifies the condition number
251*> of the computed DL. Not modified.
252*> \endverbatim
253*>
254*> \param[in,out] DR
255*> \verbatim
256*> DR is COMPLEX array, dimension (N)
257*> If MODER=0, then on entry this array specifies the diagonal
258*> entries of a diagonal matrix used as described under GRADE
259*> above. If MODER is not zero, then DR will be set according
260*> to MODER and CONDR, analogous to the way D is set according
261*> to MODE and COND (except there is no DMAX parameter for DR).
262*> Not referenced if GRADE = 'N', 'L', 'H' or 'S'.
263*> Changed on exit.
264*> \endverbatim
265*>
266*> \param[in] MODER
267*> \verbatim
268*> MODER is INTEGER
269*> This specifies how the diagonal array DR is to be computed,
270*> just as MODE specifies how D is to be computed.
271*> Not modified.
272*> \endverbatim
273*>
274*> \param[in] CONDR
275*> \verbatim
276*> CONDR is REAL
277*> When MODER is not zero, this specifies the condition number
278*> of the computed DR. Not modified.
279*> \endverbatim
280*>
281*> \param[in] PIVTNG
282*> \verbatim
283*> PIVTNG is CHARACTER*1
284*> On entry specifies pivoting permutations as follows:
285*> 'N' or ' ' => none.
286*> 'L' => left or row pivoting (matrix must be nonsymmetric).
287*> 'R' => right or column pivoting (matrix must be
288*> nonsymmetric).
289*> 'B' or 'F' => both or full pivoting, i.e., on both sides.
290*> In this case, M must equal N
291*>
292*> If two calls to CLATMR both have full bandwidth (KL = M-1
293*> and KU = N-1), and differ only in the PIVTNG and PACK
294*> parameters, then the matrices generated will differ only
295*> in the order of the rows and/or columns, and otherwise
296*> contain the same data. This consistency cannot be
297*> maintained with less than full bandwidth.
298*> \endverbatim
299*>
300*> \param[in] IPIVOT
301*> \verbatim
302*> IPIVOT is INTEGER array, dimension (N or M)
303*> This array specifies the permutation used. After the
304*> basic matrix is generated, the rows, columns, or both
305*> are permuted. If, say, row pivoting is selected, CLATMR
306*> starts with the *last* row and interchanges the M-th and
307*> IPIVOT(M)-th rows, then moves to the next-to-last row,
308*> interchanging the (M-1)-th and the IPIVOT(M-1)-th rows,
309*> and so on. In terms of "2-cycles", the permutation is
310*> (1 IPIVOT(1)) (2 IPIVOT(2)) ... (M IPIVOT(M))
311*> where the rightmost cycle is applied first. This is the
312*> *inverse* of the effect of pivoting in LINPACK. The idea
313*> is that factoring (with pivoting) an identity matrix
314*> which has been inverse-pivoted in this way should
315*> result in a pivot vector identical to IPIVOT.
316*> Not referenced if PIVTNG = 'N'. Not modified.
317*> \endverbatim
318*>
319*> \param[in] KL
320*> \verbatim
321*> KL is INTEGER
322*> On entry specifies the lower bandwidth of the matrix. For
323*> example, KL=0 implies upper triangular, KL=1 implies upper
324*> Hessenberg, and KL at least M-1 implies the matrix is not
325*> banded. Must equal KU if matrix is symmetric or Hermitian.
326*> Not modified.
327*> \endverbatim
328*>
329*> \param[in] KU
330*> \verbatim
331*> KU is INTEGER
332*> On entry specifies the upper bandwidth of the matrix. For
333*> example, KU=0 implies lower triangular, KU=1 implies lower
334*> Hessenberg, and KU at least N-1 implies the matrix is not
335*> banded. Must equal KL if matrix is symmetric or Hermitian.
336*> Not modified.
337*> \endverbatim
338*>
339*> \param[in] SPARSE
340*> \verbatim
341*> SPARSE is REAL
342*> On entry specifies the sparsity of the matrix if a sparse
343*> matrix is to be generated. SPARSE should lie between
344*> 0 and 1. To generate a sparse matrix, for each matrix entry
345*> a uniform ( 0, 1 ) random number x is generated and
346*> compared to SPARSE; if x is larger the matrix entry
347*> is unchanged and if x is smaller the entry is set
348*> to zero. Thus on the average a fraction SPARSE of the
349*> entries will be set to zero.
350*> Not modified.
351*> \endverbatim
352*>
353*> \param[in] ANORM
354*> \verbatim
355*> ANORM is REAL
356*> On entry specifies maximum entry of output matrix
357*> (output matrix will by multiplied by a constant so that
358*> its largest absolute entry equal ANORM)
359*> if ANORM is nonnegative. If ANORM is negative no scaling
360*> is done. Not modified.
361*> \endverbatim
362*>
363*> \param[in] PACK
364*> \verbatim
365*> PACK is CHARACTER*1
366*> On entry specifies packing of matrix as follows:
367*> 'N' => no packing
368*> 'U' => zero out all subdiagonal entries
369*> (if symmetric or Hermitian)
370*> 'L' => zero out all superdiagonal entries
371*> (if symmetric or Hermitian)
372*> 'C' => store the upper triangle columnwise
373*> (only if matrix symmetric or Hermitian or
374*> square upper triangular)
375*> 'R' => store the lower triangle columnwise
376*> (only if matrix symmetric or Hermitian or
377*> square lower triangular)
378*> (same as upper half rowwise if symmetric)
379*> (same as conjugate upper half rowwise if Hermitian)
380*> 'B' => store the lower triangle in band storage scheme
381*> (only if matrix symmetric or Hermitian)
382*> 'Q' => store the upper triangle in band storage scheme
383*> (only if matrix symmetric or Hermitian)
384*> 'Z' => store the entire matrix in band storage scheme
385*> (pivoting can be provided for by using this
386*> option to store A in the trailing rows of
387*> the allocated storage)
388*>
389*> Using these options, the various LAPACK packed and banded
390*> storage schemes can be obtained:
391*> GB - use 'Z'
392*> PB, HB or TB - use 'B' or 'Q'
393*> PP, HP or TP - use 'C' or 'R'
394*>
395*> If two calls to CLATMR differ only in the PACK parameter,
396*> they will generate mathematically equivalent matrices.
397*> Not modified.
398*> \endverbatim
399*>
400*> \param[in,out] A
401*> \verbatim
402*> A is COMPLEX array, dimension (LDA,N)
403*> On exit A is the desired test matrix. Only those
404*> entries of A which are significant on output
405*> will be referenced (even if A is in packed or band
406*> storage format). The 'unoccupied corners' of A in
407*> band format will be zeroed out.
408*> \endverbatim
409*>
410*> \param[in] LDA
411*> \verbatim
412*> LDA is INTEGER
413*> on entry LDA specifies the first dimension of A as
414*> declared in the calling program.
415*> If PACK='N', 'U' or 'L', LDA must be at least max ( 1, M ).
416*> If PACK='C' or 'R', LDA must be at least 1.
417*> If PACK='B', or 'Q', LDA must be MIN ( KU+1, N )
418*> If PACK='Z', LDA must be at least KUU+KLL+1, where
419*> KUU = MIN ( KU, N-1 ) and KLL = MIN ( KL, M-1 )
420*> Not modified.
421*> \endverbatim
422*>
423*> \param[out] IWORK
424*> \verbatim
425*> IWORK is INTEGER array, dimension (N or M)
426*> Workspace. Not referenced if PIVTNG = 'N'. Changed on exit.
427*> \endverbatim
428*>
429*> \param[out] INFO
430*> \verbatim
431*> INFO is INTEGER
432*> Error parameter on exit:
433*> 0 => normal return
434*> -1 => M negative or unequal to N and SYM='S' or 'H'
435*> -2 => N negative
436*> -3 => DIST illegal string
437*> -5 => SYM illegal string
438*> -7 => MODE not in range -6 to 6
439*> -8 => COND less than 1.0, and MODE neither -6, 0 nor 6
440*> -10 => MODE neither -6, 0 nor 6 and RSIGN illegal string
441*> -11 => GRADE illegal string, or GRADE='E' and
442*> M not equal to N, or GRADE='L', 'R', 'B', 'S' or 'E'
443*> and SYM = 'H', or GRADE='L', 'R', 'B', 'H' or 'E'
444*> and SYM = 'S'
445*> -12 => GRADE = 'E' and DL contains zero
446*> -13 => MODEL not in range -6 to 6 and GRADE= 'L', 'B', 'H',
447*> 'S' or 'E'
448*> -14 => CONDL less than 1.0, GRADE='L', 'B', 'H', 'S' or 'E',
449*> and MODEL neither -6, 0 nor 6
450*> -16 => MODER not in range -6 to 6 and GRADE= 'R' or 'B'
451*> -17 => CONDR less than 1.0, GRADE='R' or 'B', and
452*> MODER neither -6, 0 nor 6
453*> -18 => PIVTNG illegal string, or PIVTNG='B' or 'F' and
454*> M not equal to N, or PIVTNG='L' or 'R' and SYM='S'
455*> or 'H'
456*> -19 => IPIVOT contains out of range number and
457*> PIVTNG not equal to 'N'
458*> -20 => KL negative
459*> -21 => KU negative, or SYM='S' or 'H' and KU not equal to KL
460*> -22 => SPARSE not in range 0. to 1.
461*> -24 => PACK illegal string, or PACK='U', 'L', 'B' or 'Q'
462*> and SYM='N', or PACK='C' and SYM='N' and either KL
463*> not equal to 0 or N not equal to M, or PACK='R' and
464*> SYM='N', and either KU not equal to 0 or N not equal
465*> to M
466*> -26 => LDA too small
467*> 1 => Error return from CLATM1 (computing D)
468*> 2 => Cannot scale diagonal to DMAX (max. entry is 0)
469*> 3 => Error return from CLATM1 (computing DL)
470*> 4 => Error return from CLATM1 (computing DR)
471*> 5 => ANORM is positive, but matrix constructed prior to
472*> attempting to scale it to have norm ANORM, is zero
473*> \endverbatim
474*
475* Authors:
476* ========
477*
478*> \author Univ. of Tennessee
479*> \author Univ. of California Berkeley
480*> \author Univ. of Colorado Denver
481*> \author NAG Ltd.
482*
483*> \ingroup complex_matgen
484*
485* =====================================================================
486 SUBROUTINE clatmr( M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX,
487 $ RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER,
488 $ CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM,
489 $ PACK, A, LDA, IWORK, INFO )
490*
491* -- LAPACK computational routine --
492* -- LAPACK is a software package provided by Univ. of Tennessee, --
493* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
494*
495* .. Scalar Arguments ..
496 CHARACTER DIST, GRADE, PACK, PIVTNG, RSIGN, SYM
497 INTEGER INFO, KL, KU, LDA, M, MODE, MODEL, MODER, N
498 REAL ANORM, COND, CONDL, CONDR, SPARSE
499 COMPLEX DMAX
500* ..
501* .. Array Arguments ..
502 INTEGER IPIVOT( * ), ISEED( 4 ), IWORK( * )
503 COMPLEX A( LDA, * ), D( * ), DL( * ), DR( * )
504* ..
505*
506* =====================================================================
507*
508* .. Parameters ..
509 REAL ZERO
510 PARAMETER ( ZERO = 0.0e0 )
511 REAL ONE
512 parameter( one = 1.0e0 )
513 COMPLEX CONE
514 parameter( cone = ( 1.0e0, 0.0e0 ) )
515 COMPLEX CZERO
516 parameter( czero = ( 0.0e0, 0.0e0 ) )
517* ..
518* .. Local Scalars ..
519 LOGICAL BADPVT, DZERO, FULBND
520 INTEGER I, IDIST, IGRADE, IISUB, IPACK, IPVTNG, IRSIGN,
521 $ ISUB, ISYM, J, JJSUB, JSUB, K, KLL, KUU, MNMIN,
522 $ mnsub, mxsub, npvts
523 REAL ONORM, TEMP
524 COMPLEX CALPHA, CTEMP
525* ..
526* .. Local Arrays ..
527 REAL TEMPA( 1 )
528* ..
529* .. External Functions ..
530 LOGICAL LSAME
531 REAL CLANGB, CLANGE, CLANSB, CLANSP, CLANSY
532 COMPLEX CLATM2, CLATM3
533 EXTERNAL lsame, clangb, clange, clansb, clansp, clansy,
534 $ clatm2, clatm3
535* ..
536* .. External Subroutines ..
537 EXTERNAL clatm1, csscal, xerbla
538* ..
539* .. Intrinsic Functions ..
540 INTRINSIC abs, conjg, max, min, mod, real
541* ..
542* .. Executable Statements ..
543*
544* 1) Decode and Test the input parameters.
545* Initialize flags & seed.
546*
547 info = 0
548*
549* Quick return if possible
550*
551 IF( m.EQ.0 .OR. n.EQ.0 )
552 $ RETURN
553*
554* Decode DIST
555*
556 IF( lsame( dist, 'U' ) ) THEN
557 idist = 1
558 ELSE IF( lsame( dist, 'S' ) ) THEN
559 idist = 2
560 ELSE IF( lsame( dist, 'N' ) ) THEN
561 idist = 3
562 ELSE IF( lsame( dist, 'D' ) ) THEN
563 idist = 4
564 ELSE
565 idist = -1
566 END IF
567*
568* Decode SYM
569*
570 IF( lsame( sym, 'H' ) ) THEN
571 isym = 0
572 ELSE IF( lsame( sym, 'N' ) ) THEN
573 isym = 1
574 ELSE IF( lsame( sym, 'S' ) ) THEN
575 isym = 2
576 ELSE
577 isym = -1
578 END IF
579*
580* Decode RSIGN
581*
582 IF( lsame( rsign, 'F' ) ) THEN
583 irsign = 0
584 ELSE IF( lsame( rsign, 'T' ) ) THEN
585 irsign = 1
586 ELSE
587 irsign = -1
588 END IF
589*
590* Decode PIVTNG
591*
592 IF( lsame( pivtng, 'N' ) ) THEN
593 ipvtng = 0
594 ELSE IF( lsame( pivtng, ' ' ) ) THEN
595 ipvtng = 0
596 ELSE IF( lsame( pivtng, 'L' ) ) THEN
597 ipvtng = 1
598 npvts = m
599 ELSE IF( lsame( pivtng, 'R' ) ) THEN
600 ipvtng = 2
601 npvts = n
602 ELSE IF( lsame( pivtng, 'B' ) ) THEN
603 ipvtng = 3
604 npvts = min( n, m )
605 ELSE IF( lsame( pivtng, 'F' ) ) THEN
606 ipvtng = 3
607 npvts = min( n, m )
608 ELSE
609 ipvtng = -1
610 END IF
611*
612* Decode GRADE
613*
614 IF( lsame( grade, 'N' ) ) THEN
615 igrade = 0
616 ELSE IF( lsame( grade, 'L' ) ) THEN
617 igrade = 1
618 ELSE IF( lsame( grade, 'R' ) ) THEN
619 igrade = 2
620 ELSE IF( lsame( grade, 'B' ) ) THEN
621 igrade = 3
622 ELSE IF( lsame( grade, 'E' ) ) THEN
623 igrade = 4
624 ELSE IF( lsame( grade, 'H' ) ) THEN
625 igrade = 5
626 ELSE IF( lsame( grade, 'S' ) ) THEN
627 igrade = 6
628 ELSE
629 igrade = -1
630 END IF
631*
632* Decode PACK
633*
634 IF( lsame( pack, 'N' ) ) THEN
635 ipack = 0
636 ELSE IF( lsame( pack, 'U' ) ) THEN
637 ipack = 1
638 ELSE IF( lsame( pack, 'L' ) ) THEN
639 ipack = 2
640 ELSE IF( lsame( pack, 'C' ) ) THEN
641 ipack = 3
642 ELSE IF( lsame( pack, 'R' ) ) THEN
643 ipack = 4
644 ELSE IF( lsame( pack, 'B' ) ) THEN
645 ipack = 5
646 ELSE IF( lsame( pack, 'Q' ) ) THEN
647 ipack = 6
648 ELSE IF( lsame( pack, 'Z' ) ) THEN
649 ipack = 7
650 ELSE
651 ipack = -1
652 END IF
653*
654* Set certain internal parameters
655*
656 mnmin = min( m, n )
657 kll = min( kl, m-1 )
658 kuu = min( ku, n-1 )
659*
660* If inv(DL) is used, check to see if DL has a zero entry.
661*
662 dzero = .false.
663 IF( igrade.EQ.4 .AND. model.EQ.0 ) THEN
664 DO 10 i = 1, m
665 IF( dl( i ).EQ.czero )
666 $ dzero = .true.
667 10 CONTINUE
668 END IF
669*
670* Check values in IPIVOT
671*
672 badpvt = .false.
673 IF( ipvtng.GT.0 ) THEN
674 DO 20 j = 1, npvts
675 IF( ipivot( j ).LE.0 .OR. ipivot( j ).GT.npvts )
676 $ badpvt = .true.
677 20 CONTINUE
678 END IF
679*
680* Set INFO if an error
681*
682 IF( m.LT.0 ) THEN
683 info = -1
684 ELSE IF( m.NE.n .AND. ( isym.EQ.0 .OR. isym.EQ.2 ) ) THEN
685 info = -1
686 ELSE IF( n.LT.0 ) THEN
687 info = -2
688 ELSE IF( idist.EQ.-1 ) THEN
689 info = -3
690 ELSE IF( isym.EQ.-1 ) THEN
691 info = -5
692 ELSE IF( mode.LT.-6 .OR. mode.GT.6 ) THEN
693 info = -7
694 ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
695 $ cond.LT.one ) THEN
696 info = -8
697 ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
698 $ irsign.EQ.-1 ) THEN
699 info = -10
700 ELSE IF( igrade.EQ.-1 .OR. ( igrade.EQ.4 .AND. m.NE.n ) .OR.
701 $ ( ( igrade.EQ.1 .OR. igrade.EQ.2 .OR. igrade.EQ.3 .OR.
702 $ igrade.EQ.4 .OR. igrade.EQ.6 ) .AND. isym.EQ.0 ) .OR.
703 $ ( ( igrade.EQ.1 .OR. igrade.EQ.2 .OR. igrade.EQ.3 .OR.
704 $ igrade.EQ.4 .OR. igrade.EQ.5 ) .AND. isym.EQ.2 ) ) THEN
705 info = -11
706 ELSE IF( igrade.EQ.4 .AND. dzero ) THEN
707 info = -12
708 ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
709 $ igrade.EQ.5 .OR. igrade.EQ.6 ) .AND.
710 $ ( model.LT.-6 .OR. model.GT.6 ) ) THEN
711 info = -13
712 ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
713 $ igrade.EQ.5 .OR. igrade.EQ.6 ) .AND.
714 $ ( model.NE.-6 .AND. model.NE.0 .AND. model.NE.6 ) .AND.
715 $ condl.LT.one ) THEN
716 info = -14
717 ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
718 $ ( moder.LT.-6 .OR. moder.GT.6 ) ) THEN
719 info = -16
720 ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
721 $ ( moder.NE.-6 .AND. moder.NE.0 .AND. moder.NE.6 ) .AND.
722 $ condr.LT.one ) THEN
723 info = -17
724 ELSE IF( ipvtng.EQ.-1 .OR. ( ipvtng.EQ.3 .AND. m.NE.n ) .OR.
725 $ ( ( ipvtng.EQ.1 .OR. ipvtng.EQ.2 ) .AND. ( isym.EQ.0 .OR.
726 $ isym.EQ.2 ) ) ) THEN
727 info = -18
728 ELSE IF( ipvtng.NE.0 .AND. badpvt ) THEN
729 info = -19
730 ELSE IF( kl.LT.0 ) THEN
731 info = -20
732 ELSE IF( ku.LT.0 .OR. ( ( isym.EQ.0 .OR. isym.EQ.2 ) .AND. kl.NE.
733 $ ku ) ) THEN
734 info = -21
735 ELSE IF( sparse.LT.zero .OR. sparse.GT.one ) THEN
736 info = -22
737 ELSE IF( ipack.EQ.-1 .OR. ( ( ipack.EQ.1 .OR. ipack.EQ.2 .OR.
738 $ ipack.EQ.5 .OR. ipack.EQ.6 ) .AND. isym.EQ.1 ) .OR.
739 $ ( ipack.EQ.3 .AND. isym.EQ.1 .AND. ( kl.NE.0 .OR. m.NE.
740 $ n ) ) .OR. ( ipack.EQ.4 .AND. isym.EQ.1 .AND. ( ku.NE.
741 $ 0 .OR. m.NE.n ) ) ) THEN
742 info = -24
743 ELSE IF( ( ( ipack.EQ.0 .OR. ipack.EQ.1 .OR. ipack.EQ.2 ) .AND.
744 $ lda.LT.max( 1, m ) ) .OR. ( ( ipack.EQ.3 .OR. ipack.EQ.
745 $ 4 ) .AND. lda.LT.1 ) .OR. ( ( ipack.EQ.5 .OR. ipack.EQ.
746 $ 6 ) .AND. lda.LT.kuu+1 ) .OR.
747 $ ( ipack.EQ.7 .AND. lda.LT.kll+kuu+1 ) ) THEN
748 info = -26
749 END IF
750*
751 IF( info.NE.0 ) THEN
752 CALL xerbla( 'CLATMR', -info )
753 RETURN
754 END IF
755*
756* Decide if we can pivot consistently
757*
758 fulbnd = .false.
759 IF( kuu.EQ.n-1 .AND. kll.EQ.m-1 )
760 $ fulbnd = .true.
761*
762* Initialize random number generator
763*
764 DO 30 i = 1, 4
765 iseed( i ) = mod( abs( iseed( i ) ), 4096 )
766 30 CONTINUE
767*
768 iseed( 4 ) = 2*( iseed( 4 ) / 2 ) + 1
769*
770* 2) Set up D, DL, and DR, if indicated.
771*
772* Compute D according to COND and MODE
773*
774 CALL clatm1( mode, cond, irsign, idist, iseed, d, mnmin, info )
775 IF( info.NE.0 ) THEN
776 info = 1
777 RETURN
778 END IF
779 IF( mode.NE.0 .AND. mode.NE.-6 .AND. mode.NE.6 ) THEN
780*
781* Scale by DMAX
782*
783 temp = abs( d( 1 ) )
784 DO 40 i = 2, mnmin
785 temp = max( temp, abs( d( i ) ) )
786 40 CONTINUE
787 IF( temp.EQ.zero .AND. dmax.NE.czero ) THEN
788 info = 2
789 RETURN
790 END IF
791 IF( temp.NE.zero ) THEN
792 calpha = dmax / temp
793 ELSE
794 calpha = cone
795 END IF
796 DO 50 i = 1, mnmin
797 d( i ) = calpha*d( i )
798 50 CONTINUE
799*
800 END IF
801*
802* If matrix Hermitian, make D real
803*
804 IF( isym.EQ.0 ) THEN
805 DO 60 i = 1, mnmin
806 d( i ) = real( d( i ) )
807 60 CONTINUE
808 END IF
809*
810* Compute DL if grading set
811*
812 IF( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR. igrade.EQ.
813 $ 5 .OR. igrade.EQ.6 ) THEN
814 CALL clatm1( model, condl, 0, idist, iseed, dl, m, info )
815 IF( info.NE.0 ) THEN
816 info = 3
817 RETURN
818 END IF
819 END IF
820*
821* Compute DR if grading set
822*
823 IF( igrade.EQ.2 .OR. igrade.EQ.3 ) THEN
824 CALL clatm1( moder, condr, 0, idist, iseed, dr, n, info )
825 IF( info.NE.0 ) THEN
826 info = 4
827 RETURN
828 END IF
829 END IF
830*
831* 3) Generate IWORK if pivoting
832*
833 IF( ipvtng.GT.0 ) THEN
834 DO 70 i = 1, npvts
835 iwork( i ) = i
836 70 CONTINUE
837 IF( fulbnd ) THEN
838 DO 80 i = 1, npvts
839 k = ipivot( i )
840 j = iwork( i )
841 iwork( i ) = iwork( k )
842 iwork( k ) = j
843 80 CONTINUE
844 ELSE
845 DO 90 i = npvts, 1, -1
846 k = ipivot( i )
847 j = iwork( i )
848 iwork( i ) = iwork( k )
849 iwork( k ) = j
850 90 CONTINUE
851 END IF
852 END IF
853*
854* 4) Generate matrices for each kind of PACKing
855* Always sweep matrix columnwise (if symmetric, upper
856* half only) so that matrix generated does not depend
857* on PACK
858*
859 IF( fulbnd ) THEN
860*
861* Use CLATM3 so matrices generated with differing PIVOTing only
862* differ only in the order of their rows and/or columns.
863*
864 IF( ipack.EQ.0 ) THEN
865 IF( isym.EQ.0 ) THEN
866 DO 110 j = 1, n
867 DO 100 i = 1, j
868 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku,
869 $ idist, iseed, d, igrade, dl, dr, ipvtng,
870 $ iwork, sparse )
871 a( isub, jsub ) = ctemp
872 a( jsub, isub ) = conjg( ctemp )
873 100 CONTINUE
874 110 CONTINUE
875 ELSE IF( isym.EQ.1 ) THEN
876 DO 130 j = 1, n
877 DO 120 i = 1, m
878 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku,
879 $ idist, iseed, d, igrade, dl, dr, ipvtng,
880 $ iwork, sparse )
881 a( isub, jsub ) = ctemp
882 120 CONTINUE
883 130 CONTINUE
884 ELSE IF( isym.EQ.2 ) THEN
885 DO 150 j = 1, n
886 DO 140 i = 1, j
887 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku,
888 $ idist, iseed, d, igrade, dl, dr, ipvtng,
889 $ iwork, sparse )
890 a( isub, jsub ) = ctemp
891 a( jsub, isub ) = ctemp
892 140 CONTINUE
893 150 CONTINUE
894 END IF
895*
896 ELSE IF( ipack.EQ.1 ) THEN
897*
898 DO 170 j = 1, n
899 DO 160 i = 1, j
900 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku, idist,
901 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
902 $ sparse )
903 mnsub = min( isub, jsub )
904 mxsub = max( isub, jsub )
905 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
906 a( mnsub, mxsub ) = conjg( ctemp )
907 ELSE
908 a( mnsub, mxsub ) = ctemp
909 END IF
910 IF( mnsub.NE.mxsub )
911 $ a( mxsub, mnsub ) = czero
912 160 CONTINUE
913 170 CONTINUE
914*
915 ELSE IF( ipack.EQ.2 ) THEN
916*
917 DO 190 j = 1, n
918 DO 180 i = 1, j
919 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku, idist,
920 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
921 $ sparse )
922 mnsub = min( isub, jsub )
923 mxsub = max( isub, jsub )
924 IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
925 a( mxsub, mnsub ) = conjg( ctemp )
926 ELSE
927 a( mxsub, mnsub ) = ctemp
928 END IF
929 IF( mnsub.NE.mxsub )
930 $ a( mnsub, mxsub ) = czero
931 180 CONTINUE
932 190 CONTINUE
933*
934 ELSE IF( ipack.EQ.3 ) THEN
935*
936 DO 210 j = 1, n
937 DO 200 i = 1, j
938 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku, idist,
939 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
940 $ sparse )
941*
942* Compute K = location of (ISUB,JSUB) entry in packed
943* array
944*
945 mnsub = min( isub, jsub )
946 mxsub = max( isub, jsub )
947 k = mxsub*( mxsub-1 ) / 2 + mnsub
948*
949* Convert K to (IISUB,JJSUB) location
950*
951 jjsub = ( k-1 ) / lda + 1
952 iisub = k - lda*( jjsub-1 )
953*
954 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
955 a( iisub, jjsub ) = conjg( ctemp )
956 ELSE
957 a( iisub, jjsub ) = ctemp
958 END IF
959 200 CONTINUE
960 210 CONTINUE
961*
962 ELSE IF( ipack.EQ.4 ) THEN
963*
964 DO 230 j = 1, n
965 DO 220 i = 1, j
966 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku, idist,
967 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
968 $ sparse )
969*
970* Compute K = location of (I,J) entry in packed array
971*
972 mnsub = min( isub, jsub )
973 mxsub = max( isub, jsub )
974 IF( mnsub.EQ.1 ) THEN
975 k = mxsub
976 ELSE
977 k = n*( n+1 ) / 2 - ( n-mnsub+1 )*( n-mnsub+2 ) /
978 $ 2 + mxsub - mnsub + 1
979 END IF
980*
981* Convert K to (IISUB,JJSUB) location
982*
983 jjsub = ( k-1 ) / lda + 1
984 iisub = k - lda*( jjsub-1 )
985*
986 IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
987 a( iisub, jjsub ) = conjg( ctemp )
988 ELSE
989 a( iisub, jjsub ) = ctemp
990 END IF
991 220 CONTINUE
992 230 CONTINUE
993*
994 ELSE IF( ipack.EQ.5 ) THEN
995*
996 DO 250 j = 1, n
997 DO 240 i = j - kuu, j
998 IF( i.LT.1 ) THEN
999 a( j-i+1, i+n ) = czero
1000 ELSE
1001 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku,
1002 $ idist, iseed, d, igrade, dl, dr, ipvtng,
1003 $ iwork, sparse )
1004 mnsub = min( isub, jsub )
1005 mxsub = max( isub, jsub )
1006 IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
1007 a( mxsub-mnsub+1, mnsub ) = conjg( ctemp )
1008 ELSE
1009 a( mxsub-mnsub+1, mnsub ) = ctemp
1010 END IF
1011 END IF
1012 240 CONTINUE
1013 250 CONTINUE
1014*
1015 ELSE IF( ipack.EQ.6 ) THEN
1016*
1017 DO 270 j = 1, n
1018 DO 260 i = j - kuu, j
1019 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku, idist,
1020 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
1021 $ sparse )
1022 mnsub = min( isub, jsub )
1023 mxsub = max( isub, jsub )
1024 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
1025 a( mnsub-mxsub+kuu+1, mxsub ) = conjg( ctemp )
1026 ELSE
1027 a( mnsub-mxsub+kuu+1, mxsub ) = ctemp
1028 END IF
1029 260 CONTINUE
1030 270 CONTINUE
1031*
1032 ELSE IF( ipack.EQ.7 ) THEN
1033*
1034 IF( isym.NE.1 ) THEN
1035 DO 290 j = 1, n
1036 DO 280 i = j - kuu, j
1037 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku,
1038 $ idist, iseed, d, igrade, dl, dr, ipvtng,
1039 $ iwork, sparse )
1040 mnsub = min( isub, jsub )
1041 mxsub = max( isub, jsub )
1042 IF( i.LT.1 )
1043 $ a( j-i+1+kuu, i+n ) = czero
1044 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
1045 a( mnsub-mxsub+kuu+1, mxsub ) = conjg( ctemp )
1046 ELSE
1047 a( mnsub-mxsub+kuu+1, mxsub ) = ctemp
1048 END IF
1049 IF( i.GE.1 .AND. mnsub.NE.mxsub ) THEN
1050 IF( mnsub.EQ.isub .AND. isym.EQ.0 ) THEN
1051 a( mxsub-mnsub+1+kuu,
1052 $ mnsub ) = conjg( ctemp )
1053 ELSE
1054 a( mxsub-mnsub+1+kuu, mnsub ) = ctemp
1055 END IF
1056 END IF
1057 280 CONTINUE
1058 290 CONTINUE
1059 ELSE IF( isym.EQ.1 ) THEN
1060 DO 310 j = 1, n
1061 DO 300 i = j - kuu, j + kll
1062 ctemp = clatm3( m, n, i, j, isub, jsub, kl, ku,
1063 $ idist, iseed, d, igrade, dl, dr, ipvtng,
1064 $ iwork, sparse )
1065 a( isub-jsub+kuu+1, jsub ) = ctemp
1066 300 CONTINUE
1067 310 CONTINUE
1068 END IF
1069*
1070 END IF
1071*
1072 ELSE
1073*
1074* Use CLATM2
1075*
1076 IF( ipack.EQ.0 ) THEN
1077 IF( isym.EQ.0 ) THEN
1078 DO 330 j = 1, n
1079 DO 320 i = 1, j
1080 a( i, j ) = clatm2( m, n, i, j, kl, ku, idist,
1081 $ iseed, d, igrade, dl, dr, ipvtng,
1082 $ iwork, sparse )
1083 a( j, i ) = conjg( a( i, j ) )
1084 320 CONTINUE
1085 330 CONTINUE
1086 ELSE IF( isym.EQ.1 ) THEN
1087 DO 350 j = 1, n
1088 DO 340 i = 1, m
1089 a( i, j ) = clatm2( m, n, i, j, kl, ku, idist,
1090 $ iseed, d, igrade, dl, dr, ipvtng,
1091 $ iwork, sparse )
1092 340 CONTINUE
1093 350 CONTINUE
1094 ELSE IF( isym.EQ.2 ) THEN
1095 DO 370 j = 1, n
1096 DO 360 i = 1, j
1097 a( i, j ) = clatm2( m, n, i, j, kl, ku, idist,
1098 $ iseed, d, igrade, dl, dr, ipvtng,
1099 $ iwork, sparse )
1100 a( j, i ) = a( i, j )
1101 360 CONTINUE
1102 370 CONTINUE
1103 END IF
1104*
1105 ELSE IF( ipack.EQ.1 ) THEN
1106*
1107 DO 390 j = 1, n
1108 DO 380 i = 1, j
1109 a( i, j ) = clatm2( m, n, i, j, kl, ku, idist, iseed,
1110 $ d, igrade, dl, dr, ipvtng, iwork, sparse )
1111 IF( i.NE.j )
1112 $ a( j, i ) = czero
1113 380 CONTINUE
1114 390 CONTINUE
1115*
1116 ELSE IF( ipack.EQ.2 ) THEN
1117*
1118 DO 410 j = 1, n
1119 DO 400 i = 1, j
1120 IF( isym.EQ.0 ) THEN
1121 a( j, i ) = conjg( clatm2( m, n, i, j, kl, ku,
1122 $ idist, iseed, d, igrade, dl, dr,
1123 $ ipvtng, iwork, sparse ) )
1124 ELSE
1125 a( j, i ) = clatm2( m, n, i, j, kl, ku, idist,
1126 $ iseed, d, igrade, dl, dr, ipvtng,
1127 $ iwork, sparse )
1128 END IF
1129 IF( i.NE.j )
1130 $ a( i, j ) = czero
1131 400 CONTINUE
1132 410 CONTINUE
1133*
1134 ELSE IF( ipack.EQ.3 ) THEN
1135*
1136 isub = 0
1137 jsub = 1
1138 DO 430 j = 1, n
1139 DO 420 i = 1, j
1140 isub = isub + 1
1141 IF( isub.GT.lda ) THEN
1142 isub = 1
1143 jsub = jsub + 1
1144 END IF
1145 a( isub, jsub ) = clatm2( m, n, i, j, kl, ku, idist,
1146 $ iseed, d, igrade, dl, dr, ipvtng,
1147 $ iwork, sparse )
1148 420 CONTINUE
1149 430 CONTINUE
1150*
1151 ELSE IF( ipack.EQ.4 ) THEN
1152*
1153 IF( isym.EQ.0 .OR. isym.EQ.2 ) THEN
1154 DO 450 j = 1, n
1155 DO 440 i = 1, j
1156*
1157* Compute K = location of (I,J) entry in packed array
1158*
1159 IF( i.EQ.1 ) THEN
1160 k = j
1161 ELSE
1162 k = n*( n+1 ) / 2 - ( n-i+1 )*( n-i+2 ) / 2 +
1163 $ j - i + 1
1164 END IF
1165*
1166* Convert K to (ISUB,JSUB) location
1167*
1168 jsub = ( k-1 ) / lda + 1
1169 isub = k - lda*( jsub-1 )
1170*
1171 a( isub, jsub ) = clatm2( m, n, i, j, kl, ku,
1172 $ idist, iseed, d, igrade, dl, dr,
1173 $ ipvtng, iwork, sparse )
1174 IF( isym.EQ.0 )
1175 $ a( isub, jsub ) = conjg( a( isub, jsub ) )
1176 440 CONTINUE
1177 450 CONTINUE
1178 ELSE
1179 isub = 0
1180 jsub = 1
1181 DO 470 j = 1, n
1182 DO 460 i = j, m
1183 isub = isub + 1
1184 IF( isub.GT.lda ) THEN
1185 isub = 1
1186 jsub = jsub + 1
1187 END IF
1188 a( isub, jsub ) = clatm2( m, n, i, j, kl, ku,
1189 $ idist, iseed, d, igrade, dl, dr,
1190 $ ipvtng, iwork, sparse )
1191 460 CONTINUE
1192 470 CONTINUE
1193 END IF
1194*
1195 ELSE IF( ipack.EQ.5 ) THEN
1196*
1197 DO 490 j = 1, n
1198 DO 480 i = j - kuu, j
1199 IF( i.LT.1 ) THEN
1200 a( j-i+1, i+n ) = czero
1201 ELSE
1202 IF( isym.EQ.0 ) THEN
1203 a( j-i+1, i ) = conjg( clatm2( m, n, i, j, kl,
1204 $ ku, idist, iseed, d, igrade, dl,
1205 $ dr, ipvtng, iwork, sparse ) )
1206 ELSE
1207 a( j-i+1, i ) = clatm2( m, n, i, j, kl, ku,
1208 $ idist, iseed, d, igrade, dl, dr,
1209 $ ipvtng, iwork, sparse )
1210 END IF
1211 END IF
1212 480 CONTINUE
1213 490 CONTINUE
1214*
1215 ELSE IF( ipack.EQ.6 ) THEN
1216*
1217 DO 510 j = 1, n
1218 DO 500 i = j - kuu, j
1219 a( i-j+kuu+1, j ) = clatm2( m, n, i, j, kl, ku, idist,
1220 $ iseed, d, igrade, dl, dr, ipvtng,
1221 $ iwork, sparse )
1222 500 CONTINUE
1223 510 CONTINUE
1224*
1225 ELSE IF( ipack.EQ.7 ) THEN
1226*
1227 IF( isym.NE.1 ) THEN
1228 DO 530 j = 1, n
1229 DO 520 i = j - kuu, j
1230 a( i-j+kuu+1, j ) = clatm2( m, n, i, j, kl, ku,
1231 $ idist, iseed, d, igrade, dl,
1232 $ dr, ipvtng, iwork, sparse )
1233 IF( i.LT.1 )
1234 $ a( j-i+1+kuu, i+n ) = czero
1235 IF( i.GE.1 .AND. i.NE.j ) THEN
1236 IF( isym.EQ.0 ) THEN
1237 a( j-i+1+kuu, i ) = conjg( a( i-j+kuu+1,
1238 $ j ) )
1239 ELSE
1240 a( j-i+1+kuu, i ) = a( i-j+kuu+1, j )
1241 END IF
1242 END IF
1243 520 CONTINUE
1244 530 CONTINUE
1245 ELSE IF( isym.EQ.1 ) THEN
1246 DO 550 j = 1, n
1247 DO 540 i = j - kuu, j + kll
1248 a( i-j+kuu+1, j ) = clatm2( m, n, i, j, kl, ku,
1249 $ idist, iseed, d, igrade, dl,
1250 $ dr, ipvtng, iwork, sparse )
1251 540 CONTINUE
1252 550 CONTINUE
1253 END IF
1254*
1255 END IF
1256*
1257 END IF
1258*
1259* 5) Scaling the norm
1260*
1261 IF( ipack.EQ.0 ) THEN
1262 onorm = clange( 'M', m, n, a, lda, tempa )
1263 ELSE IF( ipack.EQ.1 ) THEN
1264 onorm = clansy( 'M', 'U', n, a, lda, tempa )
1265 ELSE IF( ipack.EQ.2 ) THEN
1266 onorm = clansy( 'M', 'L', n, a, lda, tempa )
1267 ELSE IF( ipack.EQ.3 ) THEN
1268 onorm = clansp( 'M', 'U', n, a, tempa )
1269 ELSE IF( ipack.EQ.4 ) THEN
1270 onorm = clansp( 'M', 'L', n, a, tempa )
1271 ELSE IF( ipack.EQ.5 ) THEN
1272 onorm = clansb( 'M', 'L', n, kll, a, lda, tempa )
1273 ELSE IF( ipack.EQ.6 ) THEN
1274 onorm = clansb( 'M', 'U', n, kuu, a, lda, tempa )
1275 ELSE IF( ipack.EQ.7 ) THEN
1276 onorm = clangb( 'M', n, kll, kuu, a, lda, tempa )
1277 END IF
1278*
1279 IF( anorm.GE.zero ) THEN
1280*
1281 IF( anorm.GT.zero .AND. onorm.EQ.zero ) THEN
1282*
1283* Desired scaling impossible
1284*
1285 info = 5
1286 RETURN
1287*
1288 ELSE IF( ( anorm.GT.one .AND. onorm.LT.one ) .OR.
1289 $ ( anorm.LT.one .AND. onorm.GT.one ) ) THEN
1290*
1291* Scale carefully to avoid over / underflow
1292*
1293 IF( ipack.LE.2 ) THEN
1294 DO 560 j = 1, n
1295 CALL csscal( m, one / onorm, a( 1, j ), 1 )
1296 CALL csscal( m, anorm, a( 1, j ), 1 )
1297 560 CONTINUE
1298*
1299 ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1300*
1301 CALL csscal( n*( n+1 ) / 2, one / onorm, a, 1 )
1302 CALL csscal( n*( n+1 ) / 2, anorm, a, 1 )
1303*
1304 ELSE IF( ipack.GE.5 ) THEN
1305*
1306 DO 570 j = 1, n
1307 CALL csscal( kll+kuu+1, one / onorm, a( 1, j ), 1 )
1308 CALL csscal( kll+kuu+1, anorm, a( 1, j ), 1 )
1309 570 CONTINUE
1310*
1311 END IF
1312*
1313 ELSE
1314*
1315* Scale straightforwardly
1316*
1317 IF( ipack.LE.2 ) THEN
1318 DO 580 j = 1, n
1319 CALL csscal( m, anorm / onorm, a( 1, j ), 1 )
1320 580 CONTINUE
1321*
1322 ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1323*
1324 CALL csscal( n*( n+1 ) / 2, anorm / onorm, a, 1 )
1325*
1326 ELSE IF( ipack.GE.5 ) THEN
1327*
1328 DO 590 j = 1, n
1329 CALL csscal( kll+kuu+1, anorm / onorm, a( 1, j ), 1 )
1330 590 CONTINUE
1331 END IF
1332*
1333 END IF
1334*
1335 END IF
1336*
1337* End of CLATMR
1338*
1339 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine clatm1(mode, cond, irsign, idist, iseed, d, n, info)
CLATM1
Definition clatm1.f:137
subroutine clatmr(m, n, dist, iseed, sym, d, mode, cond, dmax, rsign, grade, dl, model, condl, dr, moder, condr, pivtng, ipivot, kl, ku, sparse, anorm, pack, a, lda, iwork, info)
CLATMR
Definition clatmr.f:490
subroutine csscal(n, sa, cx, incx)
CSSCAL
Definition csscal.f:78