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