LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zlatmr.f
Go to the documentation of this file.
1*> \brief \b ZLATMR
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 ZLATMR( 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* DOUBLE PRECISION ANORM, COND, CONDL, CONDR, SPARSE
20* COMPLEX*16 DMAX
21* ..
22* .. Array Arguments ..
23* INTEGER IPIVOT( * ), ISEED( 4 ), IWORK( * )
24* COMPLEX*16 A( LDA, * ), D( * ), DL( * ), DR( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> ZLATMR generates random matrices of various types for testing
34*> LAPACK programs.
35*>
36*> ZLATMR 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 ZLATMR differ only in the PACK parameter,
80*> they will generate mathematically equivalent matrices.
81*>
82*> If two calls to ZLATMR 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 ZLATMR
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*16 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 DOUBLE PRECISION
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*16
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*16 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 DOUBLE PRECISION
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*16 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 DOUBLE PRECISION
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 ZLATMR 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, ZLATMR
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 DOUBLE PRECISION
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 DOUBLE PRECISION
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 ZLATMR 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*16 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 ZLATM1 (computing D)
468*> 2 => Cannot scale diagonal to DMAX (max. entry is 0)
469*> 3 => Error return from ZLATM1 (computing DL)
470*> 4 => Error return from ZLATM1 (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 complex16_matgen
484*
485* =====================================================================
486 SUBROUTINE zlatmr( 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 DOUBLE PRECISION ANORM, COND, CONDL, CONDR, SPARSE
499 COMPLEX*16 DMAX
500* ..
501* .. Array Arguments ..
502 INTEGER IPIVOT( * ), ISEED( 4 ), IWORK( * )
503 COMPLEX*16 A( LDA, * ), D( * ), DL( * ), DR( * )
504* ..
505*
506* =====================================================================
507*
508* .. Parameters ..
509 DOUBLE PRECISION ZERO
510 PARAMETER ( ZERO = 0.0d0 )
511 DOUBLE PRECISION ONE
512 parameter( one = 1.0d0 )
513 COMPLEX*16 CONE
514 parameter( cone = ( 1.0d0, 0.0d0 ) )
515 COMPLEX*16 CZERO
516 parameter( czero = ( 0.0d0, 0.0d0 ) )
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 DOUBLE PRECISION ONORM, TEMP
524 COMPLEX*16 CALPHA, CTEMP
525* ..
526* .. Local Arrays ..
527 DOUBLE PRECISION TEMPA( 1 )
528* ..
529* .. External Functions ..
530 LOGICAL LSAME
531 DOUBLE PRECISION ZLANGB, ZLANGE, ZLANSB,
532 $ ZLANSP, ZLANSY
533 COMPLEX*16 ZLATM2, ZLATM3
534 EXTERNAL lsame, zlangb, zlange,
535 $ zlansb, zlansp, zlansy,
536 $ zlatm2, zlatm3
537* ..
538* .. External Subroutines ..
539 EXTERNAL xerbla, zdscal, zlatm1
540* ..
541* .. Intrinsic Functions ..
542 INTRINSIC abs, dble, dconjg, max, min, mod
543* ..
544* .. Executable Statements ..
545*
546* 1) Decode and Test the input parameters.
547* Initialize flags & seed.
548*
549 info = 0
550*
551* Quick return if possible
552*
553 IF( m.EQ.0 .OR. n.EQ.0 )
554 $ RETURN
555*
556* Decode DIST
557*
558 IF( lsame( dist, 'U' ) ) THEN
559 idist = 1
560 ELSE IF( lsame( dist, 'S' ) ) THEN
561 idist = 2
562 ELSE IF( lsame( dist, 'N' ) ) THEN
563 idist = 3
564 ELSE IF( lsame( dist, 'D' ) ) THEN
565 idist = 4
566 ELSE
567 idist = -1
568 END IF
569*
570* Decode SYM
571*
572 IF( lsame( sym, 'H' ) ) THEN
573 isym = 0
574 ELSE IF( lsame( sym, 'N' ) ) THEN
575 isym = 1
576 ELSE IF( lsame( sym, 'S' ) ) THEN
577 isym = 2
578 ELSE
579 isym = -1
580 END IF
581*
582* Decode RSIGN
583*
584 IF( lsame( rsign, 'F' ) ) THEN
585 irsign = 0
586 ELSE IF( lsame( rsign, 'T' ) ) THEN
587 irsign = 1
588 ELSE
589 irsign = -1
590 END IF
591*
592* Decode PIVTNG
593*
594 IF( lsame( pivtng, 'N' ) ) THEN
595 ipvtng = 0
596 ELSE IF( lsame( pivtng, ' ' ) ) THEN
597 ipvtng = 0
598 ELSE IF( lsame( pivtng, 'L' ) ) THEN
599 ipvtng = 1
600 npvts = m
601 ELSE IF( lsame( pivtng, 'R' ) ) THEN
602 ipvtng = 2
603 npvts = n
604 ELSE IF( lsame( pivtng, 'B' ) ) THEN
605 ipvtng = 3
606 npvts = min( n, m )
607 ELSE IF( lsame( pivtng, 'F' ) ) THEN
608 ipvtng = 3
609 npvts = min( n, m )
610 ELSE
611 ipvtng = -1
612 END IF
613*
614* Decode GRADE
615*
616 IF( lsame( grade, 'N' ) ) THEN
617 igrade = 0
618 ELSE IF( lsame( grade, 'L' ) ) THEN
619 igrade = 1
620 ELSE IF( lsame( grade, 'R' ) ) THEN
621 igrade = 2
622 ELSE IF( lsame( grade, 'B' ) ) THEN
623 igrade = 3
624 ELSE IF( lsame( grade, 'E' ) ) THEN
625 igrade = 4
626 ELSE IF( lsame( grade, 'H' ) ) THEN
627 igrade = 5
628 ELSE IF( lsame( grade, 'S' ) ) THEN
629 igrade = 6
630 ELSE
631 igrade = -1
632 END IF
633*
634* Decode PACK
635*
636 IF( lsame( pack, 'N' ) ) THEN
637 ipack = 0
638 ELSE IF( lsame( pack, 'U' ) ) THEN
639 ipack = 1
640 ELSE IF( lsame( pack, 'L' ) ) THEN
641 ipack = 2
642 ELSE IF( lsame( pack, 'C' ) ) THEN
643 ipack = 3
644 ELSE IF( lsame( pack, 'R' ) ) THEN
645 ipack = 4
646 ELSE IF( lsame( pack, 'B' ) ) THEN
647 ipack = 5
648 ELSE IF( lsame( pack, 'Q' ) ) THEN
649 ipack = 6
650 ELSE IF( lsame( pack, 'Z' ) ) THEN
651 ipack = 7
652 ELSE
653 ipack = -1
654 END IF
655*
656* Set certain internal parameters
657*
658 mnmin = min( m, n )
659 kll = min( kl, m-1 )
660 kuu = min( ku, n-1 )
661*
662* If inv(DL) is used, check to see if DL has a zero entry.
663*
664 dzero = .false.
665 IF( igrade.EQ.4 .AND. model.EQ.0 ) THEN
666 DO 10 i = 1, m
667 IF( dl( i ).EQ.czero )
668 $ dzero = .true.
669 10 CONTINUE
670 END IF
671*
672* Check values in IPIVOT
673*
674 badpvt = .false.
675 IF( ipvtng.GT.0 ) THEN
676 DO 20 j = 1, npvts
677 IF( ipivot( j ).LE.0 .OR. ipivot( j ).GT.npvts )
678 $ badpvt = .true.
679 20 CONTINUE
680 END IF
681*
682* Set INFO if an error
683*
684 IF( m.LT.0 ) THEN
685 info = -1
686 ELSE IF( m.NE.n .AND. ( isym.EQ.0 .OR. isym.EQ.2 ) ) THEN
687 info = -1
688 ELSE IF( n.LT.0 ) THEN
689 info = -2
690 ELSE IF( idist.EQ.-1 ) THEN
691 info = -3
692 ELSE IF( isym.EQ.-1 ) THEN
693 info = -5
694 ELSE IF( mode.LT.-6 .OR. mode.GT.6 ) THEN
695 info = -7
696 ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
697 $ cond.LT.one ) THEN
698 info = -8
699 ELSE IF( ( mode.NE.-6 .AND. mode.NE.0 .AND. mode.NE.6 ) .AND.
700 $ irsign.EQ.-1 ) THEN
701 info = -10
702 ELSE IF( igrade.EQ.-1 .OR. ( igrade.EQ.4 .AND. m.NE.n ) .OR.
703 $ ( ( igrade.EQ.1 .OR. igrade.EQ.2 .OR. igrade.EQ.3 .OR.
704 $ igrade.EQ.4 .OR. igrade.EQ.6 ) .AND. isym.EQ.0 ) .OR.
705 $ ( ( igrade.EQ.1 .OR. igrade.EQ.2 .OR. igrade.EQ.3 .OR.
706 $ igrade.EQ.4 .OR. igrade.EQ.5 ) .AND. isym.EQ.2 ) ) THEN
707 info = -11
708 ELSE IF( igrade.EQ.4 .AND. dzero ) THEN
709 info = -12
710 ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
711 $ igrade.EQ.5 .OR. igrade.EQ.6 ) .AND.
712 $ ( model.LT.-6 .OR. model.GT.6 ) ) THEN
713 info = -13
714 ELSE IF( ( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR.
715 $ igrade.EQ.5 .OR. igrade.EQ.6 ) .AND.
716 $ ( model.NE.-6 .AND. model.NE.0 .AND. model.NE.6 ) .AND.
717 $ condl.LT.one ) THEN
718 info = -14
719 ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
720 $ ( moder.LT.-6 .OR. moder.GT.6 ) ) THEN
721 info = -16
722 ELSE IF( ( igrade.EQ.2 .OR. igrade.EQ.3 ) .AND.
723 $ ( moder.NE.-6 .AND. moder.NE.0 .AND. moder.NE.6 ) .AND.
724 $ condr.LT.one ) THEN
725 info = -17
726 ELSE IF( ipvtng.EQ.-1 .OR. ( ipvtng.EQ.3 .AND. m.NE.n ) .OR.
727 $ ( ( ipvtng.EQ.1 .OR. ipvtng.EQ.2 ) .AND. ( isym.EQ.0 .OR.
728 $ isym.EQ.2 ) ) ) THEN
729 info = -18
730 ELSE IF( ipvtng.NE.0 .AND. badpvt ) THEN
731 info = -19
732 ELSE IF( kl.LT.0 ) THEN
733 info = -20
734 ELSE IF( ku.LT.0 .OR. ( ( isym.EQ.0 .OR. isym.EQ.2 ) .AND. kl.NE.
735 $ ku ) ) THEN
736 info = -21
737 ELSE IF( sparse.LT.zero .OR. sparse.GT.one ) THEN
738 info = -22
739 ELSE IF( ipack.EQ.-1 .OR. ( ( ipack.EQ.1 .OR. ipack.EQ.2 .OR.
740 $ ipack.EQ.5 .OR. ipack.EQ.6 ) .AND. isym.EQ.1 ) .OR.
741 $ ( ipack.EQ.3 .AND. isym.EQ.1 .AND. ( kl.NE.0 .OR. m.NE.
742 $ n ) ) .OR. ( ipack.EQ.4 .AND. isym.EQ.1 .AND. ( ku.NE.
743 $ 0 .OR. m.NE.n ) ) ) THEN
744 info = -24
745 ELSE IF( ( ( ipack.EQ.0 .OR. ipack.EQ.1 .OR. ipack.EQ.2 ) .AND.
746 $ lda.LT.max( 1, m ) ) .OR. ( ( ipack.EQ.3 .OR. ipack.EQ.
747 $ 4 ) .AND. lda.LT.1 ) .OR. ( ( ipack.EQ.5 .OR. ipack.EQ.
748 $ 6 ) .AND. lda.LT.kuu+1 ) .OR.
749 $ ( ipack.EQ.7 .AND. lda.LT.kll+kuu+1 ) ) THEN
750 info = -26
751 END IF
752*
753 IF( info.NE.0 ) THEN
754 CALL xerbla( 'ZLATMR', -info )
755 RETURN
756 END IF
757*
758* Decide if we can pivot consistently
759*
760 fulbnd = .false.
761 IF( kuu.EQ.n-1 .AND. kll.EQ.m-1 )
762 $ fulbnd = .true.
763*
764* Initialize random number generator
765*
766 DO 30 i = 1, 4
767 iseed( i ) = mod( abs( iseed( i ) ), 4096 )
768 30 CONTINUE
769*
770 iseed( 4 ) = 2*( iseed( 4 ) / 2 ) + 1
771*
772* 2) Set up D, DL, and DR, if indicated.
773*
774* Compute D according to COND and MODE
775*
776 CALL zlatm1( mode, cond, irsign, idist, iseed, d, mnmin, info )
777 IF( info.NE.0 ) THEN
778 info = 1
779 RETURN
780 END IF
781 IF( mode.NE.0 .AND. mode.NE.-6 .AND. mode.NE.6 ) THEN
782*
783* Scale by DMAX
784*
785 temp = abs( d( 1 ) )
786 DO 40 i = 2, mnmin
787 temp = max( temp, abs( d( i ) ) )
788 40 CONTINUE
789 IF( temp.EQ.zero .AND. dmax.NE.czero ) THEN
790 info = 2
791 RETURN
792 END IF
793 IF( temp.NE.zero ) THEN
794 calpha = dmax / temp
795 ELSE
796 calpha = cone
797 END IF
798 DO 50 i = 1, mnmin
799 d( i ) = calpha*d( i )
800 50 CONTINUE
801*
802 END IF
803*
804* If matrix Hermitian, make D real
805*
806 IF( isym.EQ.0 ) THEN
807 DO 60 i = 1, mnmin
808 d( i ) = dble( d( i ) )
809 60 CONTINUE
810 END IF
811*
812* Compute DL if grading set
813*
814 IF( igrade.EQ.1 .OR. igrade.EQ.3 .OR. igrade.EQ.4 .OR. igrade.EQ.
815 $ 5 .OR. igrade.EQ.6 ) THEN
816 CALL zlatm1( model, condl, 0, idist, iseed, dl, m, info )
817 IF( info.NE.0 ) THEN
818 info = 3
819 RETURN
820 END IF
821 END IF
822*
823* Compute DR if grading set
824*
825 IF( igrade.EQ.2 .OR. igrade.EQ.3 ) THEN
826 CALL zlatm1( moder, condr, 0, idist, iseed, dr, n, info )
827 IF( info.NE.0 ) THEN
828 info = 4
829 RETURN
830 END IF
831 END IF
832*
833* 3) Generate IWORK if pivoting
834*
835 IF( ipvtng.GT.0 ) THEN
836 DO 70 i = 1, npvts
837 iwork( i ) = i
838 70 CONTINUE
839 IF( fulbnd ) THEN
840 DO 80 i = 1, npvts
841 k = ipivot( i )
842 j = iwork( i )
843 iwork( i ) = iwork( k )
844 iwork( k ) = j
845 80 CONTINUE
846 ELSE
847 DO 90 i = npvts, 1, -1
848 k = ipivot( i )
849 j = iwork( i )
850 iwork( i ) = iwork( k )
851 iwork( k ) = j
852 90 CONTINUE
853 END IF
854 END IF
855*
856* 4) Generate matrices for each kind of PACKing
857* Always sweep matrix columnwise (if symmetric, upper
858* half only) so that matrix generated does not depend
859* on PACK
860*
861 IF( fulbnd ) THEN
862*
863* Use ZLATM3 so matrices generated with differing PIVOTing only
864* differ only in the order of their rows and/or columns.
865*
866 IF( ipack.EQ.0 ) THEN
867 IF( isym.EQ.0 ) THEN
868 DO 110 j = 1, n
869 DO 100 i = 1, j
870 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
871 $ idist, iseed, d, igrade, dl, dr, ipvtng,
872 $ iwork, sparse )
873 a( isub, jsub ) = ctemp
874 a( jsub, isub ) = dconjg( ctemp )
875 100 CONTINUE
876 110 CONTINUE
877 ELSE IF( isym.EQ.1 ) THEN
878 DO 130 j = 1, n
879 DO 120 i = 1, m
880 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
881 $ idist, iseed, d, igrade, dl, dr, ipvtng,
882 $ iwork, sparse )
883 a( isub, jsub ) = ctemp
884 120 CONTINUE
885 130 CONTINUE
886 ELSE IF( isym.EQ.2 ) THEN
887 DO 150 j = 1, n
888 DO 140 i = 1, j
889 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
890 $ idist, iseed, d, igrade, dl, dr, ipvtng,
891 $ iwork, sparse )
892 a( isub, jsub ) = ctemp
893 a( jsub, isub ) = ctemp
894 140 CONTINUE
895 150 CONTINUE
896 END IF
897*
898 ELSE IF( ipack.EQ.1 ) THEN
899*
900 DO 170 j = 1, n
901 DO 160 i = 1, j
902 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
903 $ idist,
904 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
905 $ sparse )
906 mnsub = min( isub, jsub )
907 mxsub = max( isub, jsub )
908 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
909 a( mnsub, mxsub ) = dconjg( ctemp )
910 ELSE
911 a( mnsub, mxsub ) = ctemp
912 END IF
913 IF( mnsub.NE.mxsub )
914 $ a( mxsub, mnsub ) = czero
915 160 CONTINUE
916 170 CONTINUE
917*
918 ELSE IF( ipack.EQ.2 ) THEN
919*
920 DO 190 j = 1, n
921 DO 180 i = 1, j
922 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
923 $ idist,
924 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
925 $ sparse )
926 mnsub = min( isub, jsub )
927 mxsub = max( isub, jsub )
928 IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
929 a( mxsub, mnsub ) = dconjg( ctemp )
930 ELSE
931 a( mxsub, mnsub ) = ctemp
932 END IF
933 IF( mnsub.NE.mxsub )
934 $ a( mnsub, mxsub ) = czero
935 180 CONTINUE
936 190 CONTINUE
937*
938 ELSE IF( ipack.EQ.3 ) THEN
939*
940 DO 210 j = 1, n
941 DO 200 i = 1, j
942 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
943 $ idist,
944 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
945 $ sparse )
946*
947* Compute K = location of (ISUB,JSUB) entry in packed
948* array
949*
950 mnsub = min( isub, jsub )
951 mxsub = max( isub, jsub )
952 k = mxsub*( mxsub-1 ) / 2 + mnsub
953*
954* Convert K to (IISUB,JJSUB) location
955*
956 jjsub = ( k-1 ) / lda + 1
957 iisub = k - lda*( jjsub-1 )
958*
959 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
960 a( iisub, jjsub ) = dconjg( ctemp )
961 ELSE
962 a( iisub, jjsub ) = ctemp
963 END IF
964 200 CONTINUE
965 210 CONTINUE
966*
967 ELSE IF( ipack.EQ.4 ) THEN
968*
969 DO 230 j = 1, n
970 DO 220 i = 1, j
971 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
972 $ idist,
973 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
974 $ sparse )
975*
976* Compute K = location of (I,J) entry in packed array
977*
978 mnsub = min( isub, jsub )
979 mxsub = max( isub, jsub )
980 IF( mnsub.EQ.1 ) THEN
981 k = mxsub
982 ELSE
983 k = n*( n+1 ) / 2 - ( n-mnsub+1 )*( n-mnsub+2 ) /
984 $ 2 + mxsub - mnsub + 1
985 END IF
986*
987* Convert K to (IISUB,JJSUB) location
988*
989 jjsub = ( k-1 ) / lda + 1
990 iisub = k - lda*( jjsub-1 )
991*
992 IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
993 a( iisub, jjsub ) = dconjg( ctemp )
994 ELSE
995 a( iisub, jjsub ) = ctemp
996 END IF
997 220 CONTINUE
998 230 CONTINUE
999*
1000 ELSE IF( ipack.EQ.5 ) THEN
1001*
1002 DO 250 j = 1, n
1003 DO 240 i = j - kuu, j
1004 IF( i.LT.1 ) THEN
1005 a( j-i+1, i+n ) = czero
1006 ELSE
1007 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1008 $ idist, iseed, d, igrade, dl, dr, ipvtng,
1009 $ iwork, sparse )
1010 mnsub = min( isub, jsub )
1011 mxsub = max( isub, jsub )
1012 IF( mxsub.EQ.jsub .AND. isym.EQ.0 ) THEN
1013 a( mxsub-mnsub+1, mnsub ) = dconjg( ctemp )
1014 ELSE
1015 a( mxsub-mnsub+1, mnsub ) = ctemp
1016 END IF
1017 END IF
1018 240 CONTINUE
1019 250 CONTINUE
1020*
1021 ELSE IF( ipack.EQ.6 ) THEN
1022*
1023 DO 270 j = 1, n
1024 DO 260 i = j - kuu, j
1025 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1026 $ idist,
1027 $ iseed, d, igrade, dl, dr, ipvtng, iwork,
1028 $ sparse )
1029 mnsub = min( isub, jsub )
1030 mxsub = max( isub, jsub )
1031 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
1032 a( mnsub-mxsub+kuu+1, mxsub ) = dconjg( ctemp )
1033 ELSE
1034 a( mnsub-mxsub+kuu+1, mxsub ) = ctemp
1035 END IF
1036 260 CONTINUE
1037 270 CONTINUE
1038*
1039 ELSE IF( ipack.EQ.7 ) THEN
1040*
1041 IF( isym.NE.1 ) THEN
1042 DO 290 j = 1, n
1043 DO 280 i = j - kuu, j
1044 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1045 $ idist, iseed, d, igrade, dl, dr, ipvtng,
1046 $ iwork, sparse )
1047 mnsub = min( isub, jsub )
1048 mxsub = max( isub, jsub )
1049 IF( i.LT.1 )
1050 $ a( j-i+1+kuu, i+n ) = czero
1051 IF( mxsub.EQ.isub .AND. isym.EQ.0 ) THEN
1052 a( mnsub-mxsub+kuu+1, mxsub ) = dconjg( ctemp )
1053 ELSE
1054 a( mnsub-mxsub+kuu+1, mxsub ) = ctemp
1055 END IF
1056 IF( i.GE.1 .AND. mnsub.NE.mxsub ) THEN
1057 IF( mnsub.EQ.isub .AND. isym.EQ.0 ) THEN
1058 a( mxsub-mnsub+1+kuu,
1059 $ mnsub ) = dconjg( ctemp )
1060 ELSE
1061 a( mxsub-mnsub+1+kuu, mnsub ) = ctemp
1062 END IF
1063 END IF
1064 280 CONTINUE
1065 290 CONTINUE
1066 ELSE IF( isym.EQ.1 ) THEN
1067 DO 310 j = 1, n
1068 DO 300 i = j - kuu, j + kll
1069 ctemp = zlatm3( m, n, i, j, isub, jsub, kl, ku,
1070 $ idist, iseed, d, igrade, dl, dr, ipvtng,
1071 $ iwork, sparse )
1072 a( isub-jsub+kuu+1, jsub ) = ctemp
1073 300 CONTINUE
1074 310 CONTINUE
1075 END IF
1076*
1077 END IF
1078*
1079 ELSE
1080*
1081* Use ZLATM2
1082*
1083 IF( ipack.EQ.0 ) THEN
1084 IF( isym.EQ.0 ) THEN
1085 DO 330 j = 1, n
1086 DO 320 i = 1, j
1087 a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1088 $ iseed, d, igrade, dl, dr, ipvtng,
1089 $ iwork, sparse )
1090 a( j, i ) = dconjg( a( i, j ) )
1091 320 CONTINUE
1092 330 CONTINUE
1093 ELSE IF( isym.EQ.1 ) THEN
1094 DO 350 j = 1, n
1095 DO 340 i = 1, m
1096 a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1097 $ iseed, d, igrade, dl, dr, ipvtng,
1098 $ iwork, sparse )
1099 340 CONTINUE
1100 350 CONTINUE
1101 ELSE IF( isym.EQ.2 ) THEN
1102 DO 370 j = 1, n
1103 DO 360 i = 1, j
1104 a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1105 $ iseed, d, igrade, dl, dr, ipvtng,
1106 $ iwork, sparse )
1107 a( j, i ) = a( i, j )
1108 360 CONTINUE
1109 370 CONTINUE
1110 END IF
1111*
1112 ELSE IF( ipack.EQ.1 ) THEN
1113*
1114 DO 390 j = 1, n
1115 DO 380 i = 1, j
1116 a( i, j ) = zlatm2( m, n, i, j, kl, ku, idist,
1117 $ iseed,
1118 $ d, igrade, dl, dr, ipvtng, iwork, sparse )
1119 IF( i.NE.j )
1120 $ a( j, i ) = czero
1121 380 CONTINUE
1122 390 CONTINUE
1123*
1124 ELSE IF( ipack.EQ.2 ) THEN
1125*
1126 DO 410 j = 1, n
1127 DO 400 i = 1, j
1128 IF( isym.EQ.0 ) THEN
1129 a( j, i ) = dconjg( zlatm2( m, n, i, j, kl, ku,
1130 $ idist, iseed, d, igrade, dl, dr,
1131 $ ipvtng, iwork, sparse ) )
1132 ELSE
1133 a( j, i ) = zlatm2( m, n, i, j, kl, ku, idist,
1134 $ iseed, d, igrade, dl, dr, ipvtng,
1135 $ iwork, sparse )
1136 END IF
1137 IF( i.NE.j )
1138 $ a( i, j ) = czero
1139 400 CONTINUE
1140 410 CONTINUE
1141*
1142 ELSE IF( ipack.EQ.3 ) THEN
1143*
1144 isub = 0
1145 jsub = 1
1146 DO 430 j = 1, n
1147 DO 420 i = 1, j
1148 isub = isub + 1
1149 IF( isub.GT.lda ) THEN
1150 isub = 1
1151 jsub = jsub + 1
1152 END IF
1153 a( isub, jsub ) = zlatm2( m, n, i, j, kl, ku,
1154 $ idist,
1155 $ iseed, d, igrade, dl, dr, ipvtng,
1156 $ iwork, sparse )
1157 420 CONTINUE
1158 430 CONTINUE
1159*
1160 ELSE IF( ipack.EQ.4 ) THEN
1161*
1162 IF( isym.EQ.0 .OR. isym.EQ.2 ) THEN
1163 DO 450 j = 1, n
1164 DO 440 i = 1, j
1165*
1166* Compute K = location of (I,J) entry in packed array
1167*
1168 IF( i.EQ.1 ) THEN
1169 k = j
1170 ELSE
1171 k = n*( n+1 ) / 2 - ( n-i+1 )*( n-i+2 ) / 2 +
1172 $ j - i + 1
1173 END IF
1174*
1175* Convert K to (ISUB,JSUB) location
1176*
1177 jsub = ( k-1 ) / lda + 1
1178 isub = k - lda*( jsub-1 )
1179*
1180 a( isub, jsub ) = zlatm2( m, n, i, j, kl, ku,
1181 $ idist, iseed, d, igrade, dl, dr,
1182 $ ipvtng, iwork, sparse )
1183 IF( isym.EQ.0 )
1184 $ a( isub, jsub ) = dconjg( a( isub, jsub ) )
1185 440 CONTINUE
1186 450 CONTINUE
1187 ELSE
1188 isub = 0
1189 jsub = 1
1190 DO 470 j = 1, n
1191 DO 460 i = j, m
1192 isub = isub + 1
1193 IF( isub.GT.lda ) THEN
1194 isub = 1
1195 jsub = jsub + 1
1196 END IF
1197 a( isub, jsub ) = zlatm2( m, n, i, j, kl, ku,
1198 $ idist, iseed, d, igrade, dl, dr,
1199 $ ipvtng, iwork, sparse )
1200 460 CONTINUE
1201 470 CONTINUE
1202 END IF
1203*
1204 ELSE IF( ipack.EQ.5 ) THEN
1205*
1206 DO 490 j = 1, n
1207 DO 480 i = j - kuu, j
1208 IF( i.LT.1 ) THEN
1209 a( j-i+1, i+n ) = czero
1210 ELSE
1211 IF( isym.EQ.0 ) THEN
1212 a( j-i+1, i ) = dconjg( zlatm2( m, n, i, j,
1213 $ kl,
1214 $ ku, idist, iseed, d, igrade, dl,
1215 $ dr, ipvtng, iwork, sparse ) )
1216 ELSE
1217 a( j-i+1, i ) = zlatm2( m, n, i, j, kl, ku,
1218 $ idist, iseed, d, igrade, dl, dr,
1219 $ ipvtng, iwork, sparse )
1220 END IF
1221 END IF
1222 480 CONTINUE
1223 490 CONTINUE
1224*
1225 ELSE IF( ipack.EQ.6 ) THEN
1226*
1227 DO 510 j = 1, n
1228 DO 500 i = j - kuu, j
1229 a( i-j+kuu+1, j ) = zlatm2( m, n, i, j, kl, ku,
1230 $ idist,
1231 $ iseed, d, igrade, dl, dr, ipvtng,
1232 $ iwork, sparse )
1233 500 CONTINUE
1234 510 CONTINUE
1235*
1236 ELSE IF( ipack.EQ.7 ) THEN
1237*
1238 IF( isym.NE.1 ) THEN
1239 DO 530 j = 1, n
1240 DO 520 i = j - kuu, j
1241 a( i-j+kuu+1, j ) = zlatm2( m, n, i, j, kl, ku,
1242 $ idist, iseed, d, igrade, dl,
1243 $ dr, ipvtng, iwork, sparse )
1244 IF( i.LT.1 )
1245 $ a( j-i+1+kuu, i+n ) = czero
1246 IF( i.GE.1 .AND. i.NE.j ) THEN
1247 IF( isym.EQ.0 ) THEN
1248 a( j-i+1+kuu, i ) = dconjg( a( i-j+kuu+1,
1249 $ j ) )
1250 ELSE
1251 a( j-i+1+kuu, i ) = a( i-j+kuu+1, j )
1252 END IF
1253 END IF
1254 520 CONTINUE
1255 530 CONTINUE
1256 ELSE IF( isym.EQ.1 ) THEN
1257 DO 550 j = 1, n
1258 DO 540 i = j - kuu, j + kll
1259 a( i-j+kuu+1, j ) = zlatm2( m, n, i, j, kl, ku,
1260 $ idist, iseed, d, igrade, dl,
1261 $ dr, ipvtng, iwork, sparse )
1262 540 CONTINUE
1263 550 CONTINUE
1264 END IF
1265*
1266 END IF
1267*
1268 END IF
1269*
1270* 5) Scaling the norm
1271*
1272 IF( ipack.EQ.0 ) THEN
1273 onorm = zlange( 'M', m, n, a, lda, tempa )
1274 ELSE IF( ipack.EQ.1 ) THEN
1275 onorm = zlansy( 'M', 'U', n, a, lda, tempa )
1276 ELSE IF( ipack.EQ.2 ) THEN
1277 onorm = zlansy( 'M', 'L', n, a, lda, tempa )
1278 ELSE IF( ipack.EQ.3 ) THEN
1279 onorm = zlansp( 'M', 'U', n, a, tempa )
1280 ELSE IF( ipack.EQ.4 ) THEN
1281 onorm = zlansp( 'M', 'L', n, a, tempa )
1282 ELSE IF( ipack.EQ.5 ) THEN
1283 onorm = zlansb( 'M', 'L', n, kll, a, lda, tempa )
1284 ELSE IF( ipack.EQ.6 ) THEN
1285 onorm = zlansb( 'M', 'U', n, kuu, a, lda, tempa )
1286 ELSE IF( ipack.EQ.7 ) THEN
1287 onorm = zlangb( 'M', n, kll, kuu, a, lda, tempa )
1288 END IF
1289*
1290 IF( anorm.GE.zero ) THEN
1291*
1292 IF( anorm.GT.zero .AND. onorm.EQ.zero ) THEN
1293*
1294* Desired scaling impossible
1295*
1296 info = 5
1297 RETURN
1298*
1299 ELSE IF( ( anorm.GT.one .AND. onorm.LT.one ) .OR.
1300 $ ( anorm.LT.one .AND. onorm.GT.one ) ) THEN
1301*
1302* Scale carefully to avoid over / underflow
1303*
1304 IF( ipack.LE.2 ) THEN
1305 DO 560 j = 1, n
1306 CALL zdscal( m, one / onorm, a( 1, j ), 1 )
1307 CALL zdscal( m, anorm, a( 1, j ), 1 )
1308 560 CONTINUE
1309*
1310 ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1311*
1312 CALL zdscal( n*( n+1 ) / 2, one / onorm, a, 1 )
1313 CALL zdscal( n*( n+1 ) / 2, anorm, a, 1 )
1314*
1315 ELSE IF( ipack.GE.5 ) THEN
1316*
1317 DO 570 j = 1, n
1318 CALL zdscal( kll+kuu+1, one / onorm, a( 1, j ), 1 )
1319 CALL zdscal( kll+kuu+1, anorm, a( 1, j ), 1 )
1320 570 CONTINUE
1321*
1322 END IF
1323*
1324 ELSE
1325*
1326* Scale straightforwardly
1327*
1328 IF( ipack.LE.2 ) THEN
1329 DO 580 j = 1, n
1330 CALL zdscal( m, anorm / onorm, a( 1, j ), 1 )
1331 580 CONTINUE
1332*
1333 ELSE IF( ipack.EQ.3 .OR. ipack.EQ.4 ) THEN
1334*
1335 CALL zdscal( n*( n+1 ) / 2, anorm / onorm, a, 1 )
1336*
1337 ELSE IF( ipack.GE.5 ) THEN
1338*
1339 DO 590 j = 1, n
1340 CALL zdscal( kll+kuu+1, anorm / onorm, a( 1, j ), 1 )
1341 590 CONTINUE
1342 END IF
1343*
1344 END IF
1345*
1346 END IF
1347*
1348* End of ZLATMR
1349*
1350 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zlatm1(mode, cond, irsign, idist, iseed, d, n, info)
ZLATM1
Definition zlatm1.f:138
subroutine zlatmr(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)
ZLATMR
Definition zlatmr.f:490