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