LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
cblat2.f
Go to the documentation of this file.
1 *> \brief \b CBLAT2
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * PROGRAM CBLAT2
12 *
13 *
14 *> \par Purpose:
15 * =============
16 *>
17 *> \verbatim
18 *>
19 *> Test program for the COMPLEX Level 2 Blas.
20 *>
21 *> The program must be driven by a short data file. The first 18 records
22 *> of the file are read using list-directed input, the last 17 records
23 *> are read using the format ( A6, L2 ). An annotated example of a data
24 *> file can be obtained by deleting the first 3 characters from the
25 *> following 35 lines:
26 *> 'cblat2.out' NAME OF SUMMARY OUTPUT FILE
27 *> 6 UNIT NUMBER OF SUMMARY FILE
28 *> 'CBLA2T.SNAP' NAME OF SNAPSHOT OUTPUT FILE
29 *> -1 UNIT NUMBER OF SNAPSHOT FILE (NOT USED IF .LT. 0)
30 *> F LOGICAL FLAG, T TO REWIND SNAPSHOT FILE AFTER EACH RECORD.
31 *> F LOGICAL FLAG, T TO STOP ON FAILURES.
32 *> T LOGICAL FLAG, T TO TEST ERROR EXITS.
33 *> 16.0 THRESHOLD VALUE OF TEST RATIO
34 *> 6 NUMBER OF VALUES OF N
35 *> 0 1 2 3 5 9 VALUES OF N
36 *> 4 NUMBER OF VALUES OF K
37 *> 0 1 2 4 VALUES OF K
38 *> 4 NUMBER OF VALUES OF INCX AND INCY
39 *> 1 2 -1 -2 VALUES OF INCX AND INCY
40 *> 3 NUMBER OF VALUES OF ALPHA
41 *> (0.0,0.0) (1.0,0.0) (0.7,-0.9) VALUES OF ALPHA
42 *> 3 NUMBER OF VALUES OF BETA
43 *> (0.0,0.0) (1.0,0.0) (1.3,-1.1) VALUES OF BETA
44 *> CGEMV T PUT F FOR NO TEST. SAME COLUMNS.
45 *> CGBMV T PUT F FOR NO TEST. SAME COLUMNS.
46 *> CHEMV T PUT F FOR NO TEST. SAME COLUMNS.
47 *> CHBMV T PUT F FOR NO TEST. SAME COLUMNS.
48 *> CHPMV T PUT F FOR NO TEST. SAME COLUMNS.
49 *> CTRMV T PUT F FOR NO TEST. SAME COLUMNS.
50 *> CTBMV T PUT F FOR NO TEST. SAME COLUMNS.
51 *> CTPMV T PUT F FOR NO TEST. SAME COLUMNS.
52 *> CTRSV T PUT F FOR NO TEST. SAME COLUMNS.
53 *> CTBSV T PUT F FOR NO TEST. SAME COLUMNS.
54 *> CTPSV T PUT F FOR NO TEST. SAME COLUMNS.
55 *> CGERC T PUT F FOR NO TEST. SAME COLUMNS.
56 *> CGERU T PUT F FOR NO TEST. SAME COLUMNS.
57 *> CHER T PUT F FOR NO TEST. SAME COLUMNS.
58 *> CHPR T PUT F FOR NO TEST. SAME COLUMNS.
59 *> CHER2 T PUT F FOR NO TEST. SAME COLUMNS.
60 *> CHPR2 T PUT F FOR NO TEST. SAME COLUMNS.
61 *>
62 *> Further Details
63 *> ===============
64 *>
65 *> See:
66 *>
67 *> Dongarra J. J., Du Croz J. J., Hammarling S. and Hanson R. J..
68 *> An extended set of Fortran Basic Linear Algebra Subprograms.
69 *>
70 *> Technical Memoranda Nos. 41 (revision 3) and 81, Mathematics
71 *> and Computer Science Division, Argonne National Laboratory,
72 *> 9700 South Cass Avenue, Argonne, Illinois 60439, US.
73 *>
74 *> Or
75 *>
76 *> NAG Technical Reports TR3/87 and TR4/87, Numerical Algorithms
77 *> Group Ltd., NAG Central Office, 256 Banbury Road, Oxford
78 *> OX2 7DE, UK, and Numerical Algorithms Group Inc., 1101 31st
79 *> Street, Suite 100, Downers Grove, Illinois 60515-1263, USA.
80 *>
81 *>
82 *> -- Written on 10-August-1987.
83 *> Richard Hanson, Sandia National Labs.
84 *> Jeremy Du Croz, NAG Central Office.
85 *>
86 *> 10-9-00: Change STATUS='NEW' to 'UNKNOWN' so that the testers
87 *> can be run multiple times without deleting generated
88 *> output files (susan)
89 *> \endverbatim
90 *
91 * Authors:
92 * ========
93 *
94 *> \author Univ. of Tennessee
95 *> \author Univ. of California Berkeley
96 *> \author Univ. of Colorado Denver
97 *> \author NAG Ltd.
98 *
99 *> \date April 2012
100 *
101 *> \ingroup complex_blas_testing
102 *
103 * =====================================================================
104  PROGRAM cblat2
105 *
106 * -- Reference BLAS test routine (version 3.4.1) --
107 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
108 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
109 * April 2012
110 *
111 * =====================================================================
112 *
113 * .. Parameters ..
114  INTEGER nin
115  parameter( nin = 5 )
116  INTEGER nsubs
117  parameter( nsubs = 17 )
118  COMPLEX zero, one
119  parameter( zero = ( 0.0, 0.0 ), one = ( 1.0, 0.0 ) )
120  REAL rzero
121  parameter( rzero = 0.0 )
122  INTEGER nmax, incmax
123  parameter( nmax = 65, incmax = 2 )
124  INTEGER ninmax, nidmax, nkbmax, nalmax, nbemax
125  parameter( ninmax = 7, nidmax = 9, nkbmax = 7,
126  $ nalmax = 7, nbemax = 7 )
127 * .. Local Scalars ..
128  REAL eps, err, thresh
129  INTEGER i, isnum, j, n, nalf, nbet, nidim, ninc, nkb,
130  $ nout, ntra
131  LOGICAL fatal, ltestt, rewi, same, sfatal, trace,
132  $ tsterr
133  CHARACTER*1 trans
134  CHARACTER*6 snamet
135  CHARACTER*32 snaps, summry
136 * .. Local Arrays ..
137  COMPLEX a( nmax, nmax ), aa( nmax*nmax ),
138  $ alf( nalmax ), as( nmax*nmax ), bet( nbemax ),
139  $ x( nmax ), xs( nmax*incmax ),
140  $ xx( nmax*incmax ), y( nmax ),
141  $ ys( nmax*incmax ), yt( nmax ),
142  $ yy( nmax*incmax ), z( 2*nmax )
143  REAL g( nmax )
144  INTEGER idim( nidmax ), inc( ninmax ), kb( nkbmax )
145  LOGICAL ltest( nsubs )
146  CHARACTER*6 snames( nsubs )
147 * .. External Functions ..
148  REAL sdiff
149  LOGICAL lce
150  EXTERNAL sdiff, lce
151 * .. External Subroutines ..
152  EXTERNAL cchk1, cchk2, cchk3, cchk4, cchk5, cchk6,
153  $ cchke, cmvch
154 * .. Intrinsic Functions ..
155  INTRINSIC abs, max, min
156 * .. Scalars in Common ..
157  INTEGER infot, noutc
158  LOGICAL lerr, ok
159  CHARACTER*6 srnamt
160 * .. Common blocks ..
161  common /infoc/infot, noutc, ok, lerr
162  common /srnamc/srnamt
163 * .. Data statements ..
164  DATA snames/'CGEMV ', 'CGBMV ', 'CHEMV ', 'CHBMV ',
165  $ 'CHPMV ', 'CTRMV ', 'CTBMV ', 'CTPMV ',
166  $ 'CTRSV ', 'CTBSV ', 'CTPSV ', 'CGERC ',
167  $ 'CGERU ', 'CHER ', 'CHPR ', 'CHER2 ',
168  $ 'CHPR2 '/
169 * .. Executable Statements ..
170 *
171 * Read name and unit number for summary output file and open file.
172 *
173  READ( nin, fmt = * )summry
174  READ( nin, fmt = * )nout
175  OPEN( nout, file = summry, status = 'UNKNOWN' )
176  noutc = nout
177 *
178 * Read name and unit number for snapshot output file and open file.
179 *
180  READ( nin, fmt = * )snaps
181  READ( nin, fmt = * )ntra
182  trace = ntra.GE.0
183  IF( trace )THEN
184  OPEN( ntra, file = snaps, status = 'UNKNOWN' )
185  END IF
186 * Read the flag that directs rewinding of the snapshot file.
187  READ( nin, fmt = * )rewi
188  rewi = rewi.AND.trace
189 * Read the flag that directs stopping on any failure.
190  READ( nin, fmt = * )sfatal
191 * Read the flag that indicates whether error exits are to be tested.
192  READ( nin, fmt = * )tsterr
193 * Read the threshold value of the test ratio
194  READ( nin, fmt = * )thresh
195 *
196 * Read and check the parameter values for the tests.
197 *
198 * Values of N
199  READ( nin, fmt = * )nidim
200  IF( nidim.LT.1.OR.nidim.GT.nidmax )THEN
201  WRITE( nout, fmt = 9997 )'N', nidmax
202  go to 230
203  END IF
204  READ( nin, fmt = * )( idim( i ), i = 1, nidim )
205  DO 10 i = 1, nidim
206  IF( idim( i ).LT.0.OR.idim( i ).GT.nmax )THEN
207  WRITE( nout, fmt = 9996 )nmax
208  go to 230
209  END IF
210  10 continue
211 * Values of K
212  READ( nin, fmt = * )nkb
213  IF( nkb.LT.1.OR.nkb.GT.nkbmax )THEN
214  WRITE( nout, fmt = 9997 )'K', nkbmax
215  go to 230
216  END IF
217  READ( nin, fmt = * )( kb( i ), i = 1, nkb )
218  DO 20 i = 1, nkb
219  IF( kb( i ).LT.0 )THEN
220  WRITE( nout, fmt = 9995 )
221  go to 230
222  END IF
223  20 continue
224 * Values of INCX and INCY
225  READ( nin, fmt = * )ninc
226  IF( ninc.LT.1.OR.ninc.GT.ninmax )THEN
227  WRITE( nout, fmt = 9997 )'INCX AND INCY', ninmax
228  go to 230
229  END IF
230  READ( nin, fmt = * )( inc( i ), i = 1, ninc )
231  DO 30 i = 1, ninc
232  IF( inc( i ).EQ.0.OR.abs( inc( i ) ).GT.incmax )THEN
233  WRITE( nout, fmt = 9994 )incmax
234  go to 230
235  END IF
236  30 continue
237 * Values of ALPHA
238  READ( nin, fmt = * )nalf
239  IF( nalf.LT.1.OR.nalf.GT.nalmax )THEN
240  WRITE( nout, fmt = 9997 )'ALPHA', nalmax
241  go to 230
242  END IF
243  READ( nin, fmt = * )( alf( i ), i = 1, nalf )
244 * Values of BETA
245  READ( nin, fmt = * )nbet
246  IF( nbet.LT.1.OR.nbet.GT.nbemax )THEN
247  WRITE( nout, fmt = 9997 )'BETA', nbemax
248  go to 230
249  END IF
250  READ( nin, fmt = * )( bet( i ), i = 1, nbet )
251 *
252 * Report values of parameters.
253 *
254  WRITE( nout, fmt = 9993 )
255  WRITE( nout, fmt = 9992 )( idim( i ), i = 1, nidim )
256  WRITE( nout, fmt = 9991 )( kb( i ), i = 1, nkb )
257  WRITE( nout, fmt = 9990 )( inc( i ), i = 1, ninc )
258  WRITE( nout, fmt = 9989 )( alf( i ), i = 1, nalf )
259  WRITE( nout, fmt = 9988 )( bet( i ), i = 1, nbet )
260  IF( .NOT.tsterr )THEN
261  WRITE( nout, fmt = * )
262  WRITE( nout, fmt = 9980 )
263  END IF
264  WRITE( nout, fmt = * )
265  WRITE( nout, fmt = 9999 )thresh
266  WRITE( nout, fmt = * )
267 *
268 * Read names of subroutines and flags which indicate
269 * whether they are to be tested.
270 *
271  DO 40 i = 1, nsubs
272  ltest( i ) = .false.
273  40 continue
274  50 READ( nin, fmt = 9984, END = 80 )snamet, ltestt
275  DO 60 i = 1, nsubs
276  IF( snamet.EQ.snames( i ) )
277  $ go to 70
278  60 continue
279  WRITE( nout, fmt = 9986 )snamet
280  stop
281  70 ltest( i ) = ltestt
282  go to 50
283 *
284  80 continue
285  CLOSE ( nin )
286 *
287 * Compute EPS (the machine precision).
288 *
289  eps = epsilon(rzero)
290  WRITE( nout, fmt = 9998 )eps
291 *
292 * Check the reliability of CMVCH using exact data.
293 *
294  n = min( 32, nmax )
295  DO 120 j = 1, n
296  DO 110 i = 1, n
297  a( i, j ) = max( i - j + 1, 0 )
298  110 continue
299  x( j ) = j
300  y( j ) = zero
301  120 continue
302  DO 130 j = 1, n
303  yy( j ) = j*( ( j + 1 )*j )/2 - ( ( j + 1 )*j*( j - 1 ) )/3
304  130 continue
305 * YY holds the exact result. On exit from CMVCH YT holds
306 * the result computed by CMVCH.
307  trans = 'N'
308  CALL cmvch( trans, n, n, one, a, nmax, x, 1, zero, y, 1, yt, g,
309  $ yy, eps, err, fatal, nout, .true. )
310  same = lce( yy, yt, n )
311  IF( .NOT.same.OR.err.NE.rzero )THEN
312  WRITE( nout, fmt = 9985 )trans, same, err
313  stop
314  END IF
315  trans = 'T'
316  CALL cmvch( trans, n, n, one, a, nmax, x, -1, zero, y, -1, yt, g,
317  $ yy, eps, err, fatal, nout, .true. )
318  same = lce( yy, yt, n )
319  IF( .NOT.same.OR.err.NE.rzero )THEN
320  WRITE( nout, fmt = 9985 )trans, same, err
321  stop
322  END IF
323 *
324 * Test each subroutine in turn.
325 *
326  DO 210 isnum = 1, nsubs
327  WRITE( nout, fmt = * )
328  IF( .NOT.ltest( isnum ) )THEN
329 * Subprogram is not to be tested.
330  WRITE( nout, fmt = 9983 )snames( isnum )
331  ELSE
332  srnamt = snames( isnum )
333 * Test error exits.
334  IF( tsterr )THEN
335  CALL cchke( isnum, snames( isnum ), nout )
336  WRITE( nout, fmt = * )
337  END IF
338 * Test computations.
339  infot = 0
340  ok = .true.
341  fatal = .false.
342  go to( 140, 140, 150, 150, 150, 160, 160,
343  $ 160, 160, 160, 160, 170, 170, 180,
344  $ 180, 190, 190 )isnum
345 * Test CGEMV, 01, and CGBMV, 02.
346  140 CALL cchk1( snames( isnum ), eps, thresh, nout, ntra, trace,
347  $ rewi, fatal, nidim, idim, nkb, kb, nalf, alf,
348  $ nbet, bet, ninc, inc, nmax, incmax, a, aa, as,
349  $ x, xx, xs, y, yy, ys, yt, g )
350  go to 200
351 * Test CHEMV, 03, CHBMV, 04, and CHPMV, 05.
352  150 CALL cchk2( snames( isnum ), eps, thresh, nout, ntra, trace,
353  $ rewi, fatal, nidim, idim, nkb, kb, nalf, alf,
354  $ nbet, bet, ninc, inc, nmax, incmax, a, aa, as,
355  $ x, xx, xs, y, yy, ys, yt, g )
356  go to 200
357 * Test CTRMV, 06, CTBMV, 07, CTPMV, 08,
358 * CTRSV, 09, CTBSV, 10, and CTPSV, 11.
359  160 CALL cchk3( snames( isnum ), eps, thresh, nout, ntra, trace,
360  $ rewi, fatal, nidim, idim, nkb, kb, ninc, inc,
361  $ nmax, incmax, a, aa, as, y, yy, ys, yt, g, z )
362  go to 200
363 * Test CGERC, 12, CGERU, 13.
364  170 CALL cchk4( snames( isnum ), eps, thresh, nout, ntra, trace,
365  $ rewi, fatal, nidim, idim, nalf, alf, ninc, inc,
366  $ nmax, incmax, a, aa, as, x, xx, xs, y, yy, ys,
367  $ yt, g, z )
368  go to 200
369 * Test CHER, 14, and CHPR, 15.
370  180 CALL cchk5( snames( isnum ), eps, thresh, nout, ntra, trace,
371  $ rewi, fatal, nidim, idim, nalf, alf, ninc, inc,
372  $ nmax, incmax, a, aa, as, x, xx, xs, y, yy, ys,
373  $ yt, g, z )
374  go to 200
375 * Test CHER2, 16, and CHPR2, 17.
376  190 CALL cchk6( snames( isnum ), eps, thresh, nout, ntra, trace,
377  $ rewi, fatal, nidim, idim, nalf, alf, ninc, inc,
378  $ nmax, incmax, a, aa, as, x, xx, xs, y, yy, ys,
379  $ yt, g, z )
380 *
381  200 IF( fatal.AND.sfatal )
382  $ go to 220
383  END IF
384  210 continue
385  WRITE( nout, fmt = 9982 )
386  go to 240
387 *
388  220 continue
389  WRITE( nout, fmt = 9981 )
390  go to 240
391 *
392  230 continue
393  WRITE( nout, fmt = 9987 )
394 *
395  240 continue
396  IF( trace )
397  $ CLOSE ( ntra )
398  CLOSE ( nout )
399  stop
400 *
401  9999 format( ' ROUTINES PASS COMPUTATIONAL TESTS IF TEST RATIO IS LES',
402  $ 'S THAN', f8.2 )
403  9998 format( ' RELATIVE MACHINE PRECISION IS TAKEN TO BE', 1p, e9.1 )
404  9997 format( ' NUMBER OF VALUES OF ', a, ' IS LESS THAN 1 OR GREATER ',
405  $ 'THAN ', i2 )
406  9996 format( ' VALUE OF N IS LESS THAN 0 OR GREATER THAN ', i2 )
407  9995 format( ' VALUE OF K IS LESS THAN 0' )
408  9994 format( ' ABSOLUTE VALUE OF INCX OR INCY IS 0 OR GREATER THAN ',
409  $ i2 )
410  9993 format( ' TESTS OF THE COMPLEX LEVEL 2 BLAS', //' THE F',
411  $ 'OLLOWING PARAMETER VALUES WILL BE USED:' )
412  9992 format( ' FOR N ', 9i6 )
413  9991 format( ' FOR K ', 7i6 )
414  9990 format( ' FOR INCX AND INCY ', 7i6 )
415  9989 format( ' FOR ALPHA ',
416  $ 7( '(', f4.1, ',', f4.1, ') ', : ) )
417  9988 format( ' FOR BETA ',
418  $ 7( '(', f4.1, ',', f4.1, ') ', : ) )
419  9987 format( ' AMEND DATA FILE OR INCREASE ARRAY SIZES IN PROGRAM',
420  $ /' ******* TESTS ABANDONED *******' )
421  9986 format( ' SUBPROGRAM NAME ', a6, ' NOT RECOGNIZED', /' ******* T',
422  $ 'ESTS ABANDONED *******' )
423  9985 format( ' ERROR IN CMVCH - IN-LINE DOT PRODUCTS ARE BEING EVALU',
424  $ 'ATED WRONGLY.', /' CMVCH WAS CALLED WITH TRANS = ', a1,
425  $ ' AND RETURNED SAME = ', l1, ' AND ERR = ', f12.3, '.', /
426  $ ' THIS MAY BE DUE TO FAULTS IN THE ARITHMETIC OR THE COMPILER.'
427  $ , /' ******* TESTS ABANDONED *******' )
428  9984 format( a6, l2 )
429  9983 format( 1x, a6, ' WAS NOT TESTED' )
430  9982 format( /' END OF TESTS' )
431  9981 format( /' ******* FATAL ERROR - TESTS ABANDONED *******' )
432  9980 format( ' ERROR-EXITS WILL NOT BE TESTED' )
433 *
434 * End of CBLAT2.
435 *
436  END
437  SUBROUTINE cchk1( SNAME, EPS, THRESH, NOUT, NTRA, TRACE, REWI,
438  $ fatal, nidim, idim, nkb, kb, nalf, alf, nbet,
439  $ bet, ninc, inc, nmax, incmax, a, aa, as, x, xx,
440  $ xs, y, yy, ys, yt, g )
441 *
442 * Tests CGEMV and CGBMV.
443 *
444 * Auxiliary routine for test program for Level 2 Blas.
445 *
446 * -- Written on 10-August-1987.
447 * Richard Hanson, Sandia National Labs.
448 * Jeremy Du Croz, NAG Central Office.
449 *
450 * .. Parameters ..
451  COMPLEX zero, half
452  parameter( zero = ( 0.0, 0.0 ), half = ( 0.5, 0.0 ) )
453  REAL rzero
454  parameter( rzero = 0.0 )
455 * .. Scalar Arguments ..
456  REAL eps, thresh
457  INTEGER incmax, nalf, nbet, nidim, ninc, nkb, nmax,
458  $ nout, ntra
459  LOGICAL fatal, rewi, trace
460  CHARACTER*6 sname
461 * .. Array Arguments ..
462  COMPLEX a( nmax, nmax ), aa( nmax*nmax ), alf( nalf ),
463  $ as( nmax*nmax ), bet( nbet ), x( nmax ),
464  $ xs( nmax*incmax ), xx( nmax*incmax ),
465  $ y( nmax ), ys( nmax*incmax ), yt( nmax ),
466  $ yy( nmax*incmax )
467  REAL g( nmax )
468  INTEGER idim( nidim ), inc( ninc ), kb( nkb )
469 * .. Local Scalars ..
470  COMPLEX alpha, als, beta, bls, transl
471  REAL err, errmax
472  INTEGER i, ia, ib, ic, iku, im, in, incx, incxs, incy,
473  $ incys, ix, iy, kl, kls, ku, kus, laa, lda,
474  $ ldas, lx, ly, m, ml, ms, n, nargs, nc, nd, nk,
475  $ nl, ns
476  LOGICAL banded, full, null, reset, same, tran
477  CHARACTER*1 trans, transs
478  CHARACTER*3 ich
479 * .. Local Arrays ..
480  LOGICAL isame( 13 )
481 * .. External Functions ..
482  LOGICAL lce, lceres
483  EXTERNAL lce, lceres
484 * .. External Subroutines ..
485  EXTERNAL cgbmv, cgemv, cmake, cmvch
486 * .. Intrinsic Functions ..
487  INTRINSIC abs, max, min
488 * .. Scalars in Common ..
489  INTEGER infot, noutc
490  LOGICAL lerr, ok
491 * .. Common blocks ..
492  common /infoc/infot, noutc, ok, lerr
493 * .. Data statements ..
494  DATA ich/'NTC'/
495 * .. Executable Statements ..
496  full = sname( 3: 3 ).EQ.'E'
497  banded = sname( 3: 3 ).EQ.'B'
498 * Define the number of arguments.
499  IF( full )THEN
500  nargs = 11
501  ELSE IF( banded )THEN
502  nargs = 13
503  END IF
504 *
505  nc = 0
506  reset = .true.
507  errmax = rzero
508 *
509  DO 120 in = 1, nidim
510  n = idim( in )
511  nd = n/2 + 1
512 *
513  DO 110 im = 1, 2
514  IF( im.EQ.1 )
515  $ m = max( n - nd, 0 )
516  IF( im.EQ.2 )
517  $ m = min( n + nd, nmax )
518 *
519  IF( banded )THEN
520  nk = nkb
521  ELSE
522  nk = 1
523  END IF
524  DO 100 iku = 1, nk
525  IF( banded )THEN
526  ku = kb( iku )
527  kl = max( ku - 1, 0 )
528  ELSE
529  ku = n - 1
530  kl = m - 1
531  END IF
532 * Set LDA to 1 more than minimum value if room.
533  IF( banded )THEN
534  lda = kl + ku + 1
535  ELSE
536  lda = m
537  END IF
538  IF( lda.LT.nmax )
539  $ lda = lda + 1
540 * Skip tests if not enough room.
541  IF( lda.GT.nmax )
542  $ go to 100
543  laa = lda*n
544  null = n.LE.0.OR.m.LE.0
545 *
546 * Generate the matrix A.
547 *
548  transl = zero
549  CALL cmake( sname( 2: 3 ), ' ', ' ', m, n, a, nmax, aa,
550  $ lda, kl, ku, reset, transl )
551 *
552  DO 90 ic = 1, 3
553  trans = ich( ic: ic )
554  tran = trans.EQ.'T'.OR.trans.EQ.'C'
555 *
556  IF( tran )THEN
557  ml = n
558  nl = m
559  ELSE
560  ml = m
561  nl = n
562  END IF
563 *
564  DO 80 ix = 1, ninc
565  incx = inc( ix )
566  lx = abs( incx )*nl
567 *
568 * Generate the vector X.
569 *
570  transl = half
571  CALL cmake( 'GE', ' ', ' ', 1, nl, x, 1, xx,
572  $ abs( incx ), 0, nl - 1, reset, transl )
573  IF( nl.GT.1 )THEN
574  x( nl/2 ) = zero
575  xx( 1 + abs( incx )*( nl/2 - 1 ) ) = zero
576  END IF
577 *
578  DO 70 iy = 1, ninc
579  incy = inc( iy )
580  ly = abs( incy )*ml
581 *
582  DO 60 ia = 1, nalf
583  alpha = alf( ia )
584 *
585  DO 50 ib = 1, nbet
586  beta = bet( ib )
587 *
588 * Generate the vector Y.
589 *
590  transl = zero
591  CALL cmake( 'GE', ' ', ' ', 1, ml, y, 1,
592  $ yy, abs( incy ), 0, ml - 1,
593  $ reset, transl )
594 *
595  nc = nc + 1
596 *
597 * Save every datum before calling the
598 * subroutine.
599 *
600  transs = trans
601  ms = m
602  ns = n
603  kls = kl
604  kus = ku
605  als = alpha
606  DO 10 i = 1, laa
607  as( i ) = aa( i )
608  10 continue
609  ldas = lda
610  DO 20 i = 1, lx
611  xs( i ) = xx( i )
612  20 continue
613  incxs = incx
614  bls = beta
615  DO 30 i = 1, ly
616  ys( i ) = yy( i )
617  30 continue
618  incys = incy
619 *
620 * Call the subroutine.
621 *
622  IF( full )THEN
623  IF( trace )
624  $ WRITE( ntra, fmt = 9994 )nc, sname,
625  $ trans, m, n, alpha, lda, incx, beta,
626  $ incy
627  IF( rewi )
628  $ rewind ntra
629  CALL cgemv( trans, m, n, alpha, aa,
630  $ lda, xx, incx, beta, yy,
631  $ incy )
632  ELSE IF( banded )THEN
633  IF( trace )
634  $ WRITE( ntra, fmt = 9995 )nc, sname,
635  $ trans, m, n, kl, ku, alpha, lda,
636  $ incx, beta, incy
637  IF( rewi )
638  $ rewind ntra
639  CALL cgbmv( trans, m, n, kl, ku, alpha,
640  $ aa, lda, xx, incx, beta,
641  $ yy, incy )
642  END IF
643 *
644 * Check if error-exit was taken incorrectly.
645 *
646  IF( .NOT.ok )THEN
647  WRITE( nout, fmt = 9993 )
648  fatal = .true.
649  go to 130
650  END IF
651 *
652 * See what data changed inside subroutines.
653 *
654  isame( 1 ) = trans.EQ.transs
655  isame( 2 ) = ms.EQ.m
656  isame( 3 ) = ns.EQ.n
657  IF( full )THEN
658  isame( 4 ) = als.EQ.alpha
659  isame( 5 ) = lce( as, aa, laa )
660  isame( 6 ) = ldas.EQ.lda
661  isame( 7 ) = lce( xs, xx, lx )
662  isame( 8 ) = incxs.EQ.incx
663  isame( 9 ) = bls.EQ.beta
664  IF( null )THEN
665  isame( 10 ) = lce( ys, yy, ly )
666  ELSE
667  isame( 10 ) = lceres( 'GE', ' ', 1,
668  $ ml, ys, yy,
669  $ abs( incy ) )
670  END IF
671  isame( 11 ) = incys.EQ.incy
672  ELSE IF( banded )THEN
673  isame( 4 ) = kls.EQ.kl
674  isame( 5 ) = kus.EQ.ku
675  isame( 6 ) = als.EQ.alpha
676  isame( 7 ) = lce( as, aa, laa )
677  isame( 8 ) = ldas.EQ.lda
678  isame( 9 ) = lce( xs, xx, lx )
679  isame( 10 ) = incxs.EQ.incx
680  isame( 11 ) = bls.EQ.beta
681  IF( null )THEN
682  isame( 12 ) = lce( ys, yy, ly )
683  ELSE
684  isame( 12 ) = lceres( 'GE', ' ', 1,
685  $ ml, ys, yy,
686  $ abs( incy ) )
687  END IF
688  isame( 13 ) = incys.EQ.incy
689  END IF
690 *
691 * If data was incorrectly changed, report
692 * and return.
693 *
694  same = .true.
695  DO 40 i = 1, nargs
696  same = same.AND.isame( i )
697  IF( .NOT.isame( i ) )
698  $ WRITE( nout, fmt = 9998 )i
699  40 continue
700  IF( .NOT.same )THEN
701  fatal = .true.
702  go to 130
703  END IF
704 *
705  IF( .NOT.null )THEN
706 *
707 * Check the result.
708 *
709  CALL cmvch( trans, m, n, alpha, a,
710  $ nmax, x, incx, beta, y,
711  $ incy, yt, g, yy, eps, err,
712  $ fatal, nout, .true. )
713  errmax = max( errmax, err )
714 * If got really bad answer, report and
715 * return.
716  IF( fatal )
717  $ go to 130
718  ELSE
719 * Avoid repeating tests with M.le.0 or
720 * N.le.0.
721  go to 110
722  END IF
723 *
724  50 continue
725 *
726  60 continue
727 *
728  70 continue
729 *
730  80 continue
731 *
732  90 continue
733 *
734  100 continue
735 *
736  110 continue
737 *
738  120 continue
739 *
740 * Report result.
741 *
742  IF( errmax.LT.thresh )THEN
743  WRITE( nout, fmt = 9999 )sname, nc
744  ELSE
745  WRITE( nout, fmt = 9997 )sname, nc, errmax
746  END IF
747  go to 140
748 *
749  130 continue
750  WRITE( nout, fmt = 9996 )sname
751  IF( full )THEN
752  WRITE( nout, fmt = 9994 )nc, sname, trans, m, n, alpha, lda,
753  $ incx, beta, incy
754  ELSE IF( banded )THEN
755  WRITE( nout, fmt = 9995 )nc, sname, trans, m, n, kl, ku,
756  $ alpha, lda, incx, beta, incy
757  END IF
758 *
759  140 continue
760  return
761 *
762  9999 format( ' ', a6, ' PASSED THE COMPUTATIONAL TESTS (', i6, ' CALL',
763  $ 'S)' )
764  9998 format( ' ******* FATAL ERROR - PARAMETER NUMBER ', i2, ' WAS CH',
765  $ 'ANGED INCORRECTLY *******' )
766  9997 format( ' ', a6, ' COMPLETED THE COMPUTATIONAL TESTS (', i6, ' C',
767  $ 'ALLS)', /' ******* BUT WITH MAXIMUM TEST RATIO', f8.2,
768  $ ' - SUSPECT *******' )
769  9996 format( ' ******* ', a6, ' FAILED ON CALL NUMBER:' )
770  9995 format( 1x, i6, ': ', a6, '(''', a1, ''',', 4( i3, ',' ), '(',
771  $ f4.1, ',', f4.1, '), A,', i3, ', X,', i2, ',(', f4.1, ',',
772  $ f4.1, '), Y,', i2, ') .' )
773  9994 format( 1x, i6, ': ', a6, '(''', a1, ''',', 2( i3, ',' ), '(',
774  $ f4.1, ',', f4.1, '), A,', i3, ', X,', i2, ',(', f4.1, ',',
775  $ f4.1, '), Y,', i2, ') .' )
776  9993 format( ' ******* FATAL ERROR - ERROR-EXIT TAKEN ON VALID CALL *',
777  $ '******' )
778 *
779 * End of CCHK1.
780 *
781  END
782  SUBROUTINE cchk2( SNAME, EPS, THRESH, NOUT, NTRA, TRACE, REWI,
783  $ fatal, nidim, idim, nkb, kb, nalf, alf, nbet,
784  $ bet, ninc, inc, nmax, incmax, a, aa, as, x, xx,
785  $ xs, y, yy, ys, yt, g )
786 *
787 * Tests CHEMV, CHBMV and CHPMV.
788 *
789 * Auxiliary routine for test program for Level 2 Blas.
790 *
791 * -- Written on 10-August-1987.
792 * Richard Hanson, Sandia National Labs.
793 * Jeremy Du Croz, NAG Central Office.
794 *
795 * .. Parameters ..
796  COMPLEX zero, half
797  parameter( zero = ( 0.0, 0.0 ), half = ( 0.5, 0.0 ) )
798  REAL rzero
799  parameter( rzero = 0.0 )
800 * .. Scalar Arguments ..
801  REAL eps, thresh
802  INTEGER incmax, nalf, nbet, nidim, ninc, nkb, nmax,
803  $ nout, ntra
804  LOGICAL fatal, rewi, trace
805  CHARACTER*6 sname
806 * .. Array Arguments ..
807  COMPLEX a( nmax, nmax ), aa( nmax*nmax ), alf( nalf ),
808  $ as( nmax*nmax ), bet( nbet ), x( nmax ),
809  $ xs( nmax*incmax ), xx( nmax*incmax ),
810  $ y( nmax ), ys( nmax*incmax ), yt( nmax ),
811  $ yy( nmax*incmax )
812  REAL g( nmax )
813  INTEGER idim( nidim ), inc( ninc ), kb( nkb )
814 * .. Local Scalars ..
815  COMPLEX alpha, als, beta, bls, transl
816  REAL err, errmax
817  INTEGER i, ia, ib, ic, ik, in, incx, incxs, incy,
818  $ incys, ix, iy, k, ks, laa, lda, ldas, lx, ly,
819  $ n, nargs, nc, nk, ns
820  LOGICAL banded, full, null, packed, reset, same
821  CHARACTER*1 uplo, uplos
822  CHARACTER*2 ich
823 * .. Local Arrays ..
824  LOGICAL isame( 13 )
825 * .. External Functions ..
826  LOGICAL lce, lceres
827  EXTERNAL lce, lceres
828 * .. External Subroutines ..
829  EXTERNAL chbmv, chemv, chpmv, cmake, cmvch
830 * .. Intrinsic Functions ..
831  INTRINSIC abs, max
832 * .. Scalars in Common ..
833  INTEGER infot, noutc
834  LOGICAL lerr, ok
835 * .. Common blocks ..
836  common /infoc/infot, noutc, ok, lerr
837 * .. Data statements ..
838  DATA ich/'UL'/
839 * .. Executable Statements ..
840  full = sname( 3: 3 ).EQ.'E'
841  banded = sname( 3: 3 ).EQ.'B'
842  packed = sname( 3: 3 ).EQ.'P'
843 * Define the number of arguments.
844  IF( full )THEN
845  nargs = 10
846  ELSE IF( banded )THEN
847  nargs = 11
848  ELSE IF( packed )THEN
849  nargs = 9
850  END IF
851 *
852  nc = 0
853  reset = .true.
854  errmax = rzero
855 *
856  DO 110 in = 1, nidim
857  n = idim( in )
858 *
859  IF( banded )THEN
860  nk = nkb
861  ELSE
862  nk = 1
863  END IF
864  DO 100 ik = 1, nk
865  IF( banded )THEN
866  k = kb( ik )
867  ELSE
868  k = n - 1
869  END IF
870 * Set LDA to 1 more than minimum value if room.
871  IF( banded )THEN
872  lda = k + 1
873  ELSE
874  lda = n
875  END IF
876  IF( lda.LT.nmax )
877  $ lda = lda + 1
878 * Skip tests if not enough room.
879  IF( lda.GT.nmax )
880  $ go to 100
881  IF( packed )THEN
882  laa = ( n*( n + 1 ) )/2
883  ELSE
884  laa = lda*n
885  END IF
886  null = n.LE.0
887 *
888  DO 90 ic = 1, 2
889  uplo = ich( ic: ic )
890 *
891 * Generate the matrix A.
892 *
893  transl = zero
894  CALL cmake( sname( 2: 3 ), uplo, ' ', n, n, a, nmax, aa,
895  $ lda, k, k, reset, transl )
896 *
897  DO 80 ix = 1, ninc
898  incx = inc( ix )
899  lx = abs( incx )*n
900 *
901 * Generate the vector X.
902 *
903  transl = half
904  CALL cmake( 'GE', ' ', ' ', 1, n, x, 1, xx,
905  $ abs( incx ), 0, n - 1, reset, transl )
906  IF( n.GT.1 )THEN
907  x( n/2 ) = zero
908  xx( 1 + abs( incx )*( n/2 - 1 ) ) = zero
909  END IF
910 *
911  DO 70 iy = 1, ninc
912  incy = inc( iy )
913  ly = abs( incy )*n
914 *
915  DO 60 ia = 1, nalf
916  alpha = alf( ia )
917 *
918  DO 50 ib = 1, nbet
919  beta = bet( ib )
920 *
921 * Generate the vector Y.
922 *
923  transl = zero
924  CALL cmake( 'GE', ' ', ' ', 1, n, y, 1, yy,
925  $ abs( incy ), 0, n - 1, reset,
926  $ transl )
927 *
928  nc = nc + 1
929 *
930 * Save every datum before calling the
931 * subroutine.
932 *
933  uplos = uplo
934  ns = n
935  ks = k
936  als = alpha
937  DO 10 i = 1, laa
938  as( i ) = aa( i )
939  10 continue
940  ldas = lda
941  DO 20 i = 1, lx
942  xs( i ) = xx( i )
943  20 continue
944  incxs = incx
945  bls = beta
946  DO 30 i = 1, ly
947  ys( i ) = yy( i )
948  30 continue
949  incys = incy
950 *
951 * Call the subroutine.
952 *
953  IF( full )THEN
954  IF( trace )
955  $ WRITE( ntra, fmt = 9993 )nc, sname,
956  $ uplo, n, alpha, lda, incx, beta, incy
957  IF( rewi )
958  $ rewind ntra
959  CALL chemv( uplo, n, alpha, aa, lda, xx,
960  $ incx, beta, yy, incy )
961  ELSE IF( banded )THEN
962  IF( trace )
963  $ WRITE( ntra, fmt = 9994 )nc, sname,
964  $ uplo, n, k, alpha, lda, incx, beta,
965  $ incy
966  IF( rewi )
967  $ rewind ntra
968  CALL chbmv( uplo, n, k, alpha, aa, lda,
969  $ xx, incx, beta, yy, incy )
970  ELSE IF( packed )THEN
971  IF( trace )
972  $ WRITE( ntra, fmt = 9995 )nc, sname,
973  $ uplo, n, alpha, incx, beta, incy
974  IF( rewi )
975  $ rewind ntra
976  CALL chpmv( uplo, n, alpha, aa, xx, incx,
977  $ beta, yy, incy )
978  END IF
979 *
980 * Check if error-exit was taken incorrectly.
981 *
982  IF( .NOT.ok )THEN
983  WRITE( nout, fmt = 9992 )
984  fatal = .true.
985  go to 120
986  END IF
987 *
988 * See what data changed inside subroutines.
989 *
990  isame( 1 ) = uplo.EQ.uplos
991  isame( 2 ) = ns.EQ.n
992  IF( full )THEN
993  isame( 3 ) = als.EQ.alpha
994  isame( 4 ) = lce( as, aa, laa )
995  isame( 5 ) = ldas.EQ.lda
996  isame( 6 ) = lce( xs, xx, lx )
997  isame( 7 ) = incxs.EQ.incx
998  isame( 8 ) = bls.EQ.beta
999  IF( null )THEN
1000  isame( 9 ) = lce( ys, yy, ly )
1001  ELSE
1002  isame( 9 ) = lceres( 'GE', ' ', 1, n,
1003  $ ys, yy, abs( incy ) )
1004  END IF
1005  isame( 10 ) = incys.EQ.incy
1006  ELSE IF( banded )THEN
1007  isame( 3 ) = ks.EQ.k
1008  isame( 4 ) = als.EQ.alpha
1009  isame( 5 ) = lce( as, aa, laa )
1010  isame( 6 ) = ldas.EQ.lda
1011  isame( 7 ) = lce( xs, xx, lx )
1012  isame( 8 ) = incxs.EQ.incx
1013  isame( 9 ) = bls.EQ.beta
1014  IF( null )THEN
1015  isame( 10 ) = lce( ys, yy, ly )
1016  ELSE
1017  isame( 10 ) = lceres( 'GE', ' ', 1, n,
1018  $ ys, yy, abs( incy ) )
1019  END IF
1020  isame( 11 ) = incys.EQ.incy
1021  ELSE IF( packed )THEN
1022  isame( 3 ) = als.EQ.alpha
1023  isame( 4 ) = lce( as, aa, laa )
1024  isame( 5 ) = lce( xs, xx, lx )
1025  isame( 6 ) = incxs.EQ.incx
1026  isame( 7 ) = bls.EQ.beta
1027  IF( null )THEN
1028  isame( 8 ) = lce( ys, yy, ly )
1029  ELSE
1030  isame( 8 ) = lceres( 'GE', ' ', 1, n,
1031  $ ys, yy, abs( incy ) )
1032  END IF
1033  isame( 9 ) = incys.EQ.incy
1034  END IF
1035 *
1036 * If data was incorrectly changed, report and
1037 * return.
1038 *
1039  same = .true.
1040  DO 40 i = 1, nargs
1041  same = same.AND.isame( i )
1042  IF( .NOT.isame( i ) )
1043  $ WRITE( nout, fmt = 9998 )i
1044  40 continue
1045  IF( .NOT.same )THEN
1046  fatal = .true.
1047  go to 120
1048  END IF
1049 *
1050  IF( .NOT.null )THEN
1051 *
1052 * Check the result.
1053 *
1054  CALL cmvch( 'N', n, n, alpha, a, nmax, x,
1055  $ incx, beta, y, incy, yt, g,
1056  $ yy, eps, err, fatal, nout,
1057  $ .true. )
1058  errmax = max( errmax, err )
1059 * If got really bad answer, report and
1060 * return.
1061  IF( fatal )
1062  $ go to 120
1063  ELSE
1064 * Avoid repeating tests with N.le.0
1065  go to 110
1066  END IF
1067 *
1068  50 continue
1069 *
1070  60 continue
1071 *
1072  70 continue
1073 *
1074  80 continue
1075 *
1076  90 continue
1077 *
1078  100 continue
1079 *
1080  110 continue
1081 *
1082 * Report result.
1083 *
1084  IF( errmax.LT.thresh )THEN
1085  WRITE( nout, fmt = 9999 )sname, nc
1086  ELSE
1087  WRITE( nout, fmt = 9997 )sname, nc, errmax
1088  END IF
1089  go to 130
1090 *
1091  120 continue
1092  WRITE( nout, fmt = 9996 )sname
1093  IF( full )THEN
1094  WRITE( nout, fmt = 9993 )nc, sname, uplo, n, alpha, lda, incx,
1095  $ beta, incy
1096  ELSE IF( banded )THEN
1097  WRITE( nout, fmt = 9994 )nc, sname, uplo, n, k, alpha, lda,
1098  $ incx, beta, incy
1099  ELSE IF( packed )THEN
1100  WRITE( nout, fmt = 9995 )nc, sname, uplo, n, alpha, incx,
1101  $ beta, incy
1102  END IF
1103 *
1104  130 continue
1105  return
1106 *
1107  9999 format( ' ', a6, ' PASSED THE COMPUTATIONAL TESTS (', i6, ' CALL',
1108  $ 'S)' )
1109  9998 format( ' ******* FATAL ERROR - PARAMETER NUMBER ', i2, ' WAS CH',
1110  $ 'ANGED INCORRECTLY *******' )
1111  9997 format( ' ', a6, ' COMPLETED THE COMPUTATIONAL TESTS (', i6, ' C',
1112  $ 'ALLS)', /' ******* BUT WITH MAXIMUM TEST RATIO', f8.2,
1113  $ ' - SUSPECT *******' )
1114  9996 format( ' ******* ', a6, ' FAILED ON CALL NUMBER:' )
1115  9995 format( 1x, i6, ': ', a6, '(''', a1, ''',', i3, ',(', f4.1, ',',
1116  $ f4.1, '), AP, X,', i2, ',(', f4.1, ',', f4.1, '), Y,', i2,
1117  $ ') .' )
1118  9994 format( 1x, i6, ': ', a6, '(''', a1, ''',', 2( i3, ',' ), '(',
1119  $ f4.1, ',', f4.1, '), A,', i3, ', X,', i2, ',(', f4.1, ',',
1120  $ f4.1, '), Y,', i2, ') .' )
1121  9993 format( 1x, i6, ': ', a6, '(''', a1, ''',', i3, ',(', f4.1, ',',
1122  $ f4.1, '), A,', i3, ', X,', i2, ',(', f4.1, ',', f4.1, '), ',
1123  $ 'Y,', i2, ') .' )
1124  9992 format( ' ******* FATAL ERROR - ERROR-EXIT TAKEN ON VALID CALL *',
1125  $ '******' )
1126 *
1127 * End of CCHK2.
1128 *
1129  END
1130  SUBROUTINE cchk3( SNAME, EPS, THRESH, NOUT, NTRA, TRACE, REWI,
1131  $ fatal, nidim, idim, nkb, kb, ninc, inc, nmax,
1132  $ incmax, a, aa, as, x, xx, xs, xt, g, z )
1133 *
1134 * Tests CTRMV, CTBMV, CTPMV, CTRSV, CTBSV and CTPSV.
1135 *
1136 * Auxiliary routine for test program for Level 2 Blas.
1137 *
1138 * -- Written on 10-August-1987.
1139 * Richard Hanson, Sandia National Labs.
1140 * Jeremy Du Croz, NAG Central Office.
1141 *
1142 * .. Parameters ..
1143  COMPLEX zero, half, one
1144  parameter( zero = ( 0.0, 0.0 ), half = ( 0.5, 0.0 ),
1145  $ one = ( 1.0, 0.0 ) )
1146  REAL rzero
1147  parameter( rzero = 0.0 )
1148 * .. Scalar Arguments ..
1149  REAL eps, thresh
1150  INTEGER incmax, nidim, ninc, nkb, nmax, nout, ntra
1151  LOGICAL fatal, rewi, trace
1152  CHARACTER*6 sname
1153 * .. Array Arguments ..
1154  COMPLEX a( nmax, nmax ), aa( nmax*nmax ),
1155  $ as( nmax*nmax ), x( nmax ), xs( nmax*incmax ),
1156  $ xt( nmax ), xx( nmax*incmax ), z( nmax )
1157  REAL g( nmax )
1158  INTEGER idim( nidim ), inc( ninc ), kb( nkb )
1159 * .. Local Scalars ..
1160  COMPLEX transl
1161  REAL err, errmax
1162  INTEGER i, icd, ict, icu, ik, in, incx, incxs, ix, k,
1163  $ ks, laa, lda, ldas, lx, n, nargs, nc, nk, ns
1164  LOGICAL banded, full, null, packed, reset, same
1165  CHARACTER*1 diag, diags, trans, transs, uplo, uplos
1166  CHARACTER*2 ichd, ichu
1167  CHARACTER*3 icht
1168 * .. Local Arrays ..
1169  LOGICAL isame( 13 )
1170 * .. External Functions ..
1171  LOGICAL lce, lceres
1172  EXTERNAL lce, lceres
1173 * .. External Subroutines ..
1174  EXTERNAL cmake, cmvch, ctbmv, ctbsv, ctpmv, ctpsv,
1175  $ ctrmv, ctrsv
1176 * .. Intrinsic Functions ..
1177  INTRINSIC abs, max
1178 * .. Scalars in Common ..
1179  INTEGER infot, noutc
1180  LOGICAL lerr, ok
1181 * .. Common blocks ..
1182  common /infoc/infot, noutc, ok, lerr
1183 * .. Data statements ..
1184  DATA ichu/'UL'/, icht/'NTC'/, ichd/'UN'/
1185 * .. Executable Statements ..
1186  full = sname( 3: 3 ).EQ.'R'
1187  banded = sname( 3: 3 ).EQ.'B'
1188  packed = sname( 3: 3 ).EQ.'P'
1189 * Define the number of arguments.
1190  IF( full )THEN
1191  nargs = 8
1192  ELSE IF( banded )THEN
1193  nargs = 9
1194  ELSE IF( packed )THEN
1195  nargs = 7
1196  END IF
1197 *
1198  nc = 0
1199  reset = .true.
1200  errmax = rzero
1201 * Set up zero vector for CMVCH.
1202  DO 10 i = 1, nmax
1203  z( i ) = zero
1204  10 continue
1205 *
1206  DO 110 in = 1, nidim
1207  n = idim( in )
1208 *
1209  IF( banded )THEN
1210  nk = nkb
1211  ELSE
1212  nk = 1
1213  END IF
1214  DO 100 ik = 1, nk
1215  IF( banded )THEN
1216  k = kb( ik )
1217  ELSE
1218  k = n - 1
1219  END IF
1220 * Set LDA to 1 more than minimum value if room.
1221  IF( banded )THEN
1222  lda = k + 1
1223  ELSE
1224  lda = n
1225  END IF
1226  IF( lda.LT.nmax )
1227  $ lda = lda + 1
1228 * Skip tests if not enough room.
1229  IF( lda.GT.nmax )
1230  $ go to 100
1231  IF( packed )THEN
1232  laa = ( n*( n + 1 ) )/2
1233  ELSE
1234  laa = lda*n
1235  END IF
1236  null = n.LE.0
1237 *
1238  DO 90 icu = 1, 2
1239  uplo = ichu( icu: icu )
1240 *
1241  DO 80 ict = 1, 3
1242  trans = icht( ict: ict )
1243 *
1244  DO 70 icd = 1, 2
1245  diag = ichd( icd: icd )
1246 *
1247 * Generate the matrix A.
1248 *
1249  transl = zero
1250  CALL cmake( sname( 2: 3 ), uplo, diag, n, n, a,
1251  $ nmax, aa, lda, k, k, reset, transl )
1252 *
1253  DO 60 ix = 1, ninc
1254  incx = inc( ix )
1255  lx = abs( incx )*n
1256 *
1257 * Generate the vector X.
1258 *
1259  transl = half
1260  CALL cmake( 'GE', ' ', ' ', 1, n, x, 1, xx,
1261  $ abs( incx ), 0, n - 1, reset,
1262  $ transl )
1263  IF( n.GT.1 )THEN
1264  x( n/2 ) = zero
1265  xx( 1 + abs( incx )*( n/2 - 1 ) ) = zero
1266  END IF
1267 *
1268  nc = nc + 1
1269 *
1270 * Save every datum before calling the subroutine.
1271 *
1272  uplos = uplo
1273  transs = trans
1274  diags = diag
1275  ns = n
1276  ks = k
1277  DO 20 i = 1, laa
1278  as( i ) = aa( i )
1279  20 continue
1280  ldas = lda
1281  DO 30 i = 1, lx
1282  xs( i ) = xx( i )
1283  30 continue
1284  incxs = incx
1285 *
1286 * Call the subroutine.
1287 *
1288  IF( sname( 4: 5 ).EQ.'MV' )THEN
1289  IF( full )THEN
1290  IF( trace )
1291  $ WRITE( ntra, fmt = 9993 )nc, sname,
1292  $ uplo, trans, diag, n, lda, incx
1293  IF( rewi )
1294  $ rewind ntra
1295  CALL ctrmv( uplo, trans, diag, n, aa, lda,
1296  $ xx, incx )
1297  ELSE IF( banded )THEN
1298  IF( trace )
1299  $ WRITE( ntra, fmt = 9994 )nc, sname,
1300  $ uplo, trans, diag, n, k, lda, incx
1301  IF( rewi )
1302  $ rewind ntra
1303  CALL ctbmv( uplo, trans, diag, n, k, aa,
1304  $ lda, xx, incx )
1305  ELSE IF( packed )THEN
1306  IF( trace )
1307  $ WRITE( ntra, fmt = 9995 )nc, sname,
1308  $ uplo, trans, diag, n, incx
1309  IF( rewi )
1310  $ rewind ntra
1311  CALL ctpmv( uplo, trans, diag, n, aa, xx,
1312  $ incx )
1313  END IF
1314  ELSE IF( sname( 4: 5 ).EQ.'SV' )THEN
1315  IF( full )THEN
1316  IF( trace )
1317  $ WRITE( ntra, fmt = 9993 )nc, sname,
1318  $ uplo, trans, diag, n, lda, incx
1319  IF( rewi )
1320  $ rewind ntra
1321  CALL ctrsv( uplo, trans, diag, n, aa, lda,
1322  $ xx, incx )
1323  ELSE IF( banded )THEN
1324  IF( trace )
1325  $ WRITE( ntra, fmt = 9994 )nc, sname,
1326  $ uplo, trans, diag, n, k, lda, incx
1327  IF( rewi )
1328  $ rewind ntra
1329  CALL ctbsv( uplo, trans, diag, n, k, aa,
1330  $ lda, xx, incx )
1331  ELSE IF( packed )THEN
1332  IF( trace )
1333  $ WRITE( ntra, fmt = 9995 )nc, sname,
1334  $ uplo, trans, diag, n, incx
1335  IF( rewi )
1336  $ rewind ntra
1337  CALL ctpsv( uplo, trans, diag, n, aa, xx,
1338  $ incx )
1339  END IF
1340  END IF
1341 *
1342 * Check if error-exit was taken incorrectly.
1343 *
1344  IF( .NOT.ok )THEN
1345  WRITE( nout, fmt = 9992 )
1346  fatal = .true.
1347  go to 120
1348  END IF
1349 *
1350 * See what data changed inside subroutines.
1351 *
1352  isame( 1 ) = uplo.EQ.uplos
1353  isame( 2 ) = trans.EQ.transs
1354  isame( 3 ) = diag.EQ.diags
1355  isame( 4 ) = ns.EQ.n
1356  IF( full )THEN
1357  isame( 5 ) = lce( as, aa, laa )
1358  isame( 6 ) = ldas.EQ.lda
1359  IF( null )THEN
1360  isame( 7 ) = lce( xs, xx, lx )
1361  ELSE
1362  isame( 7 ) = lceres( 'GE', ' ', 1, n, xs,
1363  $ xx, abs( incx ) )
1364  END IF
1365  isame( 8 ) = incxs.EQ.incx
1366  ELSE IF( banded )THEN
1367  isame( 5 ) = ks.EQ.k
1368  isame( 6 ) = lce( as, aa, laa )
1369  isame( 7 ) = ldas.EQ.lda
1370  IF( null )THEN
1371  isame( 8 ) = lce( xs, xx, lx )
1372  ELSE
1373  isame( 8 ) = lceres( 'GE', ' ', 1, n, xs,
1374  $ xx, abs( incx ) )
1375  END IF
1376  isame( 9 ) = incxs.EQ.incx
1377  ELSE IF( packed )THEN
1378  isame( 5 ) = lce( as, aa, laa )
1379  IF( null )THEN
1380  isame( 6 ) = lce( xs, xx, lx )
1381  ELSE
1382  isame( 6 ) = lceres( 'GE', ' ', 1, n, xs,
1383  $ xx, abs( incx ) )
1384  END IF
1385  isame( 7 ) = incxs.EQ.incx
1386  END IF
1387 *
1388 * If data was incorrectly changed, report and
1389 * return.
1390 *
1391  same = .true.
1392  DO 40 i = 1, nargs
1393  same = same.AND.isame( i )
1394  IF( .NOT.isame( i ) )
1395  $ WRITE( nout, fmt = 9998 )i
1396  40 continue
1397  IF( .NOT.same )THEN
1398  fatal = .true.
1399  go to 120
1400  END IF
1401 *
1402  IF( .NOT.null )THEN
1403  IF( sname( 4: 5 ).EQ.'MV' )THEN
1404 *
1405 * Check the result.
1406 *
1407  CALL cmvch( trans, n, n, one, a, nmax, x,
1408  $ incx, zero, z, incx, xt, g,
1409  $ xx, eps, err, fatal, nout,
1410  $ .true. )
1411  ELSE IF( sname( 4: 5 ).EQ.'SV' )THEN
1412 *
1413 * Compute approximation to original vector.
1414 *
1415  DO 50 i = 1, n
1416  z( i ) = xx( 1 + ( i - 1 )*
1417  $ abs( incx ) )
1418  xx( 1 + ( i - 1 )*abs( incx ) )
1419  $ = x( i )
1420  50 continue
1421  CALL cmvch( trans, n, n, one, a, nmax, z,
1422  $ incx, zero, x, incx, xt, g,
1423  $ xx, eps, err, fatal, nout,
1424  $ .false. )
1425  END IF
1426  errmax = max( errmax, err )
1427 * If got really bad answer, report and return.
1428  IF( fatal )
1429  $ go to 120
1430  ELSE
1431 * Avoid repeating tests with N.le.0.
1432  go to 110
1433  END IF
1434 *
1435  60 continue
1436 *
1437  70 continue
1438 *
1439  80 continue
1440 *
1441  90 continue
1442 *
1443  100 continue
1444 *
1445  110 continue
1446 *
1447 * Report result.
1448 *
1449  IF( errmax.LT.thresh )THEN
1450  WRITE( nout, fmt = 9999 )sname, nc
1451  ELSE
1452  WRITE( nout, fmt = 9997 )sname, nc, errmax
1453  END IF
1454  go to 130
1455 *
1456  120 continue
1457  WRITE( nout, fmt = 9996 )sname
1458  IF( full )THEN
1459  WRITE( nout, fmt = 9993 )nc, sname, uplo, trans, diag, n, lda,
1460  $ incx
1461  ELSE IF( banded )THEN
1462  WRITE( nout, fmt = 9994 )nc, sname, uplo, trans, diag, n, k,
1463  $ lda, incx
1464  ELSE IF( packed )THEN
1465  WRITE( nout, fmt = 9995 )nc, sname, uplo, trans, diag, n, incx
1466  END IF
1467 *
1468  130 continue
1469  return
1470 *
1471  9999 format( ' ', a6, ' PASSED THE COMPUTATIONAL TESTS (', i6, ' CALL',
1472  $ 'S)' )
1473  9998 format( ' ******* FATAL ERROR - PARAMETER NUMBER ', i2, ' WAS CH',
1474  $ 'ANGED INCORRECTLY *******' )
1475  9997 format( ' ', a6, ' COMPLETED THE COMPUTATIONAL TESTS (', i6, ' C',
1476  $ 'ALLS)', /' ******* BUT WITH MAXIMUM TEST RATIO', f8.2,
1477  $ ' - SUSPECT *******' )
1478  9996 format( ' ******* ', a6, ' FAILED ON CALL NUMBER:' )
1479  9995 format( 1x, i6, ': ', a6, '(', 3( '''', a1, ''',' ), i3, ', AP, ',
1480  $ 'X,', i2, ') .' )
1481  9994 format( 1x, i6, ': ', a6, '(', 3( '''', a1, ''',' ), 2( i3, ',' ),
1482  $ ' A,', i3, ', X,', i2, ') .' )
1483  9993 format( 1x, i6, ': ', a6, '(', 3( '''', a1, ''',' ), i3, ', A,',
1484  $ i3, ', X,', i2, ') .' )
1485  9992 format( ' ******* FATAL ERROR - ERROR-EXIT TAKEN ON VALID CALL *',
1486  $ '******' )
1487 *
1488 * End of CCHK3.
1489 *
1490  END
1491  SUBROUTINE cchk4( SNAME, EPS, THRESH, NOUT, NTRA, TRACE, REWI,
1492  $ fatal, nidim, idim, nalf, alf, ninc, inc, nmax,
1493  $ incmax, a, aa, as, x, xx, xs, y, yy, ys, yt, g,
1494  $ z )
1495 *
1496 * Tests CGERC and CGERU.
1497 *
1498 * Auxiliary routine for test program for Level 2 Blas.
1499 *
1500 * -- Written on 10-August-1987.
1501 * Richard Hanson, Sandia National Labs.
1502 * Jeremy Du Croz, NAG Central Office.
1503 *
1504 * .. Parameters ..
1505  COMPLEX zero, half, one
1506  parameter( zero = ( 0.0, 0.0 ), half = ( 0.5, 0.0 ),
1507  $ one = ( 1.0, 0.0 ) )
1508  REAL rzero
1509  parameter( rzero = 0.0 )
1510 * .. Scalar Arguments ..
1511  REAL eps, thresh
1512  INTEGER incmax, nalf, nidim, ninc, nmax, nout, ntra
1513  LOGICAL fatal, rewi, trace
1514  CHARACTER*6 sname
1515 * .. Array Arguments ..
1516  COMPLEX a( nmax, nmax ), aa( nmax*nmax ), alf( nalf ),
1517  $ as( nmax*nmax ), x( nmax ), xs( nmax*incmax ),
1518  $ xx( nmax*incmax ), y( nmax ),
1519  $ ys( nmax*incmax ), yt( nmax ),
1520  $ yy( nmax*incmax ), z( nmax )
1521  REAL g( nmax )
1522  INTEGER idim( nidim ), inc( ninc )
1523 * .. Local Scalars ..
1524  COMPLEX alpha, als, transl
1525  REAL err, errmax
1526  INTEGER i, ia, im, in, incx, incxs, incy, incys, ix,
1527  $ iy, j, laa, lda, ldas, lx, ly, m, ms, n, nargs,
1528  $ nc, nd, ns
1529  LOGICAL conj, null, reset, same
1530 * .. Local Arrays ..
1531  COMPLEX w( 1 )
1532  LOGICAL isame( 13 )
1533 * .. External Functions ..
1534  LOGICAL lce, lceres
1535  EXTERNAL lce, lceres
1536 * .. External Subroutines ..
1537  EXTERNAL cgerc, cgeru, cmake, cmvch
1538 * .. Intrinsic Functions ..
1539  INTRINSIC abs, conjg, max, min
1540 * .. Scalars in Common ..
1541  INTEGER infot, noutc
1542  LOGICAL lerr, ok
1543 * .. Common blocks ..
1544  common /infoc/infot, noutc, ok, lerr
1545 * .. Executable Statements ..
1546  conj = sname( 5: 5 ).EQ.'C'
1547 * Define the number of arguments.
1548  nargs = 9
1549 *
1550  nc = 0
1551  reset = .true.
1552  errmax = rzero
1553 *
1554  DO 120 in = 1, nidim
1555  n = idim( in )
1556  nd = n/2 + 1
1557 *
1558  DO 110 im = 1, 2
1559  IF( im.EQ.1 )
1560  $ m = max( n - nd, 0 )
1561  IF( im.EQ.2 )
1562  $ m = min( n + nd, nmax )
1563 *
1564 * Set LDA to 1 more than minimum value if room.
1565  lda = m
1566  IF( lda.LT.nmax )
1567  $ lda = lda + 1
1568 * Skip tests if not enough room.
1569  IF( lda.GT.nmax )
1570  $ go to 110
1571  laa = lda*n
1572  null = n.LE.0.OR.m.LE.0
1573 *
1574  DO 100 ix = 1, ninc
1575  incx = inc( ix )
1576  lx = abs( incx )*m
1577 *
1578 * Generate the vector X.
1579 *
1580  transl = half
1581  CALL cmake( 'GE', ' ', ' ', 1, m, x, 1, xx, abs( incx ),
1582  $ 0, m - 1, reset, transl )
1583  IF( m.GT.1 )THEN
1584  x( m/2 ) = zero
1585  xx( 1 + abs( incx )*( m/2 - 1 ) ) = zero
1586  END IF
1587 *
1588  DO 90 iy = 1, ninc
1589  incy = inc( iy )
1590  ly = abs( incy )*n
1591 *
1592 * Generate the vector Y.
1593 *
1594  transl = zero
1595  CALL cmake( 'GE', ' ', ' ', 1, n, y, 1, yy,
1596  $ abs( incy ), 0, n - 1, reset, transl )
1597  IF( n.GT.1 )THEN
1598  y( n/2 ) = zero
1599  yy( 1 + abs( incy )*( n/2 - 1 ) ) = zero
1600  END IF
1601 *
1602  DO 80 ia = 1, nalf
1603  alpha = alf( ia )
1604 *
1605 * Generate the matrix A.
1606 *
1607  transl = zero
1608  CALL cmake( sname( 2: 3 ), ' ', ' ', m, n, a, nmax,
1609  $ aa, lda, m - 1, n - 1, reset, transl )
1610 *
1611  nc = nc + 1
1612 *
1613 * Save every datum before calling the subroutine.
1614 *
1615  ms = m
1616  ns = n
1617  als = alpha
1618  DO 10 i = 1, laa
1619  as( i ) = aa( i )
1620  10 continue
1621  ldas = lda
1622  DO 20 i = 1, lx
1623  xs( i ) = xx( i )
1624  20 continue
1625  incxs = incx
1626  DO 30 i = 1, ly
1627  ys( i ) = yy( i )
1628  30 continue
1629  incys = incy
1630 *
1631 * Call the subroutine.
1632 *
1633  IF( trace )
1634  $ WRITE( ntra, fmt = 9994 )nc, sname, m, n,
1635  $ alpha, incx, incy, lda
1636  IF( conj )THEN
1637  IF( rewi )
1638  $ rewind ntra
1639  CALL cgerc( m, n, alpha, xx, incx, yy, incy, aa,
1640  $ lda )
1641  ELSE
1642  IF( rewi )
1643  $ rewind ntra
1644  CALL cgeru( m, n, alpha, xx, incx, yy, incy, aa,
1645  $ lda )
1646  END IF
1647 *
1648 * Check if error-exit was taken incorrectly.
1649 *
1650  IF( .NOT.ok )THEN
1651  WRITE( nout, fmt = 9993 )
1652  fatal = .true.
1653  go to 140
1654  END IF
1655 *
1656 * See what data changed inside subroutine.
1657 *
1658  isame( 1 ) = ms.EQ.m
1659  isame( 2 ) = ns.EQ.n
1660  isame( 3 ) = als.EQ.alpha
1661  isame( 4 ) = lce( xs, xx, lx )
1662  isame( 5 ) = incxs.EQ.incx
1663  isame( 6 ) = lce( ys, yy, ly )
1664  isame( 7 ) = incys.EQ.incy
1665  IF( null )THEN
1666  isame( 8 ) = lce( as, aa, laa )
1667  ELSE
1668  isame( 8 ) = lceres( 'GE', ' ', m, n, as, aa,
1669  $ lda )
1670  END IF
1671  isame( 9 ) = ldas.EQ.lda
1672 *
1673 * If data was incorrectly changed, report and return.
1674 *
1675  same = .true.
1676  DO 40 i = 1, nargs
1677  same = same.AND.isame( i )
1678  IF( .NOT.isame( i ) )
1679  $ WRITE( nout, fmt = 9998 )i
1680  40 continue
1681  IF( .NOT.same )THEN
1682  fatal = .true.
1683  go to 140
1684  END IF
1685 *
1686  IF( .NOT.null )THEN
1687 *
1688 * Check the result column by column.
1689 *
1690  IF( incx.GT.0 )THEN
1691  DO 50 i = 1, m
1692  z( i ) = x( i )
1693  50 continue
1694  ELSE
1695  DO 60 i = 1, m
1696  z( i ) = x( m - i + 1 )
1697  60 continue
1698  END IF
1699  DO 70 j = 1, n
1700  IF( incy.GT.0 )THEN
1701  w( 1 ) = y( j )
1702  ELSE
1703  w( 1 ) = y( n - j + 1 )
1704  END IF
1705  IF( conj )
1706  $ w( 1 ) = conjg( w( 1 ) )
1707  CALL cmvch( 'N', m, 1, alpha, z, nmax, w, 1,
1708  $ one, a( 1, j ), 1, yt, g,
1709  $ aa( 1 + ( j - 1 )*lda ), eps,
1710  $ err, fatal, nout, .true. )
1711  errmax = max( errmax, err )
1712 * If got really bad answer, report and return.
1713  IF( fatal )
1714  $ go to 130
1715  70 continue
1716  ELSE
1717 * Avoid repeating tests with M.le.0 or N.le.0.
1718  go to 110
1719  END IF
1720 *
1721  80 continue
1722 *
1723  90 continue
1724 *
1725  100 continue
1726 *
1727  110 continue
1728 *
1729  120 continue
1730 *
1731 * Report result.
1732 *
1733  IF( errmax.LT.thresh )THEN
1734  WRITE( nout, fmt = 9999 )sname, nc
1735  ELSE
1736  WRITE( nout, fmt = 9997 )sname, nc, errmax
1737  END IF
1738  go to 150
1739 *
1740  130 continue
1741  WRITE( nout, fmt = 9995 )j
1742 *
1743  140 continue
1744  WRITE( nout, fmt = 9996 )sname
1745  WRITE( nout, fmt = 9994 )nc, sname, m, n, alpha, incx, incy, lda
1746 *
1747  150 continue
1748  return
1749 *
1750  9999 format( ' ', a6, ' PASSED THE COMPUTATIONAL TESTS (', i6, ' CALL',
1751  $ 'S)' )
1752  9998 format( ' ******* FATAL ERROR - PARAMETER NUMBER ', i2, ' WAS CH',
1753  $ 'ANGED INCORRECTLY *******' )
1754  9997 format( ' ', a6, ' COMPLETED THE COMPUTATIONAL TESTS (', i6, ' C',
1755  $ 'ALLS)', /' ******* BUT WITH MAXIMUM TEST RATIO', f8.2,
1756  $ ' - SUSPECT *******' )
1757  9996 format( ' ******* ', a6, ' FAILED ON CALL NUMBER:' )
1758  9995 format( ' THESE ARE THE RESULTS FOR COLUMN ', i3 )
1759  9994 format( 1x, i6, ': ', a6, '(', 2( i3, ',' ), '(', f4.1, ',', f4.1,
1760  $ '), X,', i2, ', Y,', i2, ', A,', i3, ') ',
1761  $ ' .' )
1762  9993 format( ' ******* FATAL ERROR - ERROR-EXIT TAKEN ON VALID CALL *',
1763  $ '******' )
1764 *
1765 * End of CCHK4.
1766 *
1767  END
1768  SUBROUTINE cchk5( SNAME, EPS, THRESH, NOUT, NTRA, TRACE, REWI,
1769  $ fatal, nidim, idim, nalf, alf, ninc, inc, nmax,
1770  $ incmax, a, aa, as, x, xx, xs, y, yy, ys, yt, g,
1771  $ z )
1772 *
1773 * Tests CHER and CHPR.
1774 *
1775 * Auxiliary routine for test program for Level 2 Blas.
1776 *
1777 * -- Written on 10-August-1987.
1778 * Richard Hanson, Sandia National Labs.
1779 * Jeremy Du Croz, NAG Central Office.
1780 *
1781 * .. Parameters ..
1782  COMPLEX zero, half, one
1783  parameter( zero = ( 0.0, 0.0 ), half = ( 0.5, 0.0 ),
1784  $ one = ( 1.0, 0.0 ) )
1785  REAL rzero
1786  parameter( rzero = 0.0 )
1787 * .. Scalar Arguments ..
1788  REAL eps, thresh
1789  INTEGER incmax, nalf, nidim, ninc, nmax, nout, ntra
1790  LOGICAL fatal, rewi, trace
1791  CHARACTER*6 sname
1792 * .. Array Arguments ..
1793  COMPLEX a( nmax, nmax ), aa( nmax*nmax ), alf( nalf ),
1794  $ as( nmax*nmax ), x( nmax ), xs( nmax*incmax ),
1795  $ xx( nmax*incmax ), y( nmax ),
1796  $ ys( nmax*incmax ), yt( nmax ),
1797  $ yy( nmax*incmax ), z( nmax )
1798  REAL g( nmax )
1799  INTEGER idim( nidim ), inc( ninc )
1800 * .. Local Scalars ..
1801  COMPLEX alpha, transl
1802  REAL err, errmax, ralpha, rals
1803  INTEGER i, ia, ic, in, incx, incxs, ix, j, ja, jj, laa,
1804  $ lda, ldas, lj, lx, n, nargs, nc, ns
1805  LOGICAL full, null, packed, reset, same, upper
1806  CHARACTER*1 uplo, uplos
1807  CHARACTER*2 ich
1808 * .. Local Arrays ..
1809  COMPLEX w( 1 )
1810  LOGICAL isame( 13 )
1811 * .. External Functions ..
1812  LOGICAL lce, lceres
1813  EXTERNAL lce, lceres
1814 * .. External Subroutines ..
1815  EXTERNAL cher, chpr, cmake, cmvch
1816 * .. Intrinsic Functions ..
1817  INTRINSIC abs, cmplx, conjg, max, real
1818 * .. Scalars in Common ..
1819  INTEGER infot, noutc
1820  LOGICAL lerr, ok
1821 * .. Common blocks ..
1822  common /infoc/infot, noutc, ok, lerr
1823 * .. Data statements ..
1824  DATA ich/'UL'/
1825 * .. Executable Statements ..
1826  full = sname( 3: 3 ).EQ.'E'
1827  packed = sname( 3: 3 ).EQ.'P'
1828 * Define the number of arguments.
1829  IF( full )THEN
1830  nargs = 7
1831  ELSE IF( packed )THEN
1832  nargs = 6
1833  END IF
1834 *
1835  nc = 0
1836  reset = .true.
1837  errmax = rzero
1838 *
1839  DO 100 in = 1, nidim
1840  n = idim( in )
1841 * Set LDA to 1 more than minimum value if room.
1842  lda = n
1843  IF( lda.LT.nmax )
1844  $ lda = lda + 1
1845 * Skip tests if not enough room.
1846  IF( lda.GT.nmax )
1847  $ go to 100
1848  IF( packed )THEN
1849  laa = ( n*( n + 1 ) )/2
1850  ELSE
1851  laa = lda*n
1852  END IF
1853 *
1854  DO 90 ic = 1, 2
1855  uplo = ich( ic: ic )
1856  upper = uplo.EQ.'U'
1857 *
1858  DO 80 ix = 1, ninc
1859  incx = inc( ix )
1860  lx = abs( incx )*n
1861 *
1862 * Generate the vector X.
1863 *
1864  transl = half
1865  CALL cmake( 'GE', ' ', ' ', 1, n, x, 1, xx, abs( incx ),
1866  $ 0, n - 1, reset, transl )
1867  IF( n.GT.1 )THEN
1868  x( n/2 ) = zero
1869  xx( 1 + abs( incx )*( n/2 - 1 ) ) = zero
1870  END IF
1871 *
1872  DO 70 ia = 1, nalf
1873  ralpha = REAL( ALF( IA ) )
1874  alpha = cmplx( ralpha, rzero )
1875  null = n.LE.0.OR.ralpha.EQ.rzero
1876 *
1877 * Generate the matrix A.
1878 *
1879  transl = zero
1880  CALL cmake( sname( 2: 3 ), uplo, ' ', n, n, a, nmax,
1881  $ aa, lda, n - 1, n - 1, reset, transl )
1882 *
1883  nc = nc + 1
1884 *
1885 * Save every datum before calling the subroutine.
1886 *
1887  uplos = uplo
1888  ns = n
1889  rals = ralpha
1890  DO 10 i = 1, laa
1891  as( i ) = aa( i )
1892  10 continue
1893  ldas = lda
1894  DO 20 i = 1, lx
1895  xs( i ) = xx( i )
1896  20 continue
1897  incxs = incx
1898 *
1899 * Call the subroutine.
1900 *
1901  IF( full )THEN
1902  IF( trace )
1903  $ WRITE( ntra, fmt = 9993 )nc, sname, uplo, n,
1904  $ ralpha, incx, lda
1905  IF( rewi )
1906  $ rewind ntra
1907  CALL cher( uplo, n, ralpha, xx, incx, aa, lda )
1908  ELSE IF( packed )THEN
1909  IF( trace )
1910  $ WRITE( ntra, fmt = 9994 )nc, sname, uplo, n,
1911  $ ralpha, incx
1912  IF( rewi )
1913  $ rewind ntra
1914  CALL chpr( uplo, n, ralpha, xx, incx, aa )
1915  END IF
1916 *
1917 * Check if error-exit was taken incorrectly.
1918 *
1919  IF( .NOT.ok )THEN
1920  WRITE( nout, fmt = 9992 )
1921  fatal = .true.
1922  go to 120
1923  END IF
1924 *
1925 * See what data changed inside subroutines.
1926 *
1927  isame( 1 ) = uplo.EQ.uplos
1928  isame( 2 ) = ns.EQ.n
1929  isame( 3 ) = rals.EQ.ralpha
1930  isame( 4 ) = lce( xs, xx, lx )
1931  isame( 5 ) = incxs.EQ.incx
1932  IF( null )THEN
1933  isame( 6 ) = lce( as, aa, laa )
1934  ELSE
1935  isame( 6 ) = lceres( sname( 2: 3 ), uplo, n, n, as,
1936  $ aa, lda )
1937  END IF
1938  IF( .NOT.packed )THEN
1939  isame( 7 ) = ldas.EQ.lda
1940  END IF
1941 *
1942 * If data was incorrectly changed, report and return.
1943 *
1944  same = .true.
1945  DO 30 i = 1, nargs
1946  same = same.AND.isame( i )
1947  IF( .NOT.isame( i ) )
1948  $ WRITE( nout, fmt = 9998 )i
1949  30 continue
1950  IF( .NOT.same )THEN
1951  fatal = .true.
1952  go to 120
1953  END IF
1954 *
1955  IF( .NOT.null )THEN
1956 *
1957 * Check the result column by column.
1958 *
1959  IF( incx.GT.0 )THEN
1960  DO 40 i = 1, n
1961  z( i ) = x( i )
1962  40 continue
1963  ELSE
1964  DO 50 i = 1, n
1965  z( i ) = x( n - i + 1 )
1966  50 continue
1967  END IF
1968  ja = 1
1969  DO 60 j = 1, n
1970  w( 1 ) = conjg( z( j ) )
1971  IF( upper )THEN
1972  jj = 1
1973  lj = j
1974  ELSE
1975  jj = j
1976  lj = n - j + 1
1977  END IF
1978  CALL cmvch( 'N', lj, 1, alpha, z( jj ), lj, w,
1979  $ 1, one, a( jj, j ), 1, yt, g,
1980  $ aa( ja ), eps, err, fatal, nout,
1981  $ .true. )
1982  IF( full )THEN
1983  IF( upper )THEN
1984  ja = ja + lda
1985  ELSE
1986  ja = ja + lda + 1
1987  END IF
1988  ELSE
1989  ja = ja + lj
1990  END IF
1991  errmax = max( errmax, err )
1992 * If got really bad answer, report and return.
1993  IF( fatal )
1994  $ go to 110
1995  60 continue
1996  ELSE
1997 * Avoid repeating tests if N.le.0.
1998  IF( n.LE.0 )
1999  $ go to 100
2000  END IF
2001 *
2002  70 continue
2003 *
2004  80 continue
2005 *
2006  90 continue
2007 *
2008  100 continue
2009 *
2010 * Report result.
2011 *
2012  IF( errmax.LT.thresh )THEN
2013  WRITE( nout, fmt = 9999 )sname, nc
2014  ELSE
2015  WRITE( nout, fmt = 9997 )sname, nc, errmax
2016  END IF
2017  go to 130
2018 *
2019  110 continue
2020  WRITE( nout, fmt = 9995 )j
2021 *
2022  120 continue
2023  WRITE( nout, fmt = 9996 )sname
2024  IF( full )THEN
2025  WRITE( nout, fmt = 9993 )nc, sname, uplo, n, ralpha, incx, lda
2026  ELSE IF( packed )THEN
2027  WRITE( nout, fmt = 9994 )nc, sname, uplo, n, ralpha, incx
2028  END IF
2029 *
2030  130 continue
2031  return
2032 *
2033  9999 format( ' ', a6, ' PASSED THE COMPUTATIONAL TESTS (', i6, ' CALL',
2034  $ 'S)' )
2035  9998 format( ' ******* FATAL ERROR - PARAMETER NUMBER ', i2, ' WAS CH',
2036  $ 'ANGED INCORRECTLY *******' )
2037  9997 format( ' ', a6, ' COMPLETED THE COMPUTATIONAL TESTS (', i6, ' C',
2038  $ 'ALLS)', /' ******* BUT WITH MAXIMUM TEST RATIO', f8.2,
2039  $ ' - SUSPECT *******' )
2040  9996 format( ' ******* ', a6, ' FAILED ON CALL NUMBER:' )
2041  9995 format( ' THESE ARE THE RESULTS FOR COLUMN ', i3 )
2042  9994 format( 1x, i6, ': ', a6, '(''', a1, ''',', i3, ',', f4.1, ', X,',
2043  $ i2, ', AP) .' )
2044  9993 format( 1x, i6, ': ', a6, '(''', a1, ''',', i3, ',', f4.1, ', X,',
2045  $ i2, ', A,', i3, ') .' )
2046  9992 format( ' ******* FATAL ERROR - ERROR-EXIT TAKEN ON VALID CALL *',
2047  $ '******' )
2048 *
2049 * End of CCHK5.
2050 *
2051  END
2052  SUBROUTINE cchk6( SNAME, EPS, THRESH, NOUT, NTRA, TRACE, REWI,
2053  $ fatal, nidim, idim, nalf, alf, ninc, inc, nmax,
2054  $ incmax, a, aa, as, x, xx, xs, y, yy, ys, yt, g,
2055  $ z )
2056 *
2057 * Tests CHER2 and CHPR2.
2058 *
2059 * Auxiliary routine for test program for Level 2 Blas.
2060 *
2061 * -- Written on 10-August-1987.
2062 * Richard Hanson, Sandia National Labs.
2063 * Jeremy Du Croz, NAG Central Office.
2064 *
2065 * .. Parameters ..
2066  COMPLEX zero, half, one
2067  parameter( zero = ( 0.0, 0.0 ), half = ( 0.5, 0.0 ),
2068  $ one = ( 1.0, 0.0 ) )
2069  REAL rzero
2070  parameter( rzero = 0.0 )
2071 * .. Scalar Arguments ..
2072  REAL eps, thresh
2073  INTEGER incmax, nalf, nidim, ninc, nmax, nout, ntra
2074  LOGICAL fatal, rewi, trace
2075  CHARACTER*6 sname
2076 * .. Array Arguments ..
2077  COMPLEX a( nmax, nmax ), aa( nmax*nmax ), alf( nalf ),
2078  $ as( nmax*nmax ), x( nmax ), xs( nmax*incmax ),
2079  $ xx( nmax*incmax ), y( nmax ),
2080  $ ys( nmax*incmax ), yt( nmax ),
2081  $ yy( nmax*incmax ), z( nmax, 2 )
2082  REAL g( nmax )
2083  INTEGER idim( nidim ), inc( ninc )
2084 * .. Local Scalars ..
2085  COMPLEX alpha, als, transl
2086  REAL err, errmax
2087  INTEGER i, ia, ic, in, incx, incxs, incy, incys, ix,
2088  $ iy, j, ja, jj, laa, lda, ldas, lj, lx, ly, n,
2089  $ nargs, nc, ns
2090  LOGICAL full, null, packed, reset, same, upper
2091  CHARACTER*1 uplo, uplos
2092  CHARACTER*2 ich
2093 * .. Local Arrays ..
2094  COMPLEX w( 2 )
2095  LOGICAL isame( 13 )
2096 * .. External Functions ..
2097  LOGICAL lce, lceres
2098  EXTERNAL lce, lceres
2099 * .. External Subroutines ..
2100  EXTERNAL cher2, chpr2, cmake, cmvch
2101 * .. Intrinsic Functions ..
2102  INTRINSIC abs, conjg, max
2103 * .. Scalars in Common ..
2104  INTEGER infot, noutc
2105  LOGICAL lerr, ok
2106 * .. Common blocks ..
2107  common /infoc/infot, noutc, ok, lerr
2108 * .. Data statements ..
2109  DATA ich/'UL'/
2110 * .. Executable Statements ..
2111  full = sname( 3: 3 ).EQ.'E'
2112  packed = sname( 3: 3 ).EQ.'P'
2113 * Define the number of arguments.
2114  IF( full )THEN
2115  nargs = 9
2116  ELSE IF( packed )THEN
2117  nargs = 8
2118  END IF
2119 *
2120  nc = 0
2121  reset = .true.
2122  errmax = rzero
2123 *
2124  DO 140 in = 1, nidim
2125  n = idim( in )
2126 * Set LDA to 1 more than minimum value if room.
2127  lda = n
2128  IF( lda.LT.nmax )
2129  $ lda = lda + 1
2130 * Skip tests if not enough room.
2131  IF( lda.GT.nmax )
2132  $ go to 140
2133  IF( packed )THEN
2134  laa = ( n*( n + 1 ) )/2
2135  ELSE
2136  laa = lda*n
2137  END IF
2138 *
2139  DO 130 ic = 1, 2
2140  uplo = ich( ic: ic )
2141  upper = uplo.EQ.'U'
2142 *
2143  DO 120 ix = 1, ninc
2144  incx = inc( ix )
2145  lx = abs( incx )*n
2146 *
2147 * Generate the vector X.
2148 *
2149  transl = half
2150  CALL cmake( 'GE', ' ', ' ', 1, n, x, 1, xx, abs( incx ),
2151  $ 0, n - 1, reset, transl )
2152  IF( n.GT.1 )THEN
2153  x( n/2 ) = zero
2154  xx( 1 + abs( incx )*( n/2 - 1 ) ) = zero
2155  END IF
2156 *
2157  DO 110 iy = 1, ninc
2158  incy = inc( iy )
2159  ly = abs( incy )*n
2160 *
2161 * Generate the vector Y.
2162 *
2163  transl = zero
2164  CALL cmake( 'GE', ' ', ' ', 1, n, y, 1, yy,
2165  $ abs( incy ), 0, n - 1, reset, transl )
2166  IF( n.GT.1 )THEN
2167  y( n/2 ) = zero
2168  yy( 1 + abs( incy )*( n/2 - 1 ) ) = zero
2169  END IF
2170 *
2171  DO 100 ia = 1, nalf
2172  alpha = alf( ia )
2173  null = n.LE.0.OR.alpha.EQ.zero
2174 *
2175 * Generate the matrix A.
2176 *
2177  transl = zero
2178  CALL cmake( sname( 2: 3 ), uplo, ' ', n, n, a,
2179  $ nmax, aa, lda, n - 1, n - 1, reset,
2180  $ transl )
2181 *
2182  nc = nc + 1
2183 *
2184 * Save every datum before calling the subroutine.
2185 *
2186  uplos = uplo
2187  ns = n
2188  als = alpha
2189  DO 10 i = 1, laa
2190  as( i ) = aa( i )
2191  10 continue
2192  ldas = lda
2193  DO 20 i = 1, lx
2194  xs( i ) = xx( i )
2195  20 continue
2196  incxs = incx
2197  DO 30 i = 1, ly
2198  ys( i ) = yy( i )
2199  30 continue
2200  incys = incy
2201 *
2202 * Call the subroutine.
2203 *
2204  IF( full )THEN
2205  IF( trace )
2206  $ WRITE( ntra, fmt = 9993 )nc, sname, uplo, n,
2207  $ alpha, incx, incy, lda
2208  IF( rewi )
2209  $ rewind ntra
2210  CALL cher2( uplo, n, alpha, xx, incx, yy, incy,
2211  $ aa, lda )
2212  ELSE IF( packed )THEN
2213  IF( trace )
2214  $ WRITE( ntra, fmt = 9994 )nc, sname, uplo, n,
2215  $ alpha, incx, incy
2216  IF( rewi )
2217  $ rewind ntra
2218  CALL chpr2( uplo, n, alpha, xx, incx, yy, incy,
2219  $ aa )
2220  END IF
2221 *
2222 * Check if error-exit was taken incorrectly.
2223 *
2224  IF( .NOT.ok )THEN
2225  WRITE( nout, fmt = 9992 )
2226  fatal = .true.
2227  go to 160
2228  END IF
2229 *
2230 * See what data changed inside subroutines.
2231 *
2232  isame( 1 ) = uplo.EQ.uplos
2233  isame( 2 ) = ns.EQ.n
2234  isame( 3 ) = als.EQ.alpha
2235  isame( 4 ) = lce( xs, xx, lx )
2236  isame( 5 ) = incxs.EQ.incx
2237  isame( 6 ) = lce( ys, yy, ly )
2238  isame( 7 ) = incys.EQ.incy
2239  IF( null )THEN
2240  isame( 8 ) = lce( as, aa, laa )
2241  ELSE
2242  isame( 8 ) = lceres( sname( 2: 3 ), uplo, n, n,
2243  $ as, aa, lda )
2244  END IF
2245  IF( .NOT.packed )THEN
2246  isame( 9 ) = ldas.EQ.lda
2247  END IF
2248 *
2249 * If data was incorrectly changed, report and return.
2250 *
2251  same = .true.
2252  DO 40 i = 1, nargs
2253  same = same.AND.isame( i )
2254  IF( .NOT.isame( i ) )
2255  $ WRITE( nout, fmt = 9998 )i
2256  40 continue
2257  IF( .NOT.same )THEN
2258  fatal = .true.
2259  go to 160
2260  END IF
2261 *
2262  IF( .NOT.null )THEN
2263 *
2264 * Check the result column by column.
2265 *
2266  IF( incx.GT.0 )THEN
2267  DO 50 i = 1, n
2268  z( i, 1 ) = x( i )
2269  50 continue
2270  ELSE
2271  DO 60 i = 1, n
2272  z( i, 1 ) = x( n - i + 1 )
2273  60 continue
2274  END IF
2275  IF( incy.GT.0 )THEN
2276  DO 70 i = 1, n
2277  z( i, 2 ) = y( i )
2278  70 continue
2279  ELSE
2280  DO 80 i = 1, n
2281  z( i, 2 ) = y( n - i + 1 )
2282  80 continue
2283  END IF
2284  ja = 1
2285  DO 90 j = 1, n
2286  w( 1 ) = alpha*conjg( z( j, 2 ) )
2287  w( 2 ) = conjg( alpha )*conjg( z( j, 1 ) )
2288  IF( upper )THEN
2289  jj = 1
2290  lj = j
2291  ELSE
2292  jj = j
2293  lj = n - j + 1
2294  END IF
2295  CALL cmvch( 'N', lj, 2, one, z( jj, 1 ),
2296  $ nmax, w, 1, one, a( jj, j ), 1,
2297  $ yt, g, aa( ja ), eps, err, fatal,
2298  $ nout, .true. )
2299  IF( full )THEN
2300  IF( upper )THEN
2301  ja = ja + lda
2302  ELSE
2303  ja = ja + lda + 1
2304  END IF
2305  ELSE
2306  ja = ja + lj
2307  END IF
2308  errmax = max( errmax, err )
2309 * If got really bad answer, report and return.
2310  IF( fatal )
2311  $ go to 150
2312  90 continue
2313  ELSE
2314 * Avoid repeating tests with N.le.0.
2315  IF( n.LE.0 )
2316  $ go to 140
2317  END IF
2318 *
2319  100 continue
2320 *
2321  110 continue
2322 *
2323  120 continue
2324 *
2325  130 continue
2326 *
2327  140 continue
2328 *
2329 * Report result.
2330 *
2331  IF( errmax.LT.thresh )THEN
2332  WRITE( nout, fmt = 9999 )sname, nc
2333  ELSE
2334  WRITE( nout, fmt = 9997 )sname, nc, errmax
2335  END IF
2336  go to 170
2337 *
2338  150 continue
2339  WRITE( nout, fmt = 9995 )j
2340 *
2341  160 continue
2342  WRITE( nout, fmt = 9996 )sname
2343  IF( full )THEN
2344  WRITE( nout, fmt = 9993 )nc, sname, uplo, n, alpha, incx,
2345  $ incy, lda
2346  ELSE IF( packed )THEN
2347  WRITE( nout, fmt = 9994 )nc, sname, uplo, n, alpha, incx, incy
2348  END IF
2349 *
2350  170 continue
2351  return
2352 *
2353  9999 format( ' ', a6, ' PASSED THE COMPUTATIONAL TESTS (', i6, ' CALL',
2354  $ 'S)' )
2355  9998 format( ' ******* FATAL ERROR - PARAMETER NUMBER ', i2, ' WAS CH',
2356  $ 'ANGED INCORRECTLY *******' )
2357  9997 format( ' ', a6, ' COMPLETED THE COMPUTATIONAL TESTS (', i6, ' C',
2358  $ 'ALLS)', /' ******* BUT WITH MAXIMUM TEST RATIO', f8.2,
2359  $ ' - SUSPECT *******' )
2360  9996 format( ' ******* ', a6, ' FAILED ON CALL NUMBER:' )
2361  9995 format( ' THESE ARE THE RESULTS FOR COLUMN ', i3 )
2362  9994 format( 1x, i6, ': ', a6, '(''', a1, ''',', i3, ',(', f4.1, ',',
2363  $ f4.1, '), X,', i2, ', Y,', i2, ', AP) ',
2364  $ ' .' )
2365  9993 format( 1x, i6, ': ', a6, '(''', a1, ''',', i3, ',(', f4.1, ',',
2366  $ f4.1, '), X,', i2, ', Y,', i2, ', A,', i3, ') ',
2367  $ ' .' )
2368  9992 format( ' ******* FATAL ERROR - ERROR-EXIT TAKEN ON VALID CALL *',
2369  $ '******' )
2370 *
2371 * End of CCHK6.
2372 *
2373  END
2374  SUBROUTINE cchke( ISNUM, SRNAMT, NOUT )
2375 *
2376 * Tests the error exits from the Level 2 Blas.
2377 * Requires a special version of the error-handling routine XERBLA.
2378 * ALPHA, RALPHA, BETA, A, X and Y should not need to be defined.
2379 *
2380 * Auxiliary routine for test program for Level 2 Blas.
2381 *
2382 * -- Written on 10-August-1987.
2383 * Richard Hanson, Sandia National Labs.
2384 * Jeremy Du Croz, NAG Central Office.
2385 *
2386 * .. Scalar Arguments ..
2387  INTEGER isnum, nout
2388  CHARACTER*6 srnamt
2389 * .. Scalars in Common ..
2390  INTEGER infot, noutc
2391  LOGICAL lerr, ok
2392 * .. Local Scalars ..
2393  COMPLEX alpha, beta
2394  REAL ralpha
2395 * .. Local Arrays ..
2396  COMPLEX a( 1, 1 ), x( 1 ), y( 1 )
2397 * .. External Subroutines ..
2398  EXTERNAL cgbmv, cgemv, cgerc, cgeru, chbmv, chemv, cher,
2399  $ cher2, chkxer, chpmv, chpr, chpr2, ctbmv,
2400  $ ctbsv, ctpmv, ctpsv, ctrmv, ctrsv
2401 * .. Common blocks ..
2402  common /infoc/infot, noutc, ok, lerr
2403 * .. Executable Statements ..
2404 * OK is set to .FALSE. by the special version of XERBLA or by CHKXER
2405 * if anything is wrong.
2406  ok = .true.
2407 * LERR is set to .TRUE. by the special version of XERBLA each time
2408 * it is called, and is then tested and re-set by CHKXER.
2409  lerr = .false.
2410  go to( 10, 20, 30, 40, 50, 60, 70, 80,
2411  $ 90, 100, 110, 120, 130, 140, 150, 160,
2412  $ 170 )isnum
2413  10 infot = 1
2414  CALL cgemv( '/', 0, 0, alpha, a, 1, x, 1, beta, y, 1 )
2415  CALL chkxer( srnamt, infot, nout, lerr, ok )
2416  infot = 2
2417  CALL cgemv( 'N', -1, 0, alpha, a, 1, x, 1, beta, y, 1 )
2418  CALL chkxer( srnamt, infot, nout, lerr, ok )
2419  infot = 3
2420  CALL cgemv( 'N', 0, -1, alpha, a, 1, x, 1, beta, y, 1 )
2421  CALL chkxer( srnamt, infot, nout, lerr, ok )
2422  infot = 6
2423  CALL cgemv( 'N', 2, 0, alpha, a, 1, x, 1, beta, y, 1 )
2424  CALL chkxer( srnamt, infot, nout, lerr, ok )
2425  infot = 8
2426  CALL cgemv( 'N', 0, 0, alpha, a, 1, x, 0, beta, y, 1 )
2427  CALL chkxer( srnamt, infot, nout, lerr, ok )
2428  infot = 11
2429  CALL cgemv( 'N', 0, 0, alpha, a, 1, x, 1, beta, y, 0 )
2430  CALL chkxer( srnamt, infot, nout, lerr, ok )
2431  go to 180
2432  20 infot = 1
2433  CALL cgbmv( '/', 0, 0, 0, 0, alpha, a, 1, x, 1, beta, y, 1 )
2434  CALL chkxer( srnamt, infot, nout, lerr, ok )
2435  infot = 2
2436  CALL cgbmv( 'N', -1, 0, 0, 0, alpha, a, 1, x, 1, beta, y, 1 )
2437  CALL chkxer( srnamt, infot, nout, lerr, ok )
2438  infot = 3
2439  CALL cgbmv( 'N', 0, -1, 0, 0, alpha, a, 1, x, 1, beta, y, 1 )
2440  CALL chkxer( srnamt, infot, nout, lerr, ok )
2441  infot = 4
2442  CALL cgbmv( 'N', 0, 0, -1, 0, alpha, a, 1, x, 1, beta, y, 1 )
2443  CALL chkxer( srnamt, infot, nout, lerr, ok )
2444  infot = 5
2445  CALL cgbmv( 'N', 2, 0, 0, -1, alpha, a, 1, x, 1, beta, y, 1 )
2446  CALL chkxer( srnamt, infot, nout, lerr, ok )
2447  infot = 8
2448  CALL cgbmv( 'N', 0, 0, 1, 0, alpha, a, 1, x, 1, beta, y, 1 )
2449  CALL chkxer( srnamt, infot, nout, lerr, ok )
2450  infot = 10
2451  CALL cgbmv( 'N', 0, 0, 0, 0, alpha, a, 1, x, 0, beta, y, 1 )
2452  CALL chkxer( srnamt, infot, nout, lerr, ok )
2453  infot = 13
2454  CALL cgbmv( 'N', 0, 0, 0, 0, alpha, a, 1, x, 1, beta, y, 0 )
2455  CALL chkxer( srnamt, infot, nout, lerr, ok )
2456  go to 180
2457  30 infot = 1
2458  CALL chemv( '/', 0, alpha, a, 1, x, 1, beta, y, 1 )
2459  CALL chkxer( srnamt, infot, nout, lerr, ok )
2460  infot = 2
2461  CALL chemv( 'U', -1, alpha, a, 1, x, 1, beta, y, 1 )
2462  CALL chkxer( srnamt, infot, nout, lerr, ok )
2463  infot = 5
2464  CALL chemv( 'U', 2, alpha, a, 1, x, 1, beta, y, 1 )
2465  CALL chkxer( srnamt, infot, nout, lerr, ok )
2466  infot = 7
2467  CALL chemv( 'U', 0, alpha, a, 1, x, 0, beta, y, 1 )
2468  CALL chkxer( srnamt, infot, nout, lerr, ok )
2469  infot = 10
2470  CALL chemv( 'U', 0, alpha, a, 1, x, 1, beta, y, 0 )
2471  CALL chkxer( srnamt, infot, nout, lerr, ok )
2472  go to 180
2473  40 infot = 1
2474  CALL chbmv( '/', 0, 0, alpha, a, 1, x, 1, beta, y, 1 )
2475  CALL chkxer( srnamt, infot, nout, lerr, ok )
2476  infot = 2
2477  CALL chbmv( 'U', -1, 0, alpha, a, 1, x, 1, beta, y, 1 )
2478  CALL chkxer( srnamt, infot, nout, lerr, ok )
2479  infot = 3
2480  CALL chbmv( 'U', 0, -1, alpha, a, 1, x, 1, beta, y, 1 )
2481  CALL chkxer( srnamt, infot, nout, lerr, ok )
2482  infot = 6
2483  CALL chbmv( 'U', 0, 1, alpha, a, 1, x, 1, beta, y, 1 )
2484  CALL chkxer( srnamt, infot, nout, lerr, ok )
2485  infot = 8
2486  CALL chbmv( 'U', 0, 0, alpha, a, 1, x, 0, beta, y, 1 )
2487  CALL chkxer( srnamt, infot, nout, lerr, ok )
2488  infot = 11
2489  CALL chbmv( 'U', 0, 0, alpha, a, 1, x, 1, beta, y, 0 )
2490  CALL chkxer( srnamt, infot, nout, lerr, ok )
2491  go to 180
2492  50 infot = 1
2493  CALL chpmv( '/', 0, alpha, a, x, 1, beta, y, 1 )
2494  CALL chkxer( srnamt, infot, nout, lerr, ok )
2495  infot = 2
2496  CALL chpmv( 'U', -1, alpha, a, x, 1, beta, y, 1 )
2497  CALL chkxer( srnamt, infot, nout, lerr, ok )
2498  infot = 6
2499  CALL chpmv( 'U', 0, alpha, a, x, 0, beta, y, 1 )
2500  CALL chkxer( srnamt, infot, nout, lerr, ok )
2501  infot = 9
2502  CALL chpmv( 'U', 0, alpha, a, x, 1, beta, y, 0 )
2503  CALL chkxer( srnamt, infot, nout, lerr, ok )
2504  go to 180
2505  60 infot = 1
2506  CALL ctrmv( '/', 'N', 'N', 0, a, 1, x, 1 )
2507  CALL chkxer( srnamt, infot, nout, lerr, ok )
2508  infot = 2
2509  CALL ctrmv( 'U', '/', 'N', 0, a, 1, x, 1 )
2510  CALL chkxer( srnamt, infot, nout, lerr, ok )
2511  infot = 3
2512  CALL ctrmv( 'U', 'N', '/', 0, a, 1, x, 1 )
2513  CALL chkxer( srnamt, infot, nout, lerr, ok )
2514  infot = 4
2515  CALL ctrmv( 'U', 'N', 'N', -1, a, 1, x, 1 )
2516  CALL chkxer( srnamt, infot, nout, lerr, ok )
2517  infot = 6
2518  CALL ctrmv( 'U', 'N', 'N', 2, a, 1, x, 1 )
2519  CALL chkxer( srnamt, infot, nout, lerr, ok )
2520  infot = 8
2521  CALL ctrmv( 'U', 'N', 'N', 0, a, 1, x, 0 )
2522  CALL chkxer( srnamt, infot, nout, lerr, ok )
2523  go to 180
2524  70 infot = 1
2525  CALL ctbmv( '/', 'N', 'N', 0, 0, a, 1, x, 1 )
2526  CALL chkxer( srnamt, infot, nout, lerr, ok )
2527  infot = 2
2528  CALL ctbmv( 'U', '/', 'N', 0, 0, a, 1, x, 1 )
2529  CALL chkxer( srnamt, infot, nout, lerr, ok )
2530  infot = 3
2531  CALL ctbmv( 'U', 'N', '/', 0, 0, a, 1, x, 1 )
2532  CALL chkxer( srnamt, infot, nout, lerr, ok )
2533  infot = 4
2534  CALL ctbmv( 'U', 'N', 'N', -1, 0, a, 1, x, 1 )
2535  CALL chkxer( srnamt, infot, nout, lerr, ok )
2536  infot = 5
2537  CALL ctbmv( 'U', 'N', 'N', 0, -1, a, 1, x, 1 )
2538  CALL chkxer( srnamt, infot, nout, lerr, ok )
2539  infot = 7
2540  CALL ctbmv( 'U', 'N', 'N', 0, 1, a, 1, x, 1 )
2541  CALL chkxer( srnamt, infot, nout, lerr, ok )
2542  infot = 9
2543  CALL ctbmv( 'U', 'N', 'N', 0, 0, a, 1, x, 0 )
2544  CALL chkxer( srnamt, infot, nout, lerr, ok )
2545  go to 180
2546  80 infot = 1
2547  CALL ctpmv( '/', 'N', 'N', 0, a, x, 1 )
2548  CALL chkxer( srnamt, infot, nout, lerr, ok )
2549  infot = 2
2550  CALL ctpmv( 'U', '/', 'N', 0, a, x, 1 )
2551  CALL chkxer( srnamt, infot, nout, lerr, ok )
2552  infot = 3
2553  CALL ctpmv( 'U', 'N', '/', 0, a, x, 1 )
2554  CALL chkxer( srnamt, infot, nout, lerr, ok )
2555  infot = 4
2556  CALL ctpmv( 'U', 'N', 'N', -1, a, x, 1 )
2557  CALL chkxer( srnamt, infot, nout, lerr, ok )
2558  infot = 7
2559  CALL ctpmv( 'U', 'N', 'N', 0, a, x, 0 )
2560  CALL chkxer( srnamt, infot, nout, lerr, ok )
2561  go to 180
2562  90 infot = 1
2563  CALL ctrsv( '/', 'N', 'N', 0, a, 1, x, 1 )
2564  CALL chkxer( srnamt, infot, nout, lerr, ok )
2565  infot = 2
2566  CALL ctrsv( 'U', '/', 'N', 0, a, 1, x, 1 )
2567  CALL chkxer( srnamt, infot, nout, lerr, ok )
2568  infot = 3
2569  CALL ctrsv( 'U', 'N', '/', 0, a, 1, x, 1 )
2570  CALL chkxer( srnamt, infot, nout, lerr, ok )
2571  infot = 4
2572  CALL ctrsv( 'U', 'N', 'N', -1, a, 1, x, 1 )
2573  CALL chkxer( srnamt, infot, nout, lerr, ok )
2574  infot = 6
2575  CALL ctrsv( 'U', 'N', 'N', 2, a, 1, x, 1 )
2576  CALL chkxer( srnamt, infot, nout, lerr, ok )
2577  infot = 8
2578  CALL ctrsv( 'U', 'N', 'N', 0, a, 1, x, 0 )
2579  CALL chkxer( srnamt, infot, nout, lerr, ok )
2580  go to 180
2581  100 infot = 1
2582  CALL ctbsv( '/', 'N', 'N', 0, 0, a, 1, x, 1 )
2583  CALL chkxer( srnamt, infot, nout, lerr, ok )
2584  infot = 2
2585  CALL ctbsv( 'U', '/', 'N', 0, 0, a, 1, x, 1 )
2586  CALL chkxer( srnamt, infot, nout, lerr, ok )
2587  infot = 3
2588  CALL ctbsv( 'U', 'N', '/', 0, 0, a, 1, x, 1 )
2589  CALL chkxer( srnamt, infot, nout, lerr, ok )
2590  infot = 4
2591  CALL ctbsv( 'U', 'N', 'N', -1, 0, a, 1, x, 1 )
2592  CALL chkxer( srnamt, infot, nout, lerr, ok )
2593  infot = 5
2594  CALL ctbsv( 'U', 'N', 'N', 0, -1, a, 1, x, 1 )
2595  CALL chkxer( srnamt, infot, nout, lerr, ok )
2596  infot = 7
2597  CALL ctbsv( 'U', 'N', 'N', 0, 1, a, 1, x, 1 )
2598  CALL chkxer( srnamt, infot, nout, lerr, ok )
2599  infot = 9
2600  CALL ctbsv( 'U', 'N', 'N', 0, 0, a, 1, x, 0 )
2601  CALL chkxer( srnamt, infot, nout, lerr, ok )
2602  go to 180
2603  110 infot = 1
2604  CALL ctpsv( '/', 'N', 'N', 0, a, x, 1 )
2605  CALL chkxer( srnamt, infot, nout, lerr, ok )
2606  infot = 2
2607  CALL ctpsv( 'U', '/', 'N', 0, a, x, 1 )
2608  CALL chkxer( srnamt, infot, nout, lerr, ok )
2609  infot = 3
2610  CALL ctpsv( 'U', 'N', '/', 0, a, x, 1 )
2611  CALL chkxer( srnamt, infot, nout, lerr, ok )
2612  infot = 4
2613  CALL ctpsv( 'U', 'N', 'N', -1, a, x, 1 )
2614  CALL chkxer( srnamt, infot, nout, lerr, ok )
2615  infot = 7
2616  CALL ctpsv( 'U', 'N', 'N', 0, a, x, 0 )
2617  CALL chkxer( srnamt, infot, nout, lerr, ok )
2618  go to 180
2619  120 infot = 1
2620  CALL cgerc( -1, 0, alpha, x, 1, y, 1, a, 1 )
2621  CALL chkxer( srnamt, infot, nout, lerr, ok )
2622  infot = 2
2623  CALL cgerc( 0, -1, alpha, x, 1, y, 1, a, 1 )
2624  CALL chkxer( srnamt, infot, nout, lerr, ok )
2625  infot = 5
2626  CALL cgerc( 0, 0, alpha, x, 0, y, 1, a, 1 )
2627  CALL chkxer( srnamt, infot, nout, lerr, ok )
2628  infot = 7
2629  CALL cgerc( 0, 0, alpha, x, 1, y, 0, a, 1 )
2630  CALL chkxer( srnamt, infot, nout, lerr, ok )
2631  infot = 9
2632  CALL cgerc( 2, 0, alpha, x, 1, y, 1, a, 1 )
2633  CALL chkxer( srnamt, infot, nout, lerr, ok )
2634  go to 180
2635  130 infot = 1
2636  CALL cgeru( -1, 0, alpha, x, 1, y, 1, a, 1 )
2637  CALL chkxer( srnamt, infot, nout, lerr, ok )
2638  infot = 2
2639  CALL cgeru( 0, -1, alpha, x, 1, y, 1, a, 1 )
2640  CALL chkxer( srnamt, infot, nout, lerr, ok )
2641  infot = 5
2642  CALL cgeru( 0, 0, alpha, x, 0, y, 1, a, 1 )
2643  CALL chkxer( srnamt, infot, nout, lerr, ok )
2644  infot = 7
2645  CALL cgeru( 0, 0, alpha, x, 1, y, 0, a, 1 )
2646  CALL chkxer( srnamt, infot, nout, lerr, ok )
2647  infot = 9
2648  CALL cgeru( 2, 0, alpha, x, 1, y, 1, a, 1 )
2649  CALL chkxer( srnamt, infot, nout, lerr, ok )
2650  go to 180
2651  140 infot = 1
2652  CALL cher( '/', 0, ralpha, x, 1, a, 1 )
2653  CALL chkxer( srnamt, infot, nout, lerr, ok )
2654  infot = 2
2655  CALL cher( 'U', -1, ralpha, x, 1, a, 1 )
2656  CALL chkxer( srnamt, infot, nout, lerr, ok )
2657  infot = 5
2658  CALL cher( 'U', 0, ralpha, x, 0, a, 1 )
2659  CALL chkxer( srnamt, infot, nout, lerr, ok )
2660  infot = 7
2661  CALL cher( 'U', 2, ralpha, x, 1, a, 1 )
2662  CALL chkxer( srnamt, infot, nout, lerr, ok )
2663  go to 180
2664  150 infot = 1
2665  CALL chpr( '/', 0, ralpha, x, 1, a )
2666  CALL chkxer( srnamt, infot, nout, lerr, ok )
2667  infot = 2
2668  CALL chpr( 'U', -1, ralpha, x, 1, a )
2669  CALL chkxer( srnamt, infot, nout, lerr, ok )
2670  infot = 5
2671  CALL chpr( 'U', 0, ralpha, x, 0, a )
2672  CALL chkxer( srnamt, infot, nout, lerr, ok )
2673  go to 180
2674  160 infot = 1
2675  CALL cher2( '/', 0, alpha, x, 1, y, 1, a, 1 )
2676  CALL chkxer( srnamt, infot, nout, lerr, ok )
2677  infot = 2
2678  CALL cher2( 'U', -1, alpha, x, 1, y, 1, a, 1 )
2679  CALL chkxer( srnamt, infot, nout, lerr, ok )
2680  infot = 5
2681  CALL cher2( 'U', 0, alpha, x, 0, y, 1, a, 1 )
2682  CALL chkxer( srnamt, infot, nout, lerr, ok )
2683  infot = 7
2684  CALL cher2( 'U', 0, alpha, x, 1, y, 0, a, 1 )
2685  CALL chkxer( srnamt, infot, nout, lerr, ok )
2686  infot = 9
2687  CALL cher2( 'U', 2, alpha, x, 1, y, 1, a, 1 )
2688  CALL chkxer( srnamt, infot, nout, lerr, ok )
2689  go to 180
2690  170 infot = 1
2691  CALL chpr2( '/', 0, alpha, x, 1, y, 1, a )
2692  CALL chkxer( srnamt, infot, nout, lerr, ok )
2693  infot = 2
2694  CALL chpr2( 'U', -1, alpha, x, 1, y, 1, a )
2695  CALL chkxer( srnamt, infot, nout, lerr, ok )
2696  infot = 5
2697  CALL chpr2( 'U', 0, alpha, x, 0, y, 1, a )
2698  CALL chkxer( srnamt, infot, nout, lerr, ok )
2699  infot = 7
2700  CALL chpr2( 'U', 0, alpha, x, 1, y, 0, a )
2701  CALL chkxer( srnamt, infot, nout, lerr, ok )
2702 *
2703  180 IF( ok )THEN
2704  WRITE( nout, fmt = 9999 )srnamt
2705  ELSE
2706  WRITE( nout, fmt = 9998 )srnamt
2707  END IF
2708  return
2709 *
2710  9999 format( ' ', a6, ' PASSED THE TESTS OF ERROR-EXITS' )
2711  9998 format( ' ******* ', a6, ' FAILED THE TESTS OF ERROR-EXITS *****',
2712  $ '**' )
2713 *
2714 * End of CCHKE.
2715 *
2716  END
2717  SUBROUTINE cmake( TYPE, UPLO, DIAG, M, N, A, NMAX, AA, LDA, KL,
2718  $ ku, reset, transl )
2719 *
2720 * Generates values for an M by N matrix A within the bandwidth
2721 * defined by KL and KU.
2722 * Stores the values in the array AA in the data structure required
2723 * by the routine, with unwanted elements set to rogue value.
2724 *
2725 * TYPE is 'GE', 'GB', 'HE', 'HB', 'HP', 'TR', 'TB' OR 'TP'.
2726 *
2727 * Auxiliary routine for test program for Level 2 Blas.
2728 *
2729 * -- Written on 10-August-1987.
2730 * Richard Hanson, Sandia National Labs.
2731 * Jeremy Du Croz, NAG Central Office.
2732 *
2733 * .. Parameters ..
2734  COMPLEX zero, one
2735  parameter( zero = ( 0.0, 0.0 ), one = ( 1.0, 0.0 ) )
2736  COMPLEX rogue
2737  parameter( rogue = ( -1.0e10, 1.0e10 ) )
2738  REAL rzero
2739  parameter( rzero = 0.0 )
2740  REAL rrogue
2741  parameter( rrogue = -1.0e10 )
2742 * .. Scalar Arguments ..
2743  COMPLEX transl
2744  INTEGER kl, ku, lda, m, n, nmax
2745  LOGICAL reset
2746  CHARACTER*1 diag, uplo
2747  CHARACTER*2 type
2748 * .. Array Arguments ..
2749  COMPLEX a( nmax, * ), aa( * )
2750 * .. Local Scalars ..
2751  INTEGER i, i1, i2, i3, ibeg, iend, ioff, j, jj, kk
2752  LOGICAL gen, lower, sym, tri, unit, upper
2753 * .. External Functions ..
2754  COMPLEX cbeg
2755  EXTERNAL cbeg
2756 * .. Intrinsic Functions ..
2757  INTRINSIC cmplx, conjg, max, min, real
2758 * .. Executable Statements ..
2759  gen = TYPE( 1: 1 ).EQ.'G'
2760  sym = TYPE( 1: 1 ).EQ.'H'
2761  tri = TYPE( 1: 1 ).EQ.'T'
2762  upper = ( sym.OR.tri ).AND.uplo.EQ.'U'
2763  lower = ( sym.OR.tri ).AND.uplo.EQ.'L'
2764  unit = tri.AND.diag.EQ.'U'
2765 *
2766 * Generate data in array A.
2767 *
2768  DO 20 j = 1, n
2769  DO 10 i = 1, m
2770  IF( gen.OR.( upper.AND.i.LE.j ).OR.( lower.AND.i.GE.j ) )
2771  $ THEN
2772  IF( ( i.LE.j.AND.j - i.LE.ku ).OR.
2773  $ ( i.GE.j.AND.i - j.LE.kl ) )THEN
2774  a( i, j ) = cbeg( reset ) + transl
2775  ELSE
2776  a( i, j ) = zero
2777  END IF
2778  IF( i.NE.j )THEN
2779  IF( sym )THEN
2780  a( j, i ) = conjg( a( i, j ) )
2781  ELSE IF( tri )THEN
2782  a( j, i ) = zero
2783  END IF
2784  END IF
2785  END IF
2786  10 continue
2787  IF( sym )
2788  $ a( j, j ) = cmplx( REAL( A( J, J ) ), rzero )
2789  IF( tri )
2790  $ a( j, j ) = a( j, j ) + one
2791  IF( unit )
2792  $ a( j, j ) = one
2793  20 continue
2794 *
2795 * Store elements in array AS in data structure required by routine.
2796 *
2797  IF( type.EQ.'GE' )THEN
2798  DO 50 j = 1, n
2799  DO 30 i = 1, m
2800  aa( i + ( j - 1 )*lda ) = a( i, j )
2801  30 continue
2802  DO 40 i = m + 1, lda
2803  aa( i + ( j - 1 )*lda ) = rogue
2804  40 continue
2805  50 continue
2806  ELSE IF( type.EQ.'GB' )THEN
2807  DO 90 j = 1, n
2808  DO 60 i1 = 1, ku + 1 - j
2809  aa( i1 + ( j - 1 )*lda ) = rogue
2810  60 continue
2811  DO 70 i2 = i1, min( kl + ku + 1, ku + 1 + m - j )
2812  aa( i2 + ( j - 1 )*lda ) = a( i2 + j - ku - 1, j )
2813  70 continue
2814  DO 80 i3 = i2, lda
2815  aa( i3 + ( j - 1 )*lda ) = rogue
2816  80 continue
2817  90 continue
2818  ELSE IF( type.EQ.'HE'.OR.type.EQ.'TR' )THEN
2819  DO 130 j = 1, n
2820  IF( upper )THEN
2821  ibeg = 1
2822  IF( unit )THEN
2823  iend = j - 1
2824  ELSE
2825  iend = j
2826  END IF
2827  ELSE
2828  IF( unit )THEN
2829  ibeg = j + 1
2830  ELSE
2831  ibeg = j
2832  END IF
2833  iend = n
2834  END IF
2835  DO 100 i = 1, ibeg - 1
2836  aa( i + ( j - 1 )*lda ) = rogue
2837  100 continue
2838  DO 110 i = ibeg, iend
2839  aa( i + ( j - 1 )*lda ) = a( i, j )
2840  110 continue
2841  DO 120 i = iend + 1, lda
2842  aa( i + ( j - 1 )*lda ) = rogue
2843  120 continue
2844  IF( sym )THEN
2845  jj = j + ( j - 1 )*lda
2846  aa( jj ) = cmplx( REAL( AA( JJ ) ), rrogue )
2847  END IF
2848  130 continue
2849  ELSE IF( type.EQ.'HB'.OR.type.EQ.'TB' )THEN
2850  DO 170 j = 1, n
2851  IF( upper )THEN
2852  kk = kl + 1
2853  ibeg = max( 1, kl + 2 - j )
2854  IF( unit )THEN
2855  iend = kl
2856  ELSE
2857  iend = kl + 1
2858  END IF
2859  ELSE
2860  kk = 1
2861  IF( unit )THEN
2862  ibeg = 2
2863  ELSE
2864  ibeg = 1
2865  END IF
2866  iend = min( kl + 1, 1 + m - j )
2867  END IF
2868  DO 140 i = 1, ibeg - 1
2869  aa( i + ( j - 1 )*lda ) = rogue
2870  140 continue
2871  DO 150 i = ibeg, iend
2872  aa( i + ( j - 1 )*lda ) = a( i + j - kk, j )
2873  150 continue
2874  DO 160 i = iend + 1, lda
2875  aa( i + ( j - 1 )*lda ) = rogue
2876  160 continue
2877  IF( sym )THEN
2878  jj = kk + ( j - 1 )*lda
2879  aa( jj ) = cmplx( REAL( AA( JJ ) ), rrogue )
2880  END IF
2881  170 continue
2882  ELSE IF( type.EQ.'HP'.OR.type.EQ.'TP' )THEN
2883  ioff = 0
2884  DO 190 j = 1, n
2885  IF( upper )THEN
2886  ibeg = 1
2887  iend = j
2888  ELSE
2889  ibeg = j
2890  iend = n
2891  END IF
2892  DO 180 i = ibeg, iend
2893  ioff = ioff + 1
2894  aa( ioff ) = a( i, j )
2895  IF( i.EQ.j )THEN
2896  IF( unit )
2897  $ aa( ioff ) = rogue
2898  IF( sym )
2899  $ aa( ioff ) = cmplx( REAL( AA( IOFF ) ), rrogue )
2900  END IF
2901  180 continue
2902  190 continue
2903  END IF
2904  return
2905 *
2906 * End of CMAKE.
2907 *
2908  END
2909  SUBROUTINE cmvch( TRANS, M, N, ALPHA, A, NMAX, X, INCX, BETA, Y,
2910  $ incy, yt, g, yy, eps, err, fatal, nout, mv )
2911 *
2912 * Checks the results of the computational tests.
2913 *
2914 * Auxiliary routine for test program for Level 2 Blas.
2915 *
2916 * -- Written on 10-August-1987.
2917 * Richard Hanson, Sandia National Labs.
2918 * Jeremy Du Croz, NAG Central Office.
2919 *
2920 * .. Parameters ..
2921  COMPLEX zero
2922  parameter( zero = ( 0.0, 0.0 ) )
2923  REAL rzero, rone
2924  parameter( rzero = 0.0, rone = 1.0 )
2925 * .. Scalar Arguments ..
2926  COMPLEX alpha, beta
2927  REAL eps, err
2928  INTEGER incx, incy, m, n, nmax, nout
2929  LOGICAL fatal, mv
2930  CHARACTER*1 trans
2931 * .. Array Arguments ..
2932  COMPLEX a( nmax, * ), x( * ), y( * ), yt( * ), yy( * )
2933  REAL g( * )
2934 * .. Local Scalars ..
2935  COMPLEX c
2936  REAL erri
2937  INTEGER i, incxl, incyl, iy, j, jx, kx, ky, ml, nl
2938  LOGICAL ctran, tran
2939 * .. Intrinsic Functions ..
2940  INTRINSIC abs, aimag, conjg, max, REAL, sqrt
2941 * .. Statement Functions ..
2942  REAL abs1
2943 * .. Statement Function definitions ..
2944  abs1( c ) = abs( REAL( C ) ) + abs( aimag( c ) )
2945 * .. Executable Statements ..
2946  tran = trans.EQ.'T'
2947  ctran = trans.EQ.'C'
2948  IF( tran.OR.ctran )THEN
2949  ml = n
2950  nl = m
2951  ELSE
2952  ml = m
2953  nl = n
2954  END IF
2955  IF( incx.LT.0 )THEN
2956  kx = nl
2957  incxl = -1
2958  ELSE
2959  kx = 1
2960  incxl = 1
2961  END IF
2962  IF( incy.LT.0 )THEN
2963  ky = ml
2964  incyl = -1
2965  ELSE
2966  ky = 1
2967  incyl = 1
2968  END IF
2969 *
2970 * Compute expected result in YT using data in A, X and Y.
2971 * Compute gauges in G.
2972 *
2973  iy = ky
2974  DO 40 i = 1, ml
2975  yt( iy ) = zero
2976  g( iy ) = rzero
2977  jx = kx
2978  IF( tran )THEN
2979  DO 10 j = 1, nl
2980  yt( iy ) = yt( iy ) + a( j, i )*x( jx )
2981  g( iy ) = g( iy ) + abs1( a( j, i ) )*abs1( x( jx ) )
2982  jx = jx + incxl
2983  10 continue
2984  ELSE IF( ctran )THEN
2985  DO 20 j = 1, nl
2986  yt( iy ) = yt( iy ) + conjg( a( j, i ) )*x( jx )
2987  g( iy ) = g( iy ) + abs1( a( j, i ) )*abs1( x( jx ) )
2988  jx = jx + incxl
2989  20 continue
2990  ELSE
2991  DO 30 j = 1, nl
2992  yt( iy ) = yt( iy ) + a( i, j )*x( jx )
2993  g( iy ) = g( iy ) + abs1( a( i, j ) )*abs1( x( jx ) )
2994  jx = jx + incxl
2995  30 continue
2996  END IF
2997  yt( iy ) = alpha*yt( iy ) + beta*y( iy )
2998  g( iy ) = abs1( alpha )*g( iy ) + abs1( beta )*abs1( y( iy ) )
2999  iy = iy + incyl
3000  40 continue
3001 *
3002 * Compute the error ratio for this result.
3003 *
3004  err = zero
3005  DO 50 i = 1, ml
3006  erri = abs( yt( i ) - yy( 1 + ( i - 1 )*abs( incy ) ) )/eps
3007  IF( g( i ).NE.rzero )
3008  $ erri = erri/g( i )
3009  err = max( err, erri )
3010  IF( err*sqrt( eps ).GE.rone )
3011  $ go to 60
3012  50 continue
3013 * If the loop completes, all results are at least half accurate.
3014  go to 80
3015 *
3016 * Report fatal error.
3017 *
3018  60 fatal = .true.
3019  WRITE( nout, fmt = 9999 )
3020  DO 70 i = 1, ml
3021  IF( mv )THEN
3022  WRITE( nout, fmt = 9998 )i, yt( i ),
3023  $ yy( 1 + ( i - 1 )*abs( incy ) )
3024  ELSE
3025  WRITE( nout, fmt = 9998 )i,
3026  $ yy( 1 + ( i - 1 )*abs( incy ) ), yt( i )
3027  END IF
3028  70 continue
3029 *
3030  80 continue
3031  return
3032 *
3033  9999 format( ' ******* FATAL ERROR - COMPUTED RESULT IS LESS THAN HAL',
3034  $ 'F ACCURATE *******', /' EXPECTED RE',
3035  $ 'SULT COMPUTED RESULT' )
3036  9998 format( 1x, i7, 2( ' (', g15.6, ',', g15.6, ')' ) )
3037 *
3038 * End of CMVCH.
3039 *
3040  END
3041  LOGICAL FUNCTION lce( RI, RJ, LR )
3042 *
3043 * Tests if two arrays are identical.
3044 *
3045 * Auxiliary routine for test program for Level 2 Blas.
3046 *
3047 * -- Written on 10-August-1987.
3048 * Richard Hanson, Sandia National Labs.
3049 * Jeremy Du Croz, NAG Central Office.
3050 *
3051 * .. Scalar Arguments ..
3052  INTEGER lr
3053 * .. Array Arguments ..
3054  COMPLEX ri( * ), rj( * )
3055 * .. Local Scalars ..
3056  INTEGER i
3057 * .. Executable Statements ..
3058  DO 10 i = 1, lr
3059  IF( ri( i ).NE.rj( i ) )
3060  $ go to 20
3061  10 continue
3062  lce = .true.
3063  go to 30
3064  20 continue
3065  lce = .false.
3066  30 return
3067 *
3068 * End of LCE.
3069 *
3070  END
3071  LOGICAL FUNCTION lceres( TYPE, UPLO, M, N, AA, AS, LDA )
3072 *
3073 * Tests if selected elements in two arrays are equal.
3074 *
3075 * TYPE is 'GE', 'HE' or 'HP'.
3076 *
3077 * Auxiliary routine for test program for Level 2 Blas.
3078 *
3079 * -- Written on 10-August-1987.
3080 * Richard Hanson, Sandia National Labs.
3081 * Jeremy Du Croz, NAG Central Office.
3082 *
3083 * .. Scalar Arguments ..
3084  INTEGER lda, m, n
3085  CHARACTER*1 uplo
3086  CHARACTER*2 type
3087 * .. Array Arguments ..
3088  COMPLEX aa( lda, * ), as( lda, * )
3089 * .. Local Scalars ..
3090  INTEGER i, ibeg, iend, j
3091  LOGICAL upper
3092 * .. Executable Statements ..
3093  upper = uplo.EQ.'U'
3094  IF( type.EQ.'GE' )THEN
3095  DO 20 j = 1, n
3096  DO 10 i = m + 1, lda
3097  IF( aa( i, j ).NE.as( i, j ) )
3098  $ go to 70
3099  10 continue
3100  20 continue
3101  ELSE IF( type.EQ.'HE' )THEN
3102  DO 50 j = 1, n
3103  IF( upper )THEN
3104  ibeg = 1
3105  iend = j
3106  ELSE
3107  ibeg = j
3108  iend = n
3109  END IF
3110  DO 30 i = 1, ibeg - 1
3111  IF( aa( i, j ).NE.as( i, j ) )
3112  $ go to 70
3113  30 continue
3114  DO 40 i = iend + 1, lda
3115  IF( aa( i, j ).NE.as( i, j ) )
3116  $ go to 70
3117  40 continue
3118  50 continue
3119  END IF
3120 *
3121  lceres = .true.
3122  go to 80
3123  70 continue
3124  lceres = .false.
3125  80 return
3126 *
3127 * End of LCERES.
3128 *
3129  END
3130  COMPLEX FUNCTION cbeg( RESET )
3131 *
3132 * Generates complex numbers as pairs of random numbers uniformly
3133 * distributed between -0.5 and 0.5.
3134 *
3135 * Auxiliary routine for test program for Level 2 Blas.
3136 *
3137 * -- Written on 10-August-1987.
3138 * Richard Hanson, Sandia National Labs.
3139 * Jeremy Du Croz, NAG Central Office.
3140 *
3141 * .. Scalar Arguments ..
3142  LOGICAL reset
3143 * .. Local Scalars ..
3144  INTEGER i, ic, j, mi, mj
3145 * .. Save statement ..
3146  SAVE i, ic, j, mi, mj
3147 * .. Intrinsic Functions ..
3148  INTRINSIC cmplx
3149 * .. Executable Statements ..
3150  IF( reset )THEN
3151 * Initialize local variables.
3152  mi = 891
3153  mj = 457
3154  i = 7
3155  j = 7
3156  ic = 0
3157  reset = .false.
3158  END IF
3159 *
3160 * The sequence of values of I or J is bounded between 1 and 999.
3161 * If initial I or J = 1,2,3,6,7 or 9, the period will be 50.
3162 * If initial I or J = 4 or 8, the period will be 25.
3163 * If initial I or J = 5, the period will be 10.
3164 * IC is used to break up the period by skipping 1 value of I or J
3165 * in 6.
3166 *
3167  ic = ic + 1
3168  10 i = i*mi
3169  j = j*mj
3170  i = i - 1000*( i/1000 )
3171  j = j - 1000*( j/1000 )
3172  IF( ic.GE.5 )THEN
3173  ic = 0
3174  go to 10
3175  END IF
3176  cbeg = cmplx( ( i - 500 )/1001.0, ( j - 500 )/1001.0 )
3177  return
3178 *
3179 * End of CBEG.
3180 *
3181  END
3182  REAL FUNCTION sdiff( X, Y )
3183 *
3184 * Auxiliary routine for test program for Level 2 Blas.
3185 *
3186 * -- Written on 10-August-1987.
3187 * Richard Hanson, Sandia National Labs.
3188 *
3189 * .. Scalar Arguments ..
3190  REAL x, y
3191 * .. Executable Statements ..
3192  sdiff = x - y
3193  return
3194 *
3195 * End of SDIFF.
3196 *
3197  END
3198  SUBROUTINE chkxer( SRNAMT, INFOT, NOUT, LERR, OK )
3199 *
3200 * Tests whether XERBLA has detected an error when it should.
3201 *
3202 * Auxiliary routine for test program for Level 2 Blas.
3203 *
3204 * -- Written on 10-August-1987.
3205 * Richard Hanson, Sandia National Labs.
3206 * Jeremy Du Croz, NAG Central Office.
3207 *
3208 * .. Scalar Arguments ..
3209  INTEGER infot, nout
3210  LOGICAL lerr, ok
3211  CHARACTER*6 srnamt
3212 * .. Executable Statements ..
3213  IF( .NOT.lerr )THEN
3214  WRITE( nout, fmt = 9999 )infot, srnamt
3215  ok = .false.
3216  END IF
3217  lerr = .false.
3218  return
3219 *
3220  9999 format( ' ***** ILLEGAL VALUE OF PARAMETER NUMBER ', i2, ' NOT D',
3221  $ 'ETECTED BY ', a6, ' *****' )
3222 *
3223 * End of CHKXER.
3224 *
3225  END
3226  SUBROUTINE xerbla( SRNAME, INFO )
3227 *
3228 * This is a special version of XERBLA to be used only as part of
3229 * the test program for testing error exits from the Level 2 BLAS
3230 * routines.
3231 *
3232 * XERBLA is an error handler for the Level 2 BLAS routines.
3233 *
3234 * It is called by the Level 2 BLAS routines if an input parameter is
3235 * invalid.
3236 *
3237 * Auxiliary routine for test program for Level 2 Blas.
3238 *
3239 * -- Written on 10-August-1987.
3240 * Richard Hanson, Sandia National Labs.
3241 * Jeremy Du Croz, NAG Central Office.
3242 *
3243 * .. Scalar Arguments ..
3244  INTEGER info
3245  CHARACTER*6 srname
3246 * .. Scalars in Common ..
3247  INTEGER infot, nout
3248  LOGICAL lerr, ok
3249  CHARACTER*6 srnamt
3250 * .. Common blocks ..
3251  common /infoc/infot, nout, ok, lerr
3252  common /srnamc/srnamt
3253 * .. Executable Statements ..
3254  lerr = .true.
3255  IF( info.NE.infot )THEN
3256  IF( infot.NE.0 )THEN
3257  WRITE( nout, fmt = 9999 )info, infot
3258  ELSE
3259  WRITE( nout, fmt = 9997 )info
3260  END IF
3261  ok = .false.
3262  END IF
3263  IF( srname.NE.srnamt )THEN
3264  WRITE( nout, fmt = 9998 )srname, srnamt
3265  ok = .false.
3266  END IF
3267  return
3268 *
3269  9999 format( ' ******* XERBLA WAS CALLED WITH INFO = ', i6, ' INSTEAD',
3270  $ ' OF ', i2, ' *******' )
3271  9998 format( ' ******* XERBLA WAS CALLED WITH SRNAME = ', a6, ' INSTE',
3272  $ 'AD OF ', a6, ' *******' )
3273  9997 format( ' ******* XERBLA WAS CALLED WITH INFO = ', i6,
3274  $ ' *******' )
3275 *
3276 * End of XERBLA
3277 *
3278  END
3279