LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
cchkaa.f
Go to the documentation of this file.
1 *> \brief \b CCHKAA
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 CCHKAA
12 *
13 *
14 *> \par Purpose:
15 * =============
16 *>
17 *> \verbatim
18 *>
19 *> CCHKAA is the main test program for the COMPLEX linear equation
20 *> routines.
21 *>
22 *> The program must be driven by a short data file. The first 15 records
23 *> (not including the first comment line) specify problem dimensions
24 *> and program options using list-directed input. The remaining lines
25 *> specify the LAPACK test paths and the number of matrix types to use
26 *> in testing. An annotated example of a data file can be obtained by
27 *> deleting the first 3 characters from the following 42 lines:
28 *> Data file for testing COMPLEX LAPACK linear equation routines
29 *> 7 Number of values of M
30 *> 0 1 2 3 5 10 16 Values of M (row dimension)
31 *> 7 Number of values of N
32 *> 0 1 2 3 5 10 16 Values of N (column dimension)
33 *> 1 Number of values of NRHS
34 *> 2 Values of NRHS (number of right hand sides)
35 *> 5 Number of values of NB
36 *> 1 3 3 3 20 Values of NB (the blocksize)
37 *> 1 0 5 9 1 Values of NX (crossover point)
38 *> 3 Number of values of RANK
39 *> 30 50 90 Values of rank (as a % of N)
40 *> 30.0 Threshold value of test ratio
41 *> T Put T to test the LAPACK routines
42 *> T Put T to test the driver routines
43 *> T Put T to test the error exits
44 *> CGE 11 List types on next line if 0 < NTYPES < 11
45 *> CGB 8 List types on next line if 0 < NTYPES < 8
46 *> CGT 12 List types on next line if 0 < NTYPES < 12
47 *> CPO 9 List types on next line if 0 < NTYPES < 9
48 *> CPO 9 List types on next line if 0 < NTYPES < 9
49 *> CPP 9 List types on next line if 0 < NTYPES < 9
50 *> CPB 8 List types on next line if 0 < NTYPES < 8
51 *> CPT 12 List types on next line if 0 < NTYPES < 12
52 *> CHE 10 List types on next line if 0 < NTYPES < 10
53 *> CHR 10 List types on next line if 0 < NTYPES < 10
54 *> CHK 10 List types on next line if 0 < NTYPES < 10
55 *> CHA 10 List types on next line if 0 < NTYPES < 10
56 *> CHP 10 List types on next line if 0 < NTYPES < 10
57 *> CSY 11 List types on next line if 0 < NTYPES < 11
58 *> CSK 11 List types on next line if 0 < NTYPES < 11
59 *> CSR 11 List types on next line if 0 < NTYPES < 11
60 *> CSP 11 List types on next line if 0 < NTYPES < 11
61 *> CTR 18 List types on next line if 0 < NTYPES < 18
62 *> CTP 18 List types on next line if 0 < NTYPES < 18
63 *> CTB 17 List types on next line if 0 < NTYPES < 17
64 *> CQR 8 List types on next line if 0 < NTYPES < 8
65 *> CRQ 8 List types on next line if 0 < NTYPES < 8
66 *> CLQ 8 List types on next line if 0 < NTYPES < 8
67 *> CQL 8 List types on next line if 0 < NTYPES < 8
68 *> CQP 6 List types on next line if 0 < NTYPES < 6
69 *> CTZ 3 List types on next line if 0 < NTYPES < 3
70 *> CLS 6 List types on next line if 0 < NTYPES < 6
71 *> CEQ
72 *> CQT
73 *> CQX
74 *> \endverbatim
75 *
76 * Parameters:
77 * ==========
78 *
79 *> \verbatim
80 *> NMAX INTEGER
81 *> The maximum allowable value for M and N.
82 *>
83 *> MAXIN INTEGER
84 *> The number of different values that can be used for each of
85 *> M, N, NRHS, NB, NX and RANK
86 *>
87 *> MAXRHS INTEGER
88 *> The maximum number of right hand sides
89 *>
90 *> MATMAX INTEGER
91 *> The maximum number of matrix types to use for testing
92 *>
93 *> NIN INTEGER
94 *> The unit number for input
95 *>
96 *> NOUT INTEGER
97 *> The unit number for output
98 *> \endverbatim
99 *
100 * Authors:
101 * ========
102 *
103 *> \author Univ. of Tennessee
104 *> \author Univ. of California Berkeley
105 *> \author Univ. of Colorado Denver
106 *> \author NAG Ltd.
107 *
108 *> \date December 2016
109 *
110 *> \ingroup complex_lin
111 *
112 * =====================================================================
113  PROGRAM cchkaa
114 *
115 * -- LAPACK test routine (version 3.7.0) --
116 * -- LAPACK is a software package provided by Univ. of Tennessee, --
117 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
118 * December 2016
119 *
120 * =====================================================================
121 *
122 * .. Parameters ..
123  INTEGER NMAX
124  parameter ( nmax = 132 )
125  INTEGER MAXIN
126  parameter ( maxin = 12 )
127  INTEGER MAXRHS
128  parameter ( maxrhs = 16 )
129  INTEGER MATMAX
130  parameter ( matmax = 30 )
131  INTEGER NIN, NOUT
132  parameter ( nin = 5, nout = 6 )
133  INTEGER KDMAX
134  parameter ( kdmax = nmax+( nmax+1 ) / 4 )
135 * ..
136 * .. Local Scalars ..
137  LOGICAL FATAL, TSTCHK, TSTDRV, TSTERR
138  CHARACTER C1
139  CHARACTER*2 C2
140  CHARACTER*3 PATH
141  CHARACTER*10 INTSTR
142  CHARACTER*72 ALINE
143  INTEGER I, IC, J, K, LA, LAFAC, LDA, NB, NM, NMATS, NN,
144  $ nnb, nnb2, nns, nrhs, ntypes, nrank,
145  $ vers_major, vers_minor, vers_patch
146  REAL EPS, S1, S2, THREQ, THRESH
147 * ..
148 * .. Local Arrays ..
149  LOGICAL DOTYPE( matmax )
150  INTEGER IWORK( 25*nmax ), MVAL( maxin ),
151  $ nbval( maxin ), nbval2( maxin ),
152  $ nsval( maxin ), nval( maxin ), nxval( maxin ),
153  $ rankval( maxin ), piv( nmax )
154  REAL RWORK( 150*nmax+2*maxrhs ), S( 2*nmax )
155  COMPLEX A( ( kdmax+1 )*nmax, 7 ), B( nmax*maxrhs, 4 ),
156  $ e( nmax ), work( nmax, nmax+maxrhs+10 )
157 * ..
158 * .. External Functions ..
159  LOGICAL LSAME, LSAMEN
160  REAL SECOND, SLAMCH
161  EXTERNAL lsame, lsamen, second, slamch
162 * ..
163 * .. External Subroutines ..
164  EXTERNAL alareq, cchkeq, cchkgb, cchkge, cchkgt, cchkhe,
174 * ..
175 * .. Scalars in Common ..
176  LOGICAL LERR, OK
177  CHARACTER*32 SRNAMT
178  INTEGER INFOT, NUNIT
179 * ..
180 * .. Arrays in Common ..
181  INTEGER IPARMS( 100 )
182 * ..
183 * .. Common blocks ..
184  COMMON / claenv / iparms
185  COMMON / infoc / infot, nunit, ok, lerr
186  COMMON / srnamc / srnamt
187 * ..
188 * .. Data statements ..
189  DATA threq / 2.0 / , intstr / '0123456789' /
190 * ..
191 * .. Executable Statements ..
192 *
193  s1 = second( )
194  lda = nmax
195  fatal = .false.
196 *
197 * Read a dummy line.
198 *
199  READ( nin, fmt = * )
200 *
201 * Report values of parameters.
202 *
203  CALL ilaver( vers_major, vers_minor, vers_patch )
204  WRITE( nout, fmt = 9994 ) vers_major, vers_minor, vers_patch
205 *
206 * Read the values of M
207 *
208  READ( nin, fmt = * )nm
209  IF( nm.LT.1 ) THEN
210  WRITE( nout, fmt = 9996 )' NM ', nm, 1
211  nm = 0
212  fatal = .true.
213  ELSE IF( nm.GT.maxin ) THEN
214  WRITE( nout, fmt = 9995 )' NM ', nm, maxin
215  nm = 0
216  fatal = .true.
217  END IF
218  READ( nin, fmt = * )( mval( i ), i = 1, nm )
219  DO 10 i = 1, nm
220  IF( mval( i ).LT.0 ) THEN
221  WRITE( nout, fmt = 9996 )' M ', mval( i ), 0
222  fatal = .true.
223  ELSE IF( mval( i ).GT.nmax ) THEN
224  WRITE( nout, fmt = 9995 )' M ', mval( i ), nmax
225  fatal = .true.
226  END IF
227  10 CONTINUE
228  IF( nm.GT.0 )
229  $ WRITE( nout, fmt = 9993 )'M ', ( mval( i ), i = 1, nm )
230 *
231 * Read the values of N
232 *
233  READ( nin, fmt = * )nn
234  IF( nn.LT.1 ) THEN
235  WRITE( nout, fmt = 9996 )' NN ', nn, 1
236  nn = 0
237  fatal = .true.
238  ELSE IF( nn.GT.maxin ) THEN
239  WRITE( nout, fmt = 9995 )' NN ', nn, maxin
240  nn = 0
241  fatal = .true.
242  END IF
243  READ( nin, fmt = * )( nval( i ), i = 1, nn )
244  DO 20 i = 1, nn
245  IF( nval( i ).LT.0 ) THEN
246  WRITE( nout, fmt = 9996 )' N ', nval( i ), 0
247  fatal = .true.
248  ELSE IF( nval( i ).GT.nmax ) THEN
249  WRITE( nout, fmt = 9995 )' N ', nval( i ), nmax
250  fatal = .true.
251  END IF
252  20 CONTINUE
253  IF( nn.GT.0 )
254  $ WRITE( nout, fmt = 9993 )'N ', ( nval( i ), i = 1, nn )
255 *
256 * Read the values of NRHS
257 *
258  READ( nin, fmt = * )nns
259  IF( nns.LT.1 ) THEN
260  WRITE( nout, fmt = 9996 )' NNS', nns, 1
261  nns = 0
262  fatal = .true.
263  ELSE IF( nns.GT.maxin ) THEN
264  WRITE( nout, fmt = 9995 )' NNS', nns, maxin
265  nns = 0
266  fatal = .true.
267  END IF
268  READ( nin, fmt = * )( nsval( i ), i = 1, nns )
269  DO 30 i = 1, nns
270  IF( nsval( i ).LT.0 ) THEN
271  WRITE( nout, fmt = 9996 )'NRHS', nsval( i ), 0
272  fatal = .true.
273  ELSE IF( nsval( i ).GT.maxrhs ) THEN
274  WRITE( nout, fmt = 9995 )'NRHS', nsval( i ), maxrhs
275  fatal = .true.
276  END IF
277  30 CONTINUE
278  IF( nns.GT.0 )
279  $ WRITE( nout, fmt = 9993 )'NRHS', ( nsval( i ), i = 1, nns )
280 *
281 * Read the values of NB
282 *
283  READ( nin, fmt = * )nnb
284  IF( nnb.LT.1 ) THEN
285  WRITE( nout, fmt = 9996 )'NNB ', nnb, 1
286  nnb = 0
287  fatal = .true.
288  ELSE IF( nnb.GT.maxin ) THEN
289  WRITE( nout, fmt = 9995 )'NNB ', nnb, maxin
290  nnb = 0
291  fatal = .true.
292  END IF
293  READ( nin, fmt = * )( nbval( i ), i = 1, nnb )
294  DO 40 i = 1, nnb
295  IF( nbval( i ).LT.0 ) THEN
296  WRITE( nout, fmt = 9996 )' NB ', nbval( i ), 0
297  fatal = .true.
298  END IF
299  40 CONTINUE
300  IF( nnb.GT.0 )
301  $ WRITE( nout, fmt = 9993 )'NB ', ( nbval( i ), i = 1, nnb )
302 *
303 * Set NBVAL2 to be the set of unique values of NB
304 *
305  nnb2 = 0
306  DO 60 i = 1, nnb
307  nb = nbval( i )
308  DO 50 j = 1, nnb2
309  IF( nb.EQ.nbval2( j ) )
310  $ GO TO 60
311  50 CONTINUE
312  nnb2 = nnb2 + 1
313  nbval2( nnb2 ) = nb
314  60 CONTINUE
315 *
316 * Read the values of NX
317 *
318  READ( nin, fmt = * )( nxval( i ), i = 1, nnb )
319  DO 70 i = 1, nnb
320  IF( nxval( i ).LT.0 ) THEN
321  WRITE( nout, fmt = 9996 )' NX ', nxval( i ), 0
322  fatal = .true.
323  END IF
324  70 CONTINUE
325  IF( nnb.GT.0 )
326  $ WRITE( nout, fmt = 9993 )'NX ', ( nxval( i ), i = 1, nnb )
327 *
328 * Read the values of RANKVAL
329 *
330  READ( nin, fmt = * )nrank
331  IF( nn.LT.1 ) THEN
332  WRITE( nout, fmt = 9996 )' NRANK ', nrank, 1
333  nrank = 0
334  fatal = .true.
335  ELSE IF( nn.GT.maxin ) THEN
336  WRITE( nout, fmt = 9995 )' NRANK ', nrank, maxin
337  nrank = 0
338  fatal = .true.
339  END IF
340  READ( nin, fmt = * )( rankval( i ), i = 1, nrank )
341  DO i = 1, nrank
342  IF( rankval( i ).LT.0 ) THEN
343  WRITE( nout, fmt = 9996 )' RANK ', rankval( i ), 0
344  fatal = .true.
345  ELSE IF( rankval( i ).GT.100 ) THEN
346  WRITE( nout, fmt = 9995 )' RANK ', rankval( i ), 100
347  fatal = .true.
348  END IF
349  END DO
350  IF( nrank.GT.0 )
351  $ WRITE( nout, fmt = 9993 )'RANK % OF N',
352  $ ( rankval( i ), i = 1, nrank )
353 *
354 * Read the threshold value for the test ratios.
355 *
356  READ( nin, fmt = * )thresh
357  WRITE( nout, fmt = 9992 )thresh
358 *
359 * Read the flag that indicates whether to test the LAPACK routines.
360 *
361  READ( nin, fmt = * )tstchk
362 *
363 * Read the flag that indicates whether to test the driver routines.
364 *
365  READ( nin, fmt = * )tstdrv
366 *
367 * Read the flag that indicates whether to test the error exits.
368 *
369  READ( nin, fmt = * )tsterr
370 *
371  IF( fatal ) THEN
372  WRITE( nout, fmt = 9999 )
373  stop
374  END IF
375 *
376 * Calculate and print the machine dependent constants.
377 *
378  eps = slamch( 'Underflow threshold' )
379  WRITE( nout, fmt = 9991 )'underflow', eps
380  eps = slamch( 'Overflow threshold' )
381  WRITE( nout, fmt = 9991 )'overflow ', eps
382  eps = slamch( 'Epsilon' )
383  WRITE( nout, fmt = 9991 )'precision', eps
384  WRITE( nout, fmt = * )
385  nrhs = nsval( 1 )
386 *
387  80 CONTINUE
388 *
389 * Read a test path and the number of matrix types to use.
390 *
391  READ( nin, fmt = '(A72)', end = 140 )aline
392  path = aline( 1: 3 )
393  nmats = matmax
394  i = 3
395  90 CONTINUE
396  i = i + 1
397  IF( i.GT.72 )
398  $ GO TO 130
399  IF( aline( i: i ).EQ.' ' )
400  $ GO TO 90
401  nmats = 0
402  100 CONTINUE
403  c1 = aline( i: i )
404  DO 110 k = 1, 10
405  IF( c1.EQ.intstr( k: k ) ) THEN
406  ic = k - 1
407  GO TO 120
408  END IF
409  110 CONTINUE
410  GO TO 130
411  120 CONTINUE
412  nmats = nmats*10 + ic
413  i = i + 1
414  IF( i.GT.72 )
415  $ GO TO 130
416  GO TO 100
417  130 CONTINUE
418  c1 = path( 1: 1 )
419  c2 = path( 2: 3 )
420 *
421 * Check first character for correct precision.
422 *
423  IF( .NOT.lsame( c1, 'Complex precision' ) ) THEN
424  WRITE( nout, fmt = 9990 )path
425 *
426  ELSE IF( nmats.LE.0 ) THEN
427 *
428 * Check for a positive number of tests requested.
429 *
430  WRITE( nout, fmt = 9989 )path
431 *
432  ELSE IF( lsamen( 2, c2, 'GE' ) ) THEN
433 *
434 * GE: general matrices
435 *
436  ntypes = 11
437  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
438 *
439  IF( tstchk ) THEN
440  CALL cchkge( dotype, nm, mval, nn, nval, nnb2, nbval2, nns,
441  $ nsval, thresh, tsterr, lda, a( 1, 1 ),
442  $ a( 1, 2 ), a( 1, 3 ), b( 1, 1 ), b( 1, 2 ),
443  $ b( 1, 3 ), work, rwork, iwork, nout )
444  ELSE
445  WRITE( nout, fmt = 9989 )path
446  END IF
447 *
448  IF( tstdrv ) THEN
449  CALL cdrvge( dotype, nn, nval, nrhs, thresh, tsterr, lda,
450  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
451  $ b( 1, 2 ), b( 1, 3 ), b( 1, 4 ), s, work,
452  $ rwork, iwork, nout )
453  ELSE
454  WRITE( nout, fmt = 9988 )path
455  END IF
456 *
457  ELSE IF( lsamen( 2, c2, 'GB' ) ) THEN
458 *
459 * GB: general banded matrices
460 *
461  la = ( 2*kdmax+1 )*nmax
462  lafac = ( 3*kdmax+1 )*nmax
463  ntypes = 8
464  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
465 *
466  IF( tstchk ) THEN
467  CALL cchkgb( dotype, nm, mval, nn, nval, nnb2, nbval2, nns,
468  $ nsval, thresh, tsterr, a( 1, 1 ), la,
469  $ a( 1, 3 ), lafac, b( 1, 1 ), b( 1, 2 ),
470  $ b( 1, 3 ), work, rwork, iwork, nout )
471  ELSE
472  WRITE( nout, fmt = 9989 )path
473  END IF
474 *
475  IF( tstdrv ) THEN
476  CALL cdrvgb( dotype, nn, nval, nrhs, thresh, tsterr,
477  $ a( 1, 1 ), la, a( 1, 3 ), lafac, a( 1, 6 ),
478  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), b( 1, 4 ), s,
479  $ work, rwork, iwork, nout )
480  ELSE
481  WRITE( nout, fmt = 9988 )path
482  END IF
483 *
484  ELSE IF( lsamen( 2, c2, 'GT' ) ) THEN
485 *
486 * GT: general tridiagonal matrices
487 *
488  ntypes = 12
489  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
490 *
491  IF( tstchk ) THEN
492  CALL cchkgt( dotype, nn, nval, nns, nsval, thresh, tsterr,
493  $ a( 1, 1 ), a( 1, 2 ), b( 1, 1 ), b( 1, 2 ),
494  $ b( 1, 3 ), work, rwork, iwork, nout )
495  ELSE
496  WRITE( nout, fmt = 9989 )path
497  END IF
498 *
499  IF( tstdrv ) THEN
500  CALL cdrvgt( dotype, nn, nval, nrhs, thresh, tsterr,
501  $ a( 1, 1 ), a( 1, 2 ), b( 1, 1 ), b( 1, 2 ),
502  $ b( 1, 3 ), work, rwork, iwork, nout )
503  ELSE
504  WRITE( nout, fmt = 9988 )path
505  END IF
506 *
507  ELSE IF( lsamen( 2, c2, 'PO' ) ) THEN
508 *
509 * PO: positive definite matrices
510 *
511  ntypes = 9
512  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
513 *
514  IF( tstchk ) THEN
515  CALL cchkpo( dotype, nn, nval, nnb2, nbval2, nns, nsval,
516  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
517  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
518  $ work, rwork, nout )
519  ELSE
520  WRITE( nout, fmt = 9989 )path
521  END IF
522 *
523  IF( tstdrv ) THEN
524  CALL cdrvpo( dotype, nn, nval, nrhs, thresh, tsterr, lda,
525  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
526  $ b( 1, 2 ), b( 1, 3 ), b( 1, 4 ), s, work,
527  $ rwork, nout )
528  ELSE
529  WRITE( nout, fmt = 9988 )path
530  END IF
531 *
532  ELSE IF( lsamen( 2, c2, 'PS' ) ) THEN
533 *
534 * PS: positive semi-definite matrices
535 *
536  ntypes = 9
537 *
538  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
539 *
540  IF( tstchk ) THEN
541  CALL cchkps( dotype, nn, nval, nnb2, nbval2, nrank,
542  $ rankval, thresh, tsterr, lda, a( 1, 1 ),
543  $ a( 1, 2 ), a( 1, 3 ), piv, work, rwork,
544  $ nout )
545  ELSE
546  WRITE( nout, fmt = 9989 )path
547  END IF
548 *
549  ELSE IF( lsamen( 2, c2, 'PP' ) ) THEN
550 *
551 * PP: positive definite packed matrices
552 *
553  ntypes = 9
554  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
555 *
556  IF( tstchk ) THEN
557  CALL cchkpp( dotype, nn, nval, nns, nsval, thresh, tsterr,
558  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
559  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work, rwork,
560  $ nout )
561  ELSE
562  WRITE( nout, fmt = 9989 )path
563  END IF
564 *
565  IF( tstdrv ) THEN
566  CALL cdrvpp( dotype, nn, nval, nrhs, thresh, tsterr, lda,
567  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
568  $ b( 1, 2 ), b( 1, 3 ), b( 1, 4 ), s, work,
569  $ rwork, nout )
570  ELSE
571  WRITE( nout, fmt = 9988 )path
572  END IF
573 *
574  ELSE IF( lsamen( 2, c2, 'PB' ) ) THEN
575 *
576 * PB: positive definite banded matrices
577 *
578  ntypes = 8
579  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
580 *
581  IF( tstchk ) THEN
582  CALL cchkpb( dotype, nn, nval, nnb2, nbval2, nns, nsval,
583  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
584  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
585  $ work, rwork, nout )
586  ELSE
587  WRITE( nout, fmt = 9989 )path
588  END IF
589 *
590  IF( tstdrv ) THEN
591  CALL cdrvpb( dotype, nn, nval, nrhs, thresh, tsterr, lda,
592  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
593  $ b( 1, 2 ), b( 1, 3 ), b( 1, 4 ), s, work,
594  $ rwork, nout )
595  ELSE
596  WRITE( nout, fmt = 9988 )path
597  END IF
598 *
599  ELSE IF( lsamen( 2, c2, 'PT' ) ) THEN
600 *
601 * PT: positive definite tridiagonal matrices
602 *
603  ntypes = 12
604  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
605 *
606  IF( tstchk ) THEN
607  CALL cchkpt( dotype, nn, nval, nns, nsval, thresh, tsterr,
608  $ a( 1, 1 ), s, a( 1, 2 ), b( 1, 1 ), b( 1, 2 ),
609  $ b( 1, 3 ), work, rwork, nout )
610  ELSE
611  WRITE( nout, fmt = 9989 )path
612  END IF
613 *
614  IF( tstdrv ) THEN
615  CALL cdrvpt( dotype, nn, nval, nrhs, thresh, tsterr,
616  $ a( 1, 1 ), s, a( 1, 2 ), b( 1, 1 ), b( 1, 2 ),
617  $ b( 1, 3 ), work, rwork, nout )
618  ELSE
619  WRITE( nout, fmt = 9988 )path
620  END IF
621 *
622  ELSE IF( lsamen( 2, c2, 'HE' ) ) THEN
623 *
624 * HE: Hermitian indefinite matrices,
625 * with partial (Bunch-Kaufman) pivoting algorithm
626 *
627  ntypes = 10
628  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
629 *
630  IF( tstchk ) THEN
631  CALL cchkhe( dotype, nn, nval, nnb2, nbval2, nns, nsval,
632  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
633  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
634  $ work, rwork, iwork, nout )
635  ELSE
636  WRITE( nout, fmt = 9989 )path
637  END IF
638 *
639  IF( tstdrv ) THEN
640  CALL cdrvhe( dotype, nn, nval, nrhs, thresh, tsterr, lda,
641  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
642  $ b( 1, 2 ), b( 1, 3 ), work, rwork, iwork,
643  $ nout )
644  ELSE
645  WRITE( nout, fmt = 9988 )path
646  END IF
647 *
648  ELSE IF( lsamen( 2, c2, 'HR' ) ) THEN
649 *
650 * HR: Hermitian indefinite matrices,
651 * with bounded Bunch-Kaufman (rook) pivoting algorithm
652 *
653  ntypes = 10
654  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
655 *
656  IF( tstchk ) THEN
657  CALL cchkhe_rook(dotype, nn, nval, nnb2, nbval2, nns, nsval,
658  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
659  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
660  $ work, rwork, iwork, nout )
661  ELSE
662  WRITE( nout, fmt = 9989 )path
663  END IF
664 *
665  IF( tstdrv ) THEN
666  CALL cdrvhe_rook( dotype, nn, nval, nrhs, thresh, tsterr,
667  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
668  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work,
669  $ rwork, iwork, nout )
670  ELSE
671  WRITE( nout, fmt = 9988 )path
672  END IF
673 *
674  ELSE IF( lsamen( 2, c2, 'HK' ) ) THEN
675 *
676 * HK: Hermitian indefinite matrices,
677 * with bounded Bunch-Kaufman (rook) pivoting algorithm,
678 * differnet matrix storage format than HR path version.
679 *
680  ntypes = 10
681  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
682 *
683  IF( tstchk ) THEN
684  CALL cchkhe_rk( dotype, nn, nval, nnb2, nbval2, nns, nsval,
685  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
686  $ e, a( 1, 3 ), b( 1, 1 ), b( 1, 2 ),
687  $ b( 1, 3 ), work, rwork, iwork, nout )
688  ELSE
689  WRITE( nout, fmt = 9989 )path
690  END IF
691 *
692  IF( tstdrv ) THEN
693  CALL cdrvhe_rk( dotype, nn, nval, nrhs, thresh, tsterr,
694  $ lda, a( 1, 1 ), a( 1, 2 ), e, a( 1, 3 ),
695  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work,
696  $ rwork, iwork, nout )
697  ELSE
698  WRITE( nout, fmt = 9988 )path
699  END IF
700 *
701  ELSE IF( lsamen( 2, c2, 'HA' ) ) THEN
702 *
703 * HA: Hermitian matrices,
704 * Aasen Algorithm
705 *
706  ntypes = 10
707  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
708 *
709  IF( tstchk ) THEN
710  CALL cchkhe_aa( dotype, nn, nval, nnb2, nbval2, nns,
711  $ nsval, thresh, tsterr, lda,
712  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
713  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
714  $ work, rwork, iwork, nout )
715  ELSE
716  WRITE( nout, fmt = 9989 )path
717  END IF
718 *
719  IF( tstdrv ) THEN
720  CALL cdrvhe_aa( dotype, nn, nval, nrhs, thresh, tsterr,
721  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
722  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
723  $ work, rwork, iwork, nout )
724  ELSE
725  WRITE( nout, fmt = 9988 )path
726  END IF
727 *
728  ELSE IF( lsamen( 2, c2, 'HP' ) ) THEN
729 *
730 * HP: Hermitian indefinite packed matrices,
731 * with partial (Bunch-Kaufman) pivoting algorithm
732 *
733  ntypes = 10
734  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
735 *
736  IF( tstchk ) THEN
737  CALL cchkhp( dotype, nn, nval, nns, nsval, thresh, tsterr,
738  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
739  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work, rwork,
740  $ iwork, nout )
741  ELSE
742  WRITE( nout, fmt = 9989 )path
743  END IF
744 *
745  IF( tstdrv ) THEN
746  CALL cdrvhp( dotype, nn, nval, nrhs, thresh, tsterr, lda,
747  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
748  $ b( 1, 2 ), b( 1, 3 ), work, rwork, iwork,
749  $ nout )
750  ELSE
751  WRITE( nout, fmt = 9988 )path
752  END IF
753 *
754  ELSE IF( lsamen( 2, c2, 'SY' ) ) THEN
755 *
756 * SY: symmetric indefinite matrices,
757 * with partial (Bunch-Kaufman) pivoting algorithm
758 *
759  ntypes = 11
760  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
761 *
762  IF( tstchk ) THEN
763  CALL cchksy( dotype, nn, nval, nnb2, nbval2, nns, nsval,
764  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
765  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
766  $ work, rwork, iwork, nout )
767  ELSE
768  WRITE( nout, fmt = 9989 )path
769  END IF
770 *
771  IF( tstdrv ) THEN
772  CALL cdrvsy( dotype, nn, nval, nrhs, thresh, tsterr, lda,
773  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
774  $ b( 1, 2 ), b( 1, 3 ), work, rwork, iwork,
775  $ nout )
776  ELSE
777  WRITE( nout, fmt = 9988 )path
778  END IF
779 *
780  ELSE IF( lsamen( 2, c2, 'SR' ) ) THEN
781 *
782 * SR: symmetric indefinite matrices,
783 * with bounded Bunch-Kaufman (rook) pivoting algorithm
784 *
785  ntypes = 11
786  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
787 *
788  IF( tstchk ) THEN
789  CALL cchksy_rook(dotype, nn, nval, nnb2, nbval2, nns, nsval,
790  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
791  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ), b( 1, 3 ),
792  $ work, rwork, iwork, nout )
793  ELSE
794  WRITE( nout, fmt = 9989 )path
795  END IF
796 *
797  IF( tstdrv ) THEN
798  CALL cdrvsy_rook( dotype, nn, nval, nrhs, thresh, tsterr,
799  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
800  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work,
801  $ rwork, iwork, nout )
802  ELSE
803  WRITE( nout, fmt = 9988 )path
804  END IF
805 *
806  ELSE IF( lsamen( 2, c2, 'SK' ) ) THEN
807 *
808 * SK: symmetric indefinite matrices,
809 * with bounded Bunch-Kaufman (rook) pivoting algorithm,
810 * differnet matrix storage format than SR path version.
811 *
812  ntypes = 11
813  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
814 *
815  IF( tstchk ) THEN
816  CALL cchksy_rk( dotype, nn, nval, nnb2, nbval2, nns, nsval,
817  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
818  $ e, a( 1, 3 ), b( 1, 1 ), b( 1, 2 ),
819  $ b( 1, 3 ), work, rwork, iwork, nout )
820  ELSE
821  WRITE( nout, fmt = 9989 )path
822  END IF
823 *
824  IF( tstdrv ) THEN
825  CALL cdrvsy_rk( dotype, nn, nval, nrhs, thresh, tsterr,
826  $ lda, a( 1, 1 ), a( 1, 2 ), e, a( 1, 3 ),
827  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work,
828  $ rwork, iwork, nout )
829  ELSE
830  WRITE( nout, fmt = 9988 )path
831  END IF
832 *
833  ELSE IF( lsamen( 2, c2, 'SA' ) ) THEN
834 *
835 * SA: symmetric indefinite matrices with Aasen's algorithm,
836 *
837  ntypes = 11
838  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
839 *
840  IF( tstchk ) THEN
841  CALL cchksy_aa( dotype, nn, nval, nnb2, nbval2, nns, nsval,
842  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
843  $ a( 1, 3 ), b( 1, 1 ), b( 1, 2 ),
844  $ b( 1, 3 ), work, rwork, iwork, nout )
845  ELSE
846  WRITE( nout, fmt = 9989 )path
847  END IF
848 *
849  IF( tstdrv ) THEN
850  CALL cdrvsy_aa( dotype, nn, nval, nrhs, thresh, tsterr,
851  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
852  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work,
853  $ rwork, iwork, nout )
854  ELSE
855  WRITE( nout, fmt = 9988 )path
856  END IF
857 *
858  ELSE IF( lsamen( 2, c2, 'SP' ) ) THEN
859 *
860 * SP: symmetric indefinite packed matrices,
861 * with partial (Bunch-Kaufman) pivoting algorithm
862 *
863  ntypes = 11
864  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
865 *
866  IF( tstchk ) THEN
867  CALL cchksp( dotype, nn, nval, nns, nsval, thresh, tsterr,
868  $ lda, a( 1, 1 ), a( 1, 2 ), a( 1, 3 ),
869  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work, rwork,
870  $ iwork, nout )
871  ELSE
872  WRITE( nout, fmt = 9989 )path
873  END IF
874 *
875  IF( tstdrv ) THEN
876  CALL cdrvsp( dotype, nn, nval, nrhs, thresh, tsterr, lda,
877  $ a( 1, 1 ), a( 1, 2 ), a( 1, 3 ), b( 1, 1 ),
878  $ b( 1, 2 ), b( 1, 3 ), work, rwork, iwork,
879  $ nout )
880  ELSE
881  WRITE( nout, fmt = 9988 )path
882  END IF
883 *
884  ELSE IF( lsamen( 2, c2, 'TR' ) ) THEN
885 *
886 * TR: triangular matrices
887 *
888  ntypes = 18
889  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
890 *
891  IF( tstchk ) THEN
892  CALL cchktr( dotype, nn, nval, nnb2, nbval2, nns, nsval,
893  $ thresh, tsterr, lda, a( 1, 1 ), a( 1, 2 ),
894  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), work, rwork,
895  $ nout )
896  ELSE
897  WRITE( nout, fmt = 9989 )path
898  END IF
899 *
900  ELSE IF( lsamen( 2, c2, 'TP' ) ) THEN
901 *
902 * TP: triangular packed matrices
903 *
904  ntypes = 18
905  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
906 *
907  IF( tstchk ) THEN
908  CALL cchktp( dotype, nn, nval, nns, nsval, thresh, tsterr,
909  $ lda, a( 1, 1 ), a( 1, 2 ), b( 1, 1 ),
910  $ b( 1, 2 ), b( 1, 3 ), work, rwork, nout )
911  ELSE
912  WRITE( nout, fmt = 9989 )path
913  END IF
914 *
915  ELSE IF( lsamen( 2, c2, 'TB' ) ) THEN
916 *
917 * TB: triangular banded matrices
918 *
919  ntypes = 17
920  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
921 *
922  IF( tstchk ) THEN
923  CALL cchktb( dotype, nn, nval, nns, nsval, thresh, tsterr,
924  $ lda, a( 1, 1 ), a( 1, 2 ), b( 1, 1 ),
925  $ b( 1, 2 ), b( 1, 3 ), work, rwork, nout )
926  ELSE
927  WRITE( nout, fmt = 9989 )path
928  END IF
929 *
930  ELSE IF( lsamen( 2, c2, 'QR' ) ) THEN
931 *
932 * QR: QR factorization
933 *
934  ntypes = 8
935  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
936 *
937  IF( tstchk ) THEN
938  CALL cchkqr( dotype, nm, mval, nn, nval, nnb, nbval, nxval,
939  $ nrhs, thresh, tsterr, nmax, a( 1, 1 ),
940  $ a( 1, 2 ), a( 1, 3 ), a( 1, 4 ), a( 1, 5 ),
941  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), b( 1, 4 ),
942  $ work, rwork, iwork, nout )
943  ELSE
944  WRITE( nout, fmt = 9989 )path
945  END IF
946 *
947  ELSE IF( lsamen( 2, c2, 'LQ' ) ) THEN
948 *
949 * LQ: LQ factorization
950 *
951  ntypes = 8
952  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
953 *
954  IF( tstchk ) THEN
955  CALL cchklq( dotype, nm, mval, nn, nval, nnb, nbval, nxval,
956  $ nrhs, thresh, tsterr, nmax, a( 1, 1 ),
957  $ a( 1, 2 ), a( 1, 3 ), a( 1, 4 ), a( 1, 5 ),
958  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), b( 1, 4 ),
959  $ work, rwork, nout )
960  ELSE
961  WRITE( nout, fmt = 9989 )path
962  END IF
963 *
964  ELSE IF( lsamen( 2, c2, 'QL' ) ) THEN
965 *
966 * QL: QL factorization
967 *
968  ntypes = 8
969  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
970 *
971  IF( tstchk ) THEN
972  CALL cchkql( dotype, nm, mval, nn, nval, nnb, nbval, nxval,
973  $ nrhs, thresh, tsterr, nmax, a( 1, 1 ),
974  $ a( 1, 2 ), a( 1, 3 ), a( 1, 4 ), a( 1, 5 ),
975  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), b( 1, 4 ),
976  $ work, rwork, nout )
977  ELSE
978  WRITE( nout, fmt = 9989 )path
979  END IF
980 *
981  ELSE IF( lsamen( 2, c2, 'RQ' ) ) THEN
982 *
983 * RQ: RQ factorization
984 *
985  ntypes = 8
986  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
987 *
988  IF( tstchk ) THEN
989  CALL cchkrq( dotype, nm, mval, nn, nval, nnb, nbval, nxval,
990  $ nrhs, thresh, tsterr, nmax, a( 1, 1 ),
991  $ a( 1, 2 ), a( 1, 3 ), a( 1, 4 ), a( 1, 5 ),
992  $ b( 1, 1 ), b( 1, 2 ), b( 1, 3 ), b( 1, 4 ),
993  $ work, rwork, iwork, nout )
994  ELSE
995  WRITE( nout, fmt = 9989 )path
996  END IF
997 *
998  ELSE IF( lsamen( 2, c2, 'EQ' ) ) THEN
999 *
1000 * EQ: Equilibration routines for general and positive definite
1001 * matrices (THREQ should be between 2 and 10)
1002 *
1003  IF( tstchk ) THEN
1004  CALL cchkeq( threq, nout )
1005  ELSE
1006  WRITE( nout, fmt = 9989 )path
1007  END IF
1008 *
1009  ELSE IF( lsamen( 2, c2, 'TZ' ) ) THEN
1010 *
1011 * TZ: Trapezoidal matrix
1012 *
1013  ntypes = 3
1014  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
1015 *
1016  IF( tstchk ) THEN
1017  CALL cchktz( dotype, nm, mval, nn, nval, thresh, tsterr,
1018  $ a( 1, 1 ), a( 1, 2 ), s( 1 ),
1019  $ b( 1, 1 ), work, rwork, nout )
1020  ELSE
1021  WRITE( nout, fmt = 9989 )path
1022  END IF
1023 *
1024  ELSE IF( lsamen( 2, c2, 'QP' ) ) THEN
1025 *
1026 * QP: QR factorization with pivoting
1027 *
1028  ntypes = 6
1029  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
1030 *
1031  IF( tstchk ) THEN
1032  CALL cchkq3( dotype, nm, mval, nn, nval, nnb, nbval, nxval,
1033  $ thresh, a( 1, 1 ), a( 1, 2 ), s( 1 ),
1034  $ b( 1, 1 ), work, rwork, iwork, nout )
1035  ELSE
1036  WRITE( nout, fmt = 9989 )path
1037  END IF
1038 *
1039  ELSE IF( lsamen( 2, c2, 'LS' ) ) THEN
1040 *
1041 * LS: Least squares drivers
1042 *
1043  ntypes = 6
1044  CALL alareq( path, nmats, dotype, ntypes, nin, nout )
1045 *
1046  IF( tstdrv ) THEN
1047  CALL cdrvls( dotype, nm, mval, nn, nval, nns, nsval, nnb,
1048  $ nbval, nxval, thresh, tsterr, a( 1, 1 ),
1049  $ a( 1, 2 ), a( 1, 3 ), a( 1, 4 ), a( 1, 5 ),
1050  $ s( 1 ), s( nmax+1 ), nout )
1051  ELSE
1052  WRITE( nout, fmt = 9989 )path
1053  END IF
1054 *
1055  ELSE IF( lsamen( 2, c2, 'QT' ) ) THEN
1056 *
1057 * QT: QRT routines for general matrices
1058 *
1059  IF( tstchk ) THEN
1060  CALL cchkqrt( thresh, tsterr, nm, mval, nn, nval, nnb,
1061  $ nbval, nout )
1062  ELSE
1063  WRITE( nout, fmt = 9989 )path
1064  END IF
1065 *
1066  ELSE IF( lsamen( 2, c2, 'QX' ) ) THEN
1067 *
1068 * QX: QRT routines for triangular-pentagonal matrices
1069 *
1070  IF( tstchk ) THEN
1071  CALL cchkqrtp( thresh, tsterr, nm, mval, nn, nval, nnb,
1072  $ nbval, nout )
1073  ELSE
1074  WRITE( nout, fmt = 9989 )path
1075  END IF
1076 *
1077  ELSE IF( lsamen( 2, c2, 'TQ' ) ) THEN
1078 *
1079 * TQ: LQT routines for general matrices
1080 *
1081  IF( tstchk ) THEN
1082  CALL cchklqt( thresh, tsterr, nm, mval, nn, nval, nnb,
1083  $ nbval, nout )
1084  ELSE
1085  WRITE( nout, fmt = 9989 )path
1086  END IF
1087 *
1088  ELSE IF( lsamen( 2, c2, 'XQ' ) ) THEN
1089 *
1090 * XQ: LQT routines for triangular-pentagonal matrices
1091 *
1092  IF( tstchk ) THEN
1093  CALL cchklqtp( thresh, tsterr, nm, mval, nn, nval, nnb,
1094  $ nbval, nout )
1095  ELSE
1096  WRITE( nout, fmt = 9989 )path
1097  END IF
1098 *
1099  ELSE IF( lsamen( 2, c2, 'TS' ) ) THEN
1100 *
1101 * TS: QR routines for tall-skinny matrices
1102 *
1103  IF( tstchk ) THEN
1104  CALL cchktsqr( thresh, tsterr, nm, mval, nn, nval, nnb,
1105  $ nbval, nout )
1106  ELSE
1107  WRITE( nout, fmt = 9989 )path
1108  END IF
1109 *
1110  ELSE
1111 *
1112  WRITE( nout, fmt = 9990 )path
1113  END IF
1114 *
1115 * Go back to get another input line.
1116 *
1117  GO TO 80
1118 *
1119 * Branch to this line when the last record is read.
1120 *
1121  140 CONTINUE
1122  CLOSE ( nin )
1123  s2 = second( )
1124  WRITE( nout, fmt = 9998 )
1125  WRITE( nout, fmt = 9997 )s2 - s1
1126 *
1127  9999 FORMAT( / ' Execution not attempted due to input errors' )
1128  9998 FORMAT( / ' End of tests' )
1129  9997 FORMAT( ' Total time used = ', f12.2, ' seconds', / )
1130  9996 FORMAT( ' Invalid input value: ', a4, '=', i6, '; must be >=',
1131  $ i6 )
1132  9995 FORMAT( ' Invalid input value: ', a4, '=', i6, '; must be <=',
1133  $ i6 )
1134  9994 FORMAT( ' Tests of the COMPLEX LAPACK routines ',
1135  $ / ' LAPACK VERSION ', i1, '.', i1, '.', i1,
1136  $ / / ' The following parameter values will be used:' )
1137  9993 FORMAT( 4x, a4, ': ', 10i6, / 11x, 10i6 )
1138  9992 FORMAT( / ' Routines pass computational tests if test ratio is ',
1139  $ 'less than', f8.2, / )
1140  9991 FORMAT( ' Relative machine ', a, ' is taken to be', e16.6 )
1141  9990 FORMAT( / 1x, a3, ': Unrecognized path name' )
1142  9989 FORMAT( / 1x, a3, ' routines were not tested' )
1143  9988 FORMAT( / 1x, a3, ' driver routines were not tested' )
1144 *
1145 * End of CCHKAA
1146 *
1147  END
subroutine cdrvgb(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, LA, AFB, LAFB, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, IWORK, NOUT)
CDRVGB
Definition: cdrvgb.f:174
subroutine cdrvhp(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVHP
Definition: cdrvhp.f:159
subroutine cchksy(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKSY
Definition: cchksy.f:173
subroutine cdrvpb(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, NOUT)
CDRVPB
Definition: cdrvpb.f:161
subroutine cdrvsy(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVSY
Definition: cdrvsy.f:155
subroutine cchklqtp(THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB, NBVAL, NOUT)
CCHKLQTP
Definition: cchklqtp.f:104
subroutine cchklqt(THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB, NBVAL, NOUT)
CCHKLQT
Definition: cchklqt.f:104
subroutine cdrvsy_aa(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVSY_AA
Definition: cdrvsy_aa.f:157
subroutine cchktsqr(THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB, NBVAL, NOUT)
CCHKQRT
Definition: cchktsqr.f:104
subroutine cdrvls(DOTYPE, NM, MVAL, NN, NVAL, NNS, NSVAL, NNB, NBVAL, NXVAL, THRESH, TSTERR, A, COPYA, B, COPYB, C, S, COPYS, NOUT)
CDRVLS
Definition: cdrvls.f:194
subroutine cdrvsy_rk(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVSY_RK
Definition: cdrvsy_rk.f:158
subroutine cdrvhe_rk(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVHE_RK
Definition: cdrvhe_rk.f:160
subroutine cdrvhe(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVHE
Definition: cdrvhe.f:155
subroutine cchkpp(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, NOUT)
CCHKPP
Definition: cchkpp.f:161
subroutine cchkeq(THRESH, NOUT)
CCHKEQ
Definition: cchkeq.f:56
subroutine cchkql(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL, NRHS, THRESH, TSTERR, NMAX, A, AF, AQ, AL, AC, B, X, XACT, TAU, WORK, RWORK, NOUT)
CCHKQL
Definition: cchkql.f:198
subroutine cchkhe_rook(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKHE_ROOK
Definition: cchkhe_rook.f:174
subroutine alareq(PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT)
ALAREQ
Definition: alareq.f:92
subroutine cdrvhe_aa(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVHE_AA
Definition: cdrvhe_aa.f:155
subroutine cchktr(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AINV, B, X, XACT, WORK, RWORK, NOUT)
CCHKTR
Definition: cchktr.f:165
subroutine cdrvhe_rook(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVHE_ROOK
Definition: cdrvhe_rook.f:155
program cchkaa
CCHKAA
Definition: cchkaa.f:113
subroutine cdrvsy_rook(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVSY_ROOK
Definition: cdrvsy_rook.f:154
subroutine cchksy_aa(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKSY_AA
Definition: cchksy_aa.f:174
subroutine cchktp(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, AP, AINVP, B, X, XACT, WORK, RWORK, NOUT)
CCHKTP
Definition: cchktp.f:153
subroutine cdrvge(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, IWORK, NOUT)
CDRVGE
Definition: cdrvge.f:166
subroutine cchkhe_aa(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKHE_AA
Definition: cchkhe_aa.f:174
subroutine cchklq(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL, NRHS, THRESH, TSTERR, NMAX, A, AF, AQ, AL, AC, B, X, XACT, TAU, WORK, RWORK, NOUT)
CCHKLQ
Definition: cchklq.f:198
subroutine cchkgb(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, A, LA, AFAC, LAFAC, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKGB
Definition: cchkgb.f:193
subroutine cchkhe_rk(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKHE_RK
Definition: cchkhe_rk.f:179
subroutine cdrvsp(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVSP
Definition: cdrvsp.f:159
subroutine cchksy_rk(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKSY_RK
Definition: cchksy_rk.f:179
subroutine cchkgt(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, A, AF, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKGT
Definition: cchkgt.f:149
subroutine cchkrq(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL, NRHS, THRESH, TSTERR, NMAX, A, AF, AQ, AR, AC, B, X, XACT, TAU, WORK, RWORK, IWORK, NOUT)
CCHKRQ
Definition: cchkrq.f:203
subroutine cchkge(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKGE
Definition: cchkge.f:188
subroutine cdrvgt(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, AF, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CDRVGT
Definition: cdrvgt.f:141
subroutine cchktb(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, AB, AINV, B, X, XACT, WORK, RWORK, NOUT)
CCHKTB
Definition: cchktb.f:151
subroutine cchkqrtp(THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB, NBVAL, NOUT)
CCHKQRTP
Definition: cchkqrtp.f:104
subroutine cdrvpt(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, D, E, B, X, XACT, WORK, RWORK, NOUT)
CDRVPT
Definition: cdrvpt.f:142
subroutine cchkpt(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, A, D, E, B, X, XACT, WORK, RWORK, NOUT)
CCHKPT
Definition: cchkpt.f:149
subroutine cchktz(DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR, A, COPYA, S, TAU, WORK, RWORK, NOUT)
CCHKTZ
Definition: cchktz.f:139
subroutine cchkhe(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKHE
Definition: cchkhe.f:173
subroutine cchkqr(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL, NRHS, THRESH, TSTERR, NMAX, A, AF, AQ, AR, AC, B, X, XACT, TAU, WORK, RWORK, IWORK, NOUT)
CCHKQR
Definition: cchkqr.f:203
subroutine ilaver(VERS_MAJOR, VERS_MINOR, VERS_PATCH)
ILAVER returns the LAPACK version.
Definition: ilaver.f:50
subroutine cchkps(DOTYPE, NN, NVAL, NNB, NBVAL, NRANK, RANKVAL, THRESH, TSTERR, NMAX, A, AFAC, PERM, PIV, WORK, RWORK, NOUT)
CCHKPS
Definition: cchkps.f:156
subroutine cdrvpo(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, NOUT)
CDRVPO
Definition: cdrvpo.f:161
subroutine cdrvpp(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, NOUT)
CDRVPP
Definition: cdrvpp.f:161
subroutine cchkpb(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, NOUT)
CCHKPB
Definition: cchkpb.f:170
subroutine cchkq3(DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL, THRESH, A, COPYA, S, TAU, WORK, RWORK, IWORK, NOUT)
CCHKQ3
Definition: cchkq3.f:160
subroutine cchksy_rook(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKSY_ROOK
Definition: cchksy_rook.f:174
subroutine cchksp(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKSP
Definition: cchksp.f:166
subroutine cchkpo(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, NOUT)
CCHKPO
Definition: cchkpo.f:170
subroutine cchkqrt(THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB, NBVAL, NOUT)
CCHKQRT
Definition: cchkqrt.f:104
subroutine cchkhp(DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
CCHKHP
Definition: cchkhp.f:166