You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

cchkaa.F 42 kB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237
  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,
  171. $ CCHKHE_ROOK, CCHKHE_RK, CCHKHE_AA, CCHKHP,
  172. $ CCHKLQ, CCHKUNHR_COL, CCHKPB, CCHKPO, CCHKPS,
  173. $ CCHKPP, CCHKPT, CCHKQ3, CCHKQL, CCHKQR, CCHKRQ,
  174. $ CCHKSP, CCHKSY, CCHKSY_ROOK, CCHKSY_RK,
  175. $ CCHKSY_AA, CCHKTB, CCHKTP, CCHKTR, CCHKTZ,
  176. $ CDRVGB, CDRVGE, CDRVGT, CDRVHE, CDRVHE_ROOK,
  177. $ CDRVHE_RK, CDRVHE_AA, CDRVHP, CDRVLS, CDRVPB,
  178. $ CDRVPO, CDRVPP, CDRVPT, CDRVSP, CDRVSY,
  179. $ CDRVSY_ROOK, CDRVSY_RK, CDRVSY_AA, ILAVER,
  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