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