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.

dchkee.F 93 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232223322342235223622372238223922402241224222432244224522462247224822492250225122522253225422552256225722582259226022612262226322642265226622672268226922702271227222732274227522762277227822792280228122822283228422852286228722882289229022912292229322942295229622972298229923002301230223032304230523062307230823092310231123122313231423152316231723182319232023212322232323242325232623272328232923302331233223332334233523362337233823392340234123422343234423452346234723482349235023512352235323542355235623572358235923602361236223632364236523662367236823692370237123722373237423752376237723782379238023812382238323842385238623872388238923902391239223932394239523962397239823992400240124022403240424052406240724082409241024112412241324142415241624172418241924202421242224232424242524262427242824292430243124322433243424352436243724382439244024412442244324442445244624472448244924502451245224532454245524562457245824592460246124622463246424652466246724682469247024712472247324742475247624772478247924802481248224832484248524862487248824892490249124922493249424952496249724982499250025012502250325042505250625072508250925102511251225132514251525162517251825192520252125222523252425252526252725282529253025312532253325342535253625372538
  1. *> \brief \b DCHKEE
  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 DCHKEE
  12. *
  13. *
  14. *> \par Purpose:
  15. * =============
  16. *>
  17. *> \verbatim
  18. *>
  19. *> DCHKEE tests the DOUBLE PRECISION LAPACK subroutines for the matrix
  20. *> eigenvalue problem. The test paths in this version are
  21. *>
  22. *> NEP (Nonsymmetric Eigenvalue Problem):
  23. *> Test DGEHRD, DORGHR, DHSEQR, DTREVC, DHSEIN, and DORMHR
  24. *>
  25. *> SEP (Symmetric Eigenvalue Problem):
  26. *> Test DSYTRD, DORGTR, DSTEQR, DSTERF, DSTEIN, DSTEDC,
  27. *> and drivers DSYEV(X), DSBEV(X), DSPEV(X), DSTEV(X),
  28. *> DSYEVD, DSBEVD, DSPEVD, DSTEVD
  29. *>
  30. *> SVD (Singular Value Decomposition):
  31. *> Test DGEBRD, DORGBR, DBDSQR, DBDSDC
  32. *> and the drivers DGESVD, DGESDD
  33. *>
  34. *> DEV (Nonsymmetric Eigenvalue/eigenvector Driver):
  35. *> Test DGEEV
  36. *>
  37. *> DES (Nonsymmetric Schur form Driver):
  38. *> Test DGEES
  39. *>
  40. *> DVX (Nonsymmetric Eigenvalue/eigenvector Expert Driver):
  41. *> Test DGEEVX
  42. *>
  43. *> DSX (Nonsymmetric Schur form Expert Driver):
  44. *> Test DGEESX
  45. *>
  46. *> DGG (Generalized Nonsymmetric Eigenvalue Problem):
  47. *> Test DGGHD3, DGGBAL, DGGBAK, DHGEQZ, and DTGEVC
  48. *>
  49. *> DGS (Generalized Nonsymmetric Schur form Driver):
  50. *> Test DGGES
  51. *>
  52. *> DGV (Generalized Nonsymmetric Eigenvalue/eigenvector Driver):
  53. *> Test DGGEV
  54. *>
  55. *> DGX (Generalized Nonsymmetric Schur form Expert Driver):
  56. *> Test DGGESX
  57. *>
  58. *> DXV (Generalized Nonsymmetric Eigenvalue/eigenvector Expert Driver):
  59. *> Test DGGEVX
  60. *>
  61. *> DSG (Symmetric Generalized Eigenvalue Problem):
  62. *> Test DSYGST, DSYGV, DSYGVD, DSYGVX, DSPGST, DSPGV, DSPGVD,
  63. *> DSPGVX, DSBGST, DSBGV, DSBGVD, and DSBGVX
  64. *>
  65. *> DSB (Symmetric Band Eigenvalue Problem):
  66. *> Test DSBTRD
  67. *>
  68. *> DBB (Band Singular Value Decomposition):
  69. *> Test DGBBRD
  70. *>
  71. *> DEC (Eigencondition estimation):
  72. *> Test DLALN2, DLASY2, DLAEQU, DLAEXC, DTRSYL, DTREXC, DTRSNA,
  73. *> DTRSEN, and DLAQTR
  74. *>
  75. *> DBL (Balancing a general matrix)
  76. *> Test DGEBAL
  77. *>
  78. *> DBK (Back transformation on a balanced matrix)
  79. *> Test DGEBAK
  80. *>
  81. *> DGL (Balancing a matrix pair)
  82. *> Test DGGBAL
  83. *>
  84. *> DGK (Back transformation on a matrix pair)
  85. *> Test DGGBAK
  86. *>
  87. *> GLM (Generalized Linear Regression Model):
  88. *> Tests DGGGLM
  89. *>
  90. *> GQR (Generalized QR and RQ factorizations):
  91. *> Tests DGGQRF and DGGRQF
  92. *>
  93. *> GSV (Generalized Singular Value Decomposition):
  94. *> Tests DGGSVD, DGGSVP, DTGSJA, DLAGS2, DLAPLL, and DLAPMT
  95. *>
  96. *> CSD (CS decomposition):
  97. *> Tests DORCSD
  98. *>
  99. *> LSE (Constrained Linear Least Squares):
  100. *> Tests DGGLSE
  101. *>
  102. *> Each test path has a different set of inputs, but the data sets for
  103. *> the driver routines xEV, xES, xVX, and xSX can be concatenated in a
  104. *> single input file. The first line of input should contain one of the
  105. *> 3-character path names in columns 1-3. The number of remaining lines
  106. *> depends on what is found on the first line.
  107. *>
  108. *> The number of matrix types used in testing is often controllable from
  109. *> the input file. The number of matrix types for each path, and the
  110. *> test routine that describes them, is as follows:
  111. *>
  112. *> Path name(s) Types Test routine
  113. *>
  114. *> DHS or NEP 21 DCHKHS
  115. *> DST or SEP 21 DCHKST (routines)
  116. *> 18 DDRVST (drivers)
  117. *> DBD or SVD 16 DCHKBD (routines)
  118. *> 5 DDRVBD (drivers)
  119. *> DEV 21 DDRVEV
  120. *> DES 21 DDRVES
  121. *> DVX 21 DDRVVX
  122. *> DSX 21 DDRVSX
  123. *> DGG 26 DCHKGG (routines)
  124. *> DGS 26 DDRGES
  125. *> DGX 5 DDRGSX
  126. *> DGV 26 DDRGEV
  127. *> DXV 2 DDRGVX
  128. *> DSG 21 DDRVSG
  129. *> DSB 15 DCHKSB
  130. *> DBB 15 DCHKBB
  131. *> DEC - DCHKEC
  132. *> DBL - DCHKBL
  133. *> DBK - DCHKBK
  134. *> DGL - DCHKGL
  135. *> DGK - DCHKGK
  136. *> GLM 8 DCKGLM
  137. *> GQR 8 DCKGQR
  138. *> GSV 8 DCKGSV
  139. *> CSD 3 DCKCSD
  140. *> LSE 8 DCKLSE
  141. *>
  142. *>-----------------------------------------------------------------------
  143. *>
  144. *> NEP input file:
  145. *>
  146. *> line 2: NN, INTEGER
  147. *> Number of values of N.
  148. *>
  149. *> line 3: NVAL, INTEGER array, dimension (NN)
  150. *> The values for the matrix dimension N.
  151. *>
  152. *> line 4: NPARMS, INTEGER
  153. *> Number of values of the parameters NB, NBMIN, NX, NS, and
  154. *> MAXB.
  155. *>
  156. *> line 5: NBVAL, INTEGER array, dimension (NPARMS)
  157. *> The values for the blocksize NB.
  158. *>
  159. *> line 6: NBMIN, INTEGER array, dimension (NPARMS)
  160. *> The values for the minimum blocksize NBMIN.
  161. *>
  162. *> line 7: NXVAL, INTEGER array, dimension (NPARMS)
  163. *> The values for the crossover point NX.
  164. *>
  165. *> line 8: INMIN, INTEGER array, dimension (NPARMS)
  166. *> LAHQR vs TTQRE crossover point, >= 11
  167. *>
  168. *> line 9: INWIN, INTEGER array, dimension (NPARMS)
  169. *> recommended deflation window size
  170. *>
  171. *> line 10: INIBL, INTEGER array, dimension (NPARMS)
  172. *> nibble crossover point
  173. *>
  174. *> line 11: ISHFTS, INTEGER array, dimension (NPARMS)
  175. *> number of simultaneous shifts)
  176. *>
  177. *> line 12: IACC22, INTEGER array, dimension (NPARMS)
  178. *> select structured matrix multiply: 0, 1 or 2)
  179. *>
  180. *> line 13: THRESH
  181. *> Threshold value for the test ratios. Information will be
  182. *> printed about each test for which the test ratio is greater
  183. *> than or equal to the threshold. To have all of the test
  184. *> ratios printed, use THRESH = 0.0 .
  185. *>
  186. *> line 14: NEWSD, INTEGER
  187. *> A code indicating how to set the random number seed.
  188. *> = 0: Set the seed to a default value before each run
  189. *> = 1: Initialize the seed to a default value only before the
  190. *> first run
  191. *> = 2: Like 1, but use the seed values on the next line
  192. *>
  193. *> If line 14 was 2:
  194. *>
  195. *> line 15: INTEGER array, dimension (4)
  196. *> Four integer values for the random number seed.
  197. *>
  198. *> lines 15-EOF: The remaining lines occur in sets of 1 or 2 and allow
  199. *> the user to specify the matrix types. Each line contains
  200. *> a 3-character path name in columns 1-3, and the number
  201. *> of matrix types must be the first nonblank item in columns
  202. *> 4-80. If the number of matrix types is at least 1 but is
  203. *> less than the maximum number of possible types, a second
  204. *> line will be read to get the numbers of the matrix types to
  205. *> be used. For example,
  206. *> NEP 21
  207. *> requests all of the matrix types for the nonsymmetric
  208. *> eigenvalue problem, while
  209. *> NEP 4
  210. *> 9 10 11 12
  211. *> requests only matrices of type 9, 10, 11, and 12.
  212. *>
  213. *> The valid 3-character path names are 'NEP' or 'SHS' for the
  214. *> nonsymmetric eigenvalue routines.
  215. *>
  216. *>-----------------------------------------------------------------------
  217. *>
  218. *> SEP or DSG input file:
  219. *>
  220. *> line 2: NN, INTEGER
  221. *> Number of values of N.
  222. *>
  223. *> line 3: NVAL, INTEGER array, dimension (NN)
  224. *> The values for the matrix dimension N.
  225. *>
  226. *> line 4: NPARMS, INTEGER
  227. *> Number of values of the parameters NB, NBMIN, and NX.
  228. *>
  229. *> line 5: NBVAL, INTEGER array, dimension (NPARMS)
  230. *> The values for the blocksize NB.
  231. *>
  232. *> line 6: NBMIN, INTEGER array, dimension (NPARMS)
  233. *> The values for the minimum blocksize NBMIN.
  234. *>
  235. *> line 7: NXVAL, INTEGER array, dimension (NPARMS)
  236. *> The values for the crossover point NX.
  237. *>
  238. *> line 8: THRESH
  239. *> Threshold value for the test ratios. Information will be
  240. *> printed about each test for which the test ratio is greater
  241. *> than or equal to the threshold.
  242. *>
  243. *> line 9: TSTCHK, LOGICAL
  244. *> Flag indicating whether or not to test the LAPACK routines.
  245. *>
  246. *> line 10: TSTDRV, LOGICAL
  247. *> Flag indicating whether or not to test the driver routines.
  248. *>
  249. *> line 11: TSTERR, LOGICAL
  250. *> Flag indicating whether or not to test the error exits for
  251. *> the LAPACK routines and driver routines.
  252. *>
  253. *> line 12: NEWSD, INTEGER
  254. *> A code indicating how to set the random number seed.
  255. *> = 0: Set the seed to a default value before each run
  256. *> = 1: Initialize the seed to a default value only before the
  257. *> first run
  258. *> = 2: Like 1, but use the seed values on the next line
  259. *>
  260. *> If line 12 was 2:
  261. *>
  262. *> line 13: INTEGER array, dimension (4)
  263. *> Four integer values for the random number seed.
  264. *>
  265. *> lines 13-EOF: Lines specifying matrix types, as for NEP.
  266. *> The 3-character path names are 'SEP' or 'SST' for the
  267. *> symmetric eigenvalue routines and driver routines, and
  268. *> 'DSG' for the routines for the symmetric generalized
  269. *> eigenvalue problem.
  270. *>
  271. *>-----------------------------------------------------------------------
  272. *>
  273. *> SVD input file:
  274. *>
  275. *> line 2: NN, INTEGER
  276. *> Number of values of M and N.
  277. *>
  278. *> line 3: MVAL, INTEGER array, dimension (NN)
  279. *> The values for the matrix row dimension M.
  280. *>
  281. *> line 4: NVAL, INTEGER array, dimension (NN)
  282. *> The values for the matrix column dimension N.
  283. *>
  284. *> line 5: NPARMS, INTEGER
  285. *> Number of values of the parameter NB, NBMIN, NX, and NRHS.
  286. *>
  287. *> line 6: NBVAL, INTEGER array, dimension (NPARMS)
  288. *> The values for the blocksize NB.
  289. *>
  290. *> line 7: NBMIN, INTEGER array, dimension (NPARMS)
  291. *> The values for the minimum blocksize NBMIN.
  292. *>
  293. *> line 8: NXVAL, INTEGER array, dimension (NPARMS)
  294. *> The values for the crossover point NX.
  295. *>
  296. *> line 9: NSVAL, INTEGER array, dimension (NPARMS)
  297. *> The values for the number of right hand sides NRHS.
  298. *>
  299. *> line 10: THRESH
  300. *> Threshold value for the test ratios. Information will be
  301. *> printed about each test for which the test ratio is greater
  302. *> than or equal to the threshold.
  303. *>
  304. *> line 11: TSTCHK, LOGICAL
  305. *> Flag indicating whether or not to test the LAPACK routines.
  306. *>
  307. *> line 12: TSTDRV, LOGICAL
  308. *> Flag indicating whether or not to test the driver routines.
  309. *>
  310. *> line 13: TSTERR, LOGICAL
  311. *> Flag indicating whether or not to test the error exits for
  312. *> the LAPACK routines and driver routines.
  313. *>
  314. *> line 14: NEWSD, INTEGER
  315. *> A code indicating how to set the random number seed.
  316. *> = 0: Set the seed to a default value before each run
  317. *> = 1: Initialize the seed to a default value only before the
  318. *> first run
  319. *> = 2: Like 1, but use the seed values on the next line
  320. *>
  321. *> If line 14 was 2:
  322. *>
  323. *> line 15: INTEGER array, dimension (4)
  324. *> Four integer values for the random number seed.
  325. *>
  326. *> lines 15-EOF: Lines specifying matrix types, as for NEP.
  327. *> The 3-character path names are 'SVD' or 'SBD' for both the
  328. *> SVD routines and the SVD driver routines.
  329. *>
  330. *>-----------------------------------------------------------------------
  331. *>
  332. *> DEV and DES data files:
  333. *>
  334. *> line 1: 'DEV' or 'DES' in columns 1 to 3.
  335. *>
  336. *> line 2: NSIZES, INTEGER
  337. *> Number of sizes of matrices to use. Should be at least 0
  338. *> and at most 20. If NSIZES = 0, no testing is done
  339. *> (although the remaining 3 lines are still read).
  340. *>
  341. *> line 3: NN, INTEGER array, dimension(NSIZES)
  342. *> Dimensions of matrices to be tested.
  343. *>
  344. *> line 4: NB, NBMIN, NX, NS, NBCOL, INTEGERs
  345. *> These integer parameters determine how blocking is done
  346. *> (see ILAENV for details)
  347. *> NB : block size
  348. *> NBMIN : minimum block size
  349. *> NX : minimum dimension for blocking
  350. *> NS : number of shifts in xHSEQR
  351. *> NBCOL : minimum column dimension for blocking
  352. *>
  353. *> line 5: THRESH, REAL
  354. *> The test threshold against which computed residuals are
  355. *> compared. Should generally be in the range from 10. to 20.
  356. *> If it is 0., all test case data will be printed.
  357. *>
  358. *> line 6: TSTERR, LOGICAL
  359. *> Flag indicating whether or not to test the error exits.
  360. *>
  361. *> line 7: NEWSD, INTEGER
  362. *> A code indicating how to set the random number seed.
  363. *> = 0: Set the seed to a default value before each run
  364. *> = 1: Initialize the seed to a default value only before the
  365. *> first run
  366. *> = 2: Like 1, but use the seed values on the next line
  367. *>
  368. *> If line 7 was 2:
  369. *>
  370. *> line 8: INTEGER array, dimension (4)
  371. *> Four integer values for the random number seed.
  372. *>
  373. *> lines 9 and following: Lines specifying matrix types, as for NEP.
  374. *> The 3-character path name is 'DEV' to test SGEEV, or
  375. *> 'DES' to test SGEES.
  376. *>
  377. *>-----------------------------------------------------------------------
  378. *>
  379. *> The DVX data has two parts. The first part is identical to DEV,
  380. *> and the second part consists of test matrices with precomputed
  381. *> solutions.
  382. *>
  383. *> line 1: 'DVX' in columns 1-3.
  384. *>
  385. *> line 2: NSIZES, INTEGER
  386. *> If NSIZES = 0, no testing of randomly generated examples
  387. *> is done, but any precomputed examples are tested.
  388. *>
  389. *> line 3: NN, INTEGER array, dimension(NSIZES)
  390. *>
  391. *> line 4: NB, NBMIN, NX, NS, NBCOL, INTEGERs
  392. *>
  393. *> line 5: THRESH, REAL
  394. *>
  395. *> line 6: TSTERR, LOGICAL
  396. *>
  397. *> line 7: NEWSD, INTEGER
  398. *>
  399. *> If line 7 was 2:
  400. *>
  401. *> line 8: INTEGER array, dimension (4)
  402. *>
  403. *> lines 9 and following: The first line contains 'DVX' in columns 1-3
  404. *> followed by the number of matrix types, possibly with
  405. *> a second line to specify certain matrix types.
  406. *> If the number of matrix types = 0, no testing of randomly
  407. *> generated examples is done, but any precomputed examples
  408. *> are tested.
  409. *>
  410. *> remaining lines : Each matrix is stored on 1+2*N lines, where N is
  411. *> its dimension. The first line contains the dimension (a
  412. *> single integer). The next N lines contain the matrix, one
  413. *> row per line. The last N lines correspond to each
  414. *> eigenvalue. Each of these last N lines contains 4 real
  415. *> values: the real part of the eigenvalue, the imaginary
  416. *> part of the eigenvalue, the reciprocal condition number of
  417. *> the eigenvalues, and the reciprocal condition number of the
  418. *> eigenvector. The end of data is indicated by dimension N=0.
  419. *> Even if no data is to be tested, there must be at least one
  420. *> line containing N=0.
  421. *>
  422. *>-----------------------------------------------------------------------
  423. *>
  424. *> The DSX data is like DVX. The first part is identical to DEV, and the
  425. *> second part consists of test matrices with precomputed solutions.
  426. *>
  427. *> line 1: 'DSX' in columns 1-3.
  428. *>
  429. *> line 2: NSIZES, INTEGER
  430. *> If NSIZES = 0, no testing of randomly generated examples
  431. *> is done, but any precomputed examples are tested.
  432. *>
  433. *> line 3: NN, INTEGER array, dimension(NSIZES)
  434. *>
  435. *> line 4: NB, NBMIN, NX, NS, NBCOL, INTEGERs
  436. *>
  437. *> line 5: THRESH, REAL
  438. *>
  439. *> line 6: TSTERR, LOGICAL
  440. *>
  441. *> line 7: NEWSD, INTEGER
  442. *>
  443. *> If line 7 was 2:
  444. *>
  445. *> line 8: INTEGER array, dimension (4)
  446. *>
  447. *> lines 9 and following: The first line contains 'DSX' in columns 1-3
  448. *> followed by the number of matrix types, possibly with
  449. *> a second line to specify certain matrix types.
  450. *> If the number of matrix types = 0, no testing of randomly
  451. *> generated examples is done, but any precomputed examples
  452. *> are tested.
  453. *>
  454. *> remaining lines : Each matrix is stored on 3+N lines, where N is its
  455. *> dimension. The first line contains the dimension N and the
  456. *> dimension M of an invariant subspace. The second line
  457. *> contains M integers, identifying the eigenvalues in the
  458. *> invariant subspace (by their position in a list of
  459. *> eigenvalues ordered by increasing real part). The next N
  460. *> lines contain the matrix. The last line contains the
  461. *> reciprocal condition number for the average of the selected
  462. *> eigenvalues, and the reciprocal condition number for the
  463. *> corresponding right invariant subspace. The end of data is
  464. *> indicated by a line containing N=0 and M=0. Even if no data
  465. *> is to be tested, there must be at least one line containing
  466. *> N=0 and M=0.
  467. *>
  468. *>-----------------------------------------------------------------------
  469. *>
  470. *> DGG input file:
  471. *>
  472. *> line 2: NN, INTEGER
  473. *> Number of values of N.
  474. *>
  475. *> line 3: NVAL, INTEGER array, dimension (NN)
  476. *> The values for the matrix dimension N.
  477. *>
  478. *> line 4: NPARMS, INTEGER
  479. *> Number of values of the parameters NB, NBMIN, NS, MAXB, and
  480. *> NBCOL.
  481. *>
  482. *> line 5: NBVAL, INTEGER array, dimension (NPARMS)
  483. *> The values for the blocksize NB.
  484. *>
  485. *> line 6: NBMIN, INTEGER array, dimension (NPARMS)
  486. *> The values for NBMIN, the minimum row dimension for blocks.
  487. *>
  488. *> line 7: NSVAL, INTEGER array, dimension (NPARMS)
  489. *> The values for the number of shifts.
  490. *>
  491. *> line 8: MXBVAL, INTEGER array, dimension (NPARMS)
  492. *> The values for MAXB, used in determining minimum blocksize.
  493. *>
  494. *> line 9: IACC22, INTEGER array, dimension (NPARMS)
  495. *> select structured matrix multiply: 1 or 2)
  496. *>
  497. *> line 10: NBCOL, INTEGER array, dimension (NPARMS)
  498. *> The values for NBCOL, the minimum column dimension for
  499. *> blocks.
  500. *>
  501. *> line 11: THRESH
  502. *> Threshold value for the test ratios. Information will be
  503. *> printed about each test for which the test ratio is greater
  504. *> than or equal to the threshold.
  505. *>
  506. *> line 12: TSTCHK, LOGICAL
  507. *> Flag indicating whether or not to test the LAPACK routines.
  508. *>
  509. *> line 13: TSTDRV, LOGICAL
  510. *> Flag indicating whether or not to test the driver routines.
  511. *>
  512. *> line 14: TSTERR, LOGICAL
  513. *> Flag indicating whether or not to test the error exits for
  514. *> the LAPACK routines and driver routines.
  515. *>
  516. *> line 15: NEWSD, INTEGER
  517. *> A code indicating how to set the random number seed.
  518. *> = 0: Set the seed to a default value before each run
  519. *> = 1: Initialize the seed to a default value only before the
  520. *> first run
  521. *> = 2: Like 1, but use the seed values on the next line
  522. *>
  523. *> If line 15 was 2:
  524. *>
  525. *> line 16: INTEGER array, dimension (4)
  526. *> Four integer values for the random number seed.
  527. *>
  528. *> lines 17-EOF: Lines specifying matrix types, as for NEP.
  529. *> The 3-character path name is 'DGG' for the generalized
  530. *> eigenvalue problem routines and driver routines.
  531. *>
  532. *>-----------------------------------------------------------------------
  533. *>
  534. *> DGS and DGV input files:
  535. *>
  536. *> line 1: 'DGS' or 'DGV' in columns 1 to 3.
  537. *>
  538. *> line 2: NN, INTEGER
  539. *> Number of values of N.
  540. *>
  541. *> line 3: NVAL, INTEGER array, dimension(NN)
  542. *> Dimensions of matrices to be tested.
  543. *>
  544. *> line 4: NB, NBMIN, NX, NS, NBCOL, INTEGERs
  545. *> These integer parameters determine how blocking is done
  546. *> (see ILAENV for details)
  547. *> NB : block size
  548. *> NBMIN : minimum block size
  549. *> NX : minimum dimension for blocking
  550. *> NS : number of shifts in xHGEQR
  551. *> NBCOL : minimum column dimension for blocking
  552. *>
  553. *> line 5: THRESH, REAL
  554. *> The test threshold against which computed residuals are
  555. *> compared. Should generally be in the range from 10. to 20.
  556. *> If it is 0., all test case data will be printed.
  557. *>
  558. *> line 6: TSTERR, LOGICAL
  559. *> Flag indicating whether or not to test the error exits.
  560. *>
  561. *> line 7: NEWSD, INTEGER
  562. *> A code indicating how to set the random number seed.
  563. *> = 0: Set the seed to a default value before each run
  564. *> = 1: Initialize the seed to a default value only before the
  565. *> first run
  566. *> = 2: Like 1, but use the seed values on the next line
  567. *>
  568. *> If line 17 was 2:
  569. *>
  570. *> line 7: INTEGER array, dimension (4)
  571. *> Four integer values for the random number seed.
  572. *>
  573. *> lines 7-EOF: Lines specifying matrix types, as for NEP.
  574. *> The 3-character path name is 'DGS' for the generalized
  575. *> eigenvalue problem routines and driver routines.
  576. *>
  577. *>-----------------------------------------------------------------------
  578. *>
  579. *> DXV input files:
  580. *>
  581. *> line 1: 'DXV' in columns 1 to 3.
  582. *>
  583. *> line 2: N, INTEGER
  584. *> Value of N.
  585. *>
  586. *> line 3: NB, NBMIN, NX, NS, NBCOL, INTEGERs
  587. *> These integer parameters determine how blocking is done
  588. *> (see ILAENV for details)
  589. *> NB : block size
  590. *> NBMIN : minimum block size
  591. *> NX : minimum dimension for blocking
  592. *> NS : number of shifts in xHGEQR
  593. *> NBCOL : minimum column dimension for blocking
  594. *>
  595. *> line 4: THRESH, REAL
  596. *> The test threshold against which computed residuals are
  597. *> compared. Should generally be in the range from 10. to 20.
  598. *> Information will be printed about each test for which the
  599. *> test ratio is greater than or equal to the threshold.
  600. *>
  601. *> line 5: TSTERR, LOGICAL
  602. *> Flag indicating whether or not to test the error exits for
  603. *> the LAPACK routines and driver routines.
  604. *>
  605. *> line 6: NEWSD, INTEGER
  606. *> A code indicating how to set the random number seed.
  607. *> = 0: Set the seed to a default value before each run
  608. *> = 1: Initialize the seed to a default value only before the
  609. *> first run
  610. *> = 2: Like 1, but use the seed values on the next line
  611. *>
  612. *> If line 6 was 2:
  613. *>
  614. *> line 7: INTEGER array, dimension (4)
  615. *> Four integer values for the random number seed.
  616. *>
  617. *> If line 2 was 0:
  618. *>
  619. *> line 7-EOF: Precomputed examples are tested.
  620. *>
  621. *> remaining lines : Each example is stored on 3+2*N lines, where N is
  622. *> its dimension. The first line contains the dimension (a
  623. *> single integer). The next N lines contain the matrix A, one
  624. *> row per line. The next N lines contain the matrix B. The
  625. *> next line contains the reciprocals of the eigenvalue
  626. *> condition numbers. The last line contains the reciprocals of
  627. *> the eigenvector condition numbers. The end of data is
  628. *> indicated by dimension N=0. Even if no data is to be tested,
  629. *> there must be at least one line containing N=0.
  630. *>
  631. *>-----------------------------------------------------------------------
  632. *>
  633. *> DGX input files:
  634. *>
  635. *> line 1: 'DGX' in columns 1 to 3.
  636. *>
  637. *> line 2: N, INTEGER
  638. *> Value of N.
  639. *>
  640. *> line 3: NB, NBMIN, NX, NS, NBCOL, INTEGERs
  641. *> These integer parameters determine how blocking is done
  642. *> (see ILAENV for details)
  643. *> NB : block size
  644. *> NBMIN : minimum block size
  645. *> NX : minimum dimension for blocking
  646. *> NS : number of shifts in xHGEQR
  647. *> NBCOL : minimum column dimension for blocking
  648. *>
  649. *> line 4: THRESH, REAL
  650. *> The test threshold against which computed residuals are
  651. *> compared. Should generally be in the range from 10. to 20.
  652. *> Information will be printed about each test for which the
  653. *> test ratio is greater than or equal to the threshold.
  654. *>
  655. *> line 5: TSTERR, LOGICAL
  656. *> Flag indicating whether or not to test the error exits for
  657. *> the LAPACK routines and driver routines.
  658. *>
  659. *> line 6: NEWSD, INTEGER
  660. *> A code indicating how to set the random number seed.
  661. *> = 0: Set the seed to a default value before each run
  662. *> = 1: Initialize the seed to a default value only before the
  663. *> first run
  664. *> = 2: Like 1, but use the seed values on the next line
  665. *>
  666. *> If line 6 was 2:
  667. *>
  668. *> line 7: INTEGER array, dimension (4)
  669. *> Four integer values for the random number seed.
  670. *>
  671. *> If line 2 was 0:
  672. *>
  673. *> line 7-EOF: Precomputed examples are tested.
  674. *>
  675. *> remaining lines : Each example is stored on 3+2*N lines, where N is
  676. *> its dimension. The first line contains the dimension (a
  677. *> single integer). The next line contains an integer k such
  678. *> that only the last k eigenvalues will be selected and appear
  679. *> in the leading diagonal blocks of $A$ and $B$. The next N
  680. *> lines contain the matrix A, one row per line. The next N
  681. *> lines contain the matrix B. The last line contains the
  682. *> reciprocal of the eigenvalue cluster condition number and the
  683. *> reciprocal of the deflating subspace (associated with the
  684. *> selected eigencluster) condition number. The end of data is
  685. *> indicated by dimension N=0. Even if no data is to be tested,
  686. *> there must be at least one line containing N=0.
  687. *>
  688. *>-----------------------------------------------------------------------
  689. *>
  690. *> DSB input file:
  691. *>
  692. *> line 2: NN, INTEGER
  693. *> Number of values of N.
  694. *>
  695. *> line 3: NVAL, INTEGER array, dimension (NN)
  696. *> The values for the matrix dimension N.
  697. *>
  698. *> line 4: NK, INTEGER
  699. *> Number of values of K.
  700. *>
  701. *> line 5: KVAL, INTEGER array, dimension (NK)
  702. *> The values for the matrix dimension K.
  703. *>
  704. *> line 6: THRESH
  705. *> Threshold value for the test ratios. Information will be
  706. *> printed about each test for which the test ratio is greater
  707. *> than or equal to the threshold.
  708. *>
  709. *> line 7: NEWSD, INTEGER
  710. *> A code indicating how to set the random number seed.
  711. *> = 0: Set the seed to a default value before each run
  712. *> = 1: Initialize the seed to a default value only before the
  713. *> first run
  714. *> = 2: Like 1, but use the seed values on the next line
  715. *>
  716. *> If line 7 was 2:
  717. *>
  718. *> line 8: INTEGER array, dimension (4)
  719. *> Four integer values for the random number seed.
  720. *>
  721. *> lines 8-EOF: Lines specifying matrix types, as for NEP.
  722. *> The 3-character path name is 'DSB'.
  723. *>
  724. *>-----------------------------------------------------------------------
  725. *>
  726. *> DBB input file:
  727. *>
  728. *> line 2: NN, INTEGER
  729. *> Number of values of M and N.
  730. *>
  731. *> line 3: MVAL, INTEGER array, dimension (NN)
  732. *> The values for the matrix row dimension M.
  733. *>
  734. *> line 4: NVAL, INTEGER array, dimension (NN)
  735. *> The values for the matrix column dimension N.
  736. *>
  737. *> line 4: NK, INTEGER
  738. *> Number of values of K.
  739. *>
  740. *> line 5: KVAL, INTEGER array, dimension (NK)
  741. *> The values for the matrix bandwidth K.
  742. *>
  743. *> line 6: NPARMS, INTEGER
  744. *> Number of values of the parameter NRHS
  745. *>
  746. *> line 7: NSVAL, INTEGER array, dimension (NPARMS)
  747. *> The values for the number of right hand sides NRHS.
  748. *>
  749. *> line 8: THRESH
  750. *> Threshold value for the test ratios. Information will be
  751. *> printed about each test for which the test ratio is greater
  752. *> than or equal to the threshold.
  753. *>
  754. *> line 9: NEWSD, INTEGER
  755. *> A code indicating how to set the random number seed.
  756. *> = 0: Set the seed to a default value before each run
  757. *> = 1: Initialize the seed to a default value only before the
  758. *> first run
  759. *> = 2: Like 1, but use the seed values on the next line
  760. *>
  761. *> If line 9 was 2:
  762. *>
  763. *> line 10: INTEGER array, dimension (4)
  764. *> Four integer values for the random number seed.
  765. *>
  766. *> lines 10-EOF: Lines specifying matrix types, as for SVD.
  767. *> The 3-character path name is 'DBB'.
  768. *>
  769. *>-----------------------------------------------------------------------
  770. *>
  771. *> DEC input file:
  772. *>
  773. *> line 2: THRESH, REAL
  774. *> Threshold value for the test ratios. Information will be
  775. *> printed about each test for which the test ratio is greater
  776. *> than or equal to the threshold.
  777. *>
  778. *> lines 3-EOF:
  779. *>
  780. *> Input for testing the eigencondition routines consists of a set of
  781. *> specially constructed test cases and their solutions. The data
  782. *> format is not intended to be modified by the user.
  783. *>
  784. *>-----------------------------------------------------------------------
  785. *>
  786. *> DBL and DBK input files:
  787. *>
  788. *> line 1: 'DBL' in columns 1-3 to test SGEBAL, or 'DBK' in
  789. *> columns 1-3 to test SGEBAK.
  790. *>
  791. *> The remaining lines consist of specially constructed test cases.
  792. *>
  793. *>-----------------------------------------------------------------------
  794. *>
  795. *> DGL and DGK input files:
  796. *>
  797. *> line 1: 'DGL' in columns 1-3 to test DGGBAL, or 'DGK' in
  798. *> columns 1-3 to test DGGBAK.
  799. *>
  800. *> The remaining lines consist of specially constructed test cases.
  801. *>
  802. *>-----------------------------------------------------------------------
  803. *>
  804. *> GLM data file:
  805. *>
  806. *> line 1: 'GLM' in columns 1 to 3.
  807. *>
  808. *> line 2: NN, INTEGER
  809. *> Number of values of M, P, and N.
  810. *>
  811. *> line 3: MVAL, INTEGER array, dimension(NN)
  812. *> Values of M (row dimension).
  813. *>
  814. *> line 4: PVAL, INTEGER array, dimension(NN)
  815. *> Values of P (row dimension).
  816. *>
  817. *> line 5: NVAL, INTEGER array, dimension(NN)
  818. *> Values of N (column dimension), note M <= N <= M+P.
  819. *>
  820. *> line 6: THRESH, REAL
  821. *> Threshold value for the test ratios. Information will be
  822. *> printed about each test for which the test ratio is greater
  823. *> than or equal to the threshold.
  824. *>
  825. *> line 7: TSTERR, LOGICAL
  826. *> Flag indicating whether or not to test the error exits for
  827. *> the LAPACK routines and driver routines.
  828. *>
  829. *> line 8: NEWSD, INTEGER
  830. *> A code indicating how to set the random number seed.
  831. *> = 0: Set the seed to a default value before each run
  832. *> = 1: Initialize the seed to a default value only before the
  833. *> first run
  834. *> = 2: Like 1, but use the seed values on the next line
  835. *>
  836. *> If line 8 was 2:
  837. *>
  838. *> line 9: INTEGER array, dimension (4)
  839. *> Four integer values for the random number seed.
  840. *>
  841. *> lines 9-EOF: Lines specifying matrix types, as for NEP.
  842. *> The 3-character path name is 'GLM' for the generalized
  843. *> linear regression model routines.
  844. *>
  845. *>-----------------------------------------------------------------------
  846. *>
  847. *> GQR data file:
  848. *>
  849. *> line 1: 'GQR' in columns 1 to 3.
  850. *>
  851. *> line 2: NN, INTEGER
  852. *> Number of values of M, P, and N.
  853. *>
  854. *> line 3: MVAL, INTEGER array, dimension(NN)
  855. *> Values of M.
  856. *>
  857. *> line 4: PVAL, INTEGER array, dimension(NN)
  858. *> Values of P.
  859. *>
  860. *> line 5: NVAL, INTEGER array, dimension(NN)
  861. *> Values of N.
  862. *>
  863. *> line 6: THRESH, REAL
  864. *> Threshold value for the test ratios. Information will be
  865. *> printed about each test for which the test ratio is greater
  866. *> than or equal to the threshold.
  867. *>
  868. *> line 7: TSTERR, LOGICAL
  869. *> Flag indicating whether or not to test the error exits for
  870. *> the LAPACK routines and driver routines.
  871. *>
  872. *> line 8: NEWSD, INTEGER
  873. *> A code indicating how to set the random number seed.
  874. *> = 0: Set the seed to a default value before each run
  875. *> = 1: Initialize the seed to a default value only before the
  876. *> first run
  877. *> = 2: Like 1, but use the seed values on the next line
  878. *>
  879. *> If line 8 was 2:
  880. *>
  881. *> line 9: INTEGER array, dimension (4)
  882. *> Four integer values for the random number seed.
  883. *>
  884. *> lines 9-EOF: Lines specifying matrix types, as for NEP.
  885. *> The 3-character path name is 'GQR' for the generalized
  886. *> QR and RQ routines.
  887. *>
  888. *>-----------------------------------------------------------------------
  889. *>
  890. *> GSV data file:
  891. *>
  892. *> line 1: 'GSV' in columns 1 to 3.
  893. *>
  894. *> line 2: NN, INTEGER
  895. *> Number of values of M, P, and N.
  896. *>
  897. *> line 3: MVAL, INTEGER array, dimension(NN)
  898. *> Values of M (row dimension).
  899. *>
  900. *> line 4: PVAL, INTEGER array, dimension(NN)
  901. *> Values of P (row dimension).
  902. *>
  903. *> line 5: NVAL, INTEGER array, dimension(NN)
  904. *> Values of N (column dimension).
  905. *>
  906. *> line 6: THRESH, REAL
  907. *> Threshold value for the test ratios. Information will be
  908. *> printed about each test for which the test ratio is greater
  909. *> than or equal to the threshold.
  910. *>
  911. *> line 7: TSTERR, LOGICAL
  912. *> Flag indicating whether or not to test the error exits for
  913. *> the LAPACK routines and driver routines.
  914. *>
  915. *> line 8: NEWSD, INTEGER
  916. *> A code indicating how to set the random number seed.
  917. *> = 0: Set the seed to a default value before each run
  918. *> = 1: Initialize the seed to a default value only before the
  919. *> first run
  920. *> = 2: Like 1, but use the seed values on the next line
  921. *>
  922. *> If line 8 was 2:
  923. *>
  924. *> line 9: INTEGER array, dimension (4)
  925. *> Four integer values for the random number seed.
  926. *>
  927. *> lines 9-EOF: Lines specifying matrix types, as for NEP.
  928. *> The 3-character path name is 'GSV' for the generalized
  929. *> SVD routines.
  930. *>
  931. *>-----------------------------------------------------------------------
  932. *>
  933. *> CSD data file:
  934. *>
  935. *> line 1: 'CSD' in columns 1 to 3.
  936. *>
  937. *> line 2: NM, INTEGER
  938. *> Number of values of M, P, and N.
  939. *>
  940. *> line 3: MVAL, INTEGER array, dimension(NM)
  941. *> Values of M (row and column dimension of orthogonal matrix).
  942. *>
  943. *> line 4: PVAL, INTEGER array, dimension(NM)
  944. *> Values of P (row dimension of top-left block).
  945. *>
  946. *> line 5: NVAL, INTEGER array, dimension(NM)
  947. *> Values of N (column dimension of top-left block).
  948. *>
  949. *> line 6: THRESH, REAL
  950. *> Threshold value for the test ratios. Information will be
  951. *> printed about each test for which the test ratio is greater
  952. *> than or equal to the threshold.
  953. *>
  954. *> line 7: TSTERR, LOGICAL
  955. *> Flag indicating whether or not to test the error exits for
  956. *> the LAPACK routines and driver routines.
  957. *>
  958. *> line 8: NEWSD, INTEGER
  959. *> A code indicating how to set the random number seed.
  960. *> = 0: Set the seed to a default value before each run
  961. *> = 1: Initialize the seed to a default value only before the
  962. *> first run
  963. *> = 2: Like 1, but use the seed values on the next line
  964. *>
  965. *> If line 8 was 2:
  966. *>
  967. *> line 9: INTEGER array, dimension (4)
  968. *> Four integer values for the random number seed.
  969. *>
  970. *> lines 9-EOF: Lines specifying matrix types, as for NEP.
  971. *> The 3-character path name is 'CSD' for the CSD routine.
  972. *>
  973. *>-----------------------------------------------------------------------
  974. *>
  975. *> LSE data file:
  976. *>
  977. *> line 1: 'LSE' in columns 1 to 3.
  978. *>
  979. *> line 2: NN, INTEGER
  980. *> Number of values of M, P, and N.
  981. *>
  982. *> line 3: MVAL, INTEGER array, dimension(NN)
  983. *> Values of M.
  984. *>
  985. *> line 4: PVAL, INTEGER array, dimension(NN)
  986. *> Values of P.
  987. *>
  988. *> line 5: NVAL, INTEGER array, dimension(NN)
  989. *> Values of N, note P <= N <= P+M.
  990. *>
  991. *> line 6: THRESH, REAL
  992. *> Threshold value for the test ratios. Information will be
  993. *> printed about each test for which the test ratio is greater
  994. *> than or equal to the threshold.
  995. *>
  996. *> line 7: TSTERR, LOGICAL
  997. *> Flag indicating whether or not to test the error exits for
  998. *> the LAPACK routines and driver routines.
  999. *>
  1000. *> line 8: NEWSD, INTEGER
  1001. *> A code indicating how to set the random number seed.
  1002. *> = 0: Set the seed to a default value before each run
  1003. *> = 1: Initialize the seed to a default value only before the
  1004. *> first run
  1005. *> = 2: Like 1, but use the seed values on the next line
  1006. *>
  1007. *> If line 8 was 2:
  1008. *>
  1009. *> line 9: INTEGER array, dimension (4)
  1010. *> Four integer values for the random number seed.
  1011. *>
  1012. *> lines 9-EOF: Lines specifying matrix types, as for NEP.
  1013. *> The 3-character path name is 'GSV' for the generalized
  1014. *> SVD routines.
  1015. *>
  1016. *>-----------------------------------------------------------------------
  1017. *>
  1018. *> NMAX is currently set to 132 and must be at least 12 for some of the
  1019. *> precomputed examples, and LWORK = NMAX*(5*NMAX+5)+1 in the parameter
  1020. *> statements below. For SVD, we assume NRHS may be as big as N. The
  1021. *> parameter NEED is set to 14 to allow for 14 N-by-N matrices for DGG.
  1022. *> \endverbatim
  1023. *
  1024. * Arguments:
  1025. * ==========
  1026. *
  1027. *
  1028. * Authors:
  1029. * ========
  1030. *
  1031. *> \author Univ. of Tennessee
  1032. *> \author Univ. of California Berkeley
  1033. *> \author Univ. of Colorado Denver
  1034. *> \author NAG Ltd.
  1035. *
  1036. *> \date June 2016
  1037. *
  1038. *> \ingroup double_eig
  1039. *
  1040. * =====================================================================
  1041. PROGRAM DCHKEE
  1042. *
  1043. #if defined(_OPENMP)
  1044. use omp_lib
  1045. #endif
  1046. *
  1047. * -- LAPACK test routine (version 3.7.0) --
  1048. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  1049. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  1050. * June 2016
  1051. *
  1052. * =====================================================================
  1053. *
  1054. * .. Parameters ..
  1055. INTEGER NMAX
  1056. PARAMETER ( NMAX = 132 )
  1057. INTEGER NCMAX
  1058. PARAMETER ( NCMAX = 20 )
  1059. INTEGER NEED
  1060. PARAMETER ( NEED = 14 )
  1061. INTEGER LWORK
  1062. PARAMETER ( LWORK = NMAX*( 5*NMAX+5 )+1 )
  1063. INTEGER LIWORK
  1064. PARAMETER ( LIWORK = NMAX*( 5*NMAX+20 ) )
  1065. INTEGER MAXIN
  1066. PARAMETER ( MAXIN = 20 )
  1067. INTEGER MAXT
  1068. PARAMETER ( MAXT = 30 )
  1069. INTEGER NIN, NOUT
  1070. PARAMETER ( NIN = 5, NOUT = 6 )
  1071. * ..
  1072. * .. Local Scalars ..
  1073. LOGICAL CSD, DBB, DGG, DSB, FATAL, GLM, GQR, GSV, LSE,
  1074. $ NEP, DBK, DBL, SEP, DES, DEV, DGK, DGL, DGS,
  1075. $ DGV, DGX, DSX, SVD, DVX, DXV, TSTCHK, TSTDIF,
  1076. $ TSTDRV, TSTERR
  1077. CHARACTER C1
  1078. CHARACTER*3 C3, PATH
  1079. CHARACTER*32 VNAME
  1080. CHARACTER*10 INTSTR
  1081. CHARACTER*80 LINE
  1082. INTEGER I, I1, IC, INFO, ITMP, K, LENP, MAXTYP, NEWSD,
  1083. $ NK, NN, NPARMS, NRHS, NTYPES,
  1084. $ VERS_MAJOR, VERS_MINOR, VERS_PATCH, N_THREADS
  1085. DOUBLE PRECISION EPS, S1, S2, THRESH, THRSHN
  1086. * ..
  1087. * .. Local Arrays ..
  1088. LOGICAL DOTYPE( MAXT ), LOGWRK( NMAX )
  1089. INTEGER IOLDSD( 4 ), ISEED( 4 ), IWORK( LIWORK ),
  1090. $ KVAL( MAXIN ), MVAL( MAXIN ), MXBVAL( MAXIN ),
  1091. $ NBCOL( MAXIN ), NBMIN( MAXIN ), NBVAL( MAXIN ),
  1092. $ NSVAL( MAXIN ), NVAL( MAXIN ), NXVAL( MAXIN ),
  1093. $ PVAL( MAXIN )
  1094. INTEGER INMIN( MAXIN ), INWIN( MAXIN ), INIBL( MAXIN ),
  1095. $ ISHFTS( MAXIN ), IACC22( MAXIN )
  1096. DOUBLE PRECISION D( NMAX, 12 ), RESULT( 500 ), TAUA( NMAX ),
  1097. $ TAUB( NMAX ), X( 5*NMAX )
  1098. * ..
  1099. * .. Allocatable Arrays ..
  1100. INTEGER AllocateStatus
  1101. DOUBLE PRECISION, DIMENSION(:), ALLOCATABLE :: WORK
  1102. DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE :: A, B, C
  1103. * ..
  1104. * .. External Functions ..
  1105. LOGICAL LSAMEN
  1106. DOUBLE PRECISION DLAMCH, DSECND
  1107. EXTERNAL LSAMEN, DLAMCH, DSECND
  1108. * ..
  1109. * .. External Subroutines ..
  1110. EXTERNAL ALAREQ, DCHKBB, DCHKBD, DCHKBK, DCHKBL, DCHKEC,
  1111. $ DCHKGG, DCHKGK, DCHKGL, DCHKHS, DCHKSB, DCHKST,
  1112. $ DCKCSD, DCKGLM, DCKGQR, DCKGSV, DCKLSE, DDRGES,
  1113. $ DDRGEV, DDRGSX, DDRGVX, DDRVBD, DDRVES, DDRVEV,
  1114. $ DDRVSG, DDRVST, DDRVSX, DDRVVX, DERRBD,
  1115. $ DERRED, DERRGG, DERRHS, DERRST, ILAVER, XLAENV,
  1116. $ DDRGES3, DDRGEV3,
  1117. $ DCHKST2STG, DDRVST2STG, DCHKSB2STG, DDRVSG2STG
  1118. * ..
  1119. * .. Intrinsic Functions ..
  1120. INTRINSIC LEN, MIN
  1121. * ..
  1122. * .. Scalars in Common ..
  1123. LOGICAL LERR, OK
  1124. CHARACTER*32 SRNAMT
  1125. INTEGER INFOT, MAXB, NPROC, NSHIFT, NUNIT, SELDIM,
  1126. $ SELOPT
  1127. * ..
  1128. * .. Arrays in Common ..
  1129. LOGICAL SELVAL( 20 )
  1130. INTEGER IPARMS( 100 )
  1131. DOUBLE PRECISION SELWI( 20 ), SELWR( 20 )
  1132. * ..
  1133. * .. Common blocks ..
  1134. COMMON / CENVIR / NPROC, NSHIFT, MAXB
  1135. COMMON / INFOC / INFOT, NUNIT, OK, LERR
  1136. COMMON / SRNAMC / SRNAMT
  1137. COMMON / SSLCT / SELOPT, SELDIM, SELVAL, SELWR, SELWI
  1138. COMMON / CLAENV / IPARMS
  1139. * ..
  1140. * .. Data statements ..
  1141. DATA INTSTR / '0123456789' /
  1142. DATA IOLDSD / 0, 0, 0, 1 /
  1143. * ..
  1144. * .. Allocate memory dynamically ..
  1145. *
  1146. ALLOCATE ( A(NMAX*NMAX,NEED), STAT = AllocateStatus )
  1147. IF (AllocateStatus /= 0) STOP "*** Not enough memory ***"
  1148. ALLOCATE ( B(NMAX*NMAX,5), STAT = AllocateStatus )
  1149. IF (AllocateStatus /= 0) STOP "*** Not enough memory ***"
  1150. ALLOCATE ( C(NCMAX*NCMAX,NCMAX*NCMAX), STAT = AllocateStatus )
  1151. IF (AllocateStatus /= 0) STOP "*** Not enough memory ***"
  1152. ALLOCATE ( WORK(LWORK), STAT = AllocateStatus )
  1153. IF (AllocateStatus /= 0) STOP "*** Not enough memory ***"
  1154. * ..
  1155. * .. Executable Statements ..
  1156. *
  1157. A = 0.0
  1158. B = 0.0
  1159. C = 0.0
  1160. D = 0.0
  1161. S1 = DSECND( )
  1162. FATAL = .FALSE.
  1163. NUNIT = NOUT
  1164. *
  1165. * Return to here to read multiple sets of data
  1166. *
  1167. 10 CONTINUE
  1168. *
  1169. * Read the first line and set the 3-character test path
  1170. *
  1171. READ( NIN, FMT = '(A80)', END = 380 )LINE
  1172. PATH = LINE( 1: 3 )
  1173. NEP = LSAMEN( 3, PATH, 'NEP' ) .OR. LSAMEN( 3, PATH, 'DHS' )
  1174. SEP = LSAMEN( 3, PATH, 'SEP' ) .OR. LSAMEN( 3, PATH, 'DST' ) .OR.
  1175. $ LSAMEN( 3, PATH, 'DSG' ) .OR. LSAMEN( 3, PATH, 'SE2' )
  1176. SVD = LSAMEN( 3, PATH, 'SVD' ) .OR. LSAMEN( 3, PATH, 'DBD' )
  1177. DEV = LSAMEN( 3, PATH, 'DEV' )
  1178. DES = LSAMEN( 3, PATH, 'DES' )
  1179. DVX = LSAMEN( 3, PATH, 'DVX' )
  1180. DSX = LSAMEN( 3, PATH, 'DSX' )
  1181. DGG = LSAMEN( 3, PATH, 'DGG' )
  1182. DGS = LSAMEN( 3, PATH, 'DGS' )
  1183. DGX = LSAMEN( 3, PATH, 'DGX' )
  1184. DGV = LSAMEN( 3, PATH, 'DGV' )
  1185. DXV = LSAMEN( 3, PATH, 'DXV' )
  1186. DSB = LSAMEN( 3, PATH, 'DSB' )
  1187. DBB = LSAMEN( 3, PATH, 'DBB' )
  1188. GLM = LSAMEN( 3, PATH, 'GLM' )
  1189. GQR = LSAMEN( 3, PATH, 'GQR' ) .OR. LSAMEN( 3, PATH, 'GRQ' )
  1190. GSV = LSAMEN( 3, PATH, 'GSV' )
  1191. CSD = LSAMEN( 3, PATH, 'CSD' )
  1192. LSE = LSAMEN( 3, PATH, 'LSE' )
  1193. DBL = LSAMEN( 3, PATH, 'DBL' )
  1194. DBK = LSAMEN( 3, PATH, 'DBK' )
  1195. DGL = LSAMEN( 3, PATH, 'DGL' )
  1196. DGK = LSAMEN( 3, PATH, 'DGK' )
  1197. *
  1198. * Report values of parameters.
  1199. *
  1200. IF( PATH.EQ.' ' ) THEN
  1201. GO TO 10
  1202. ELSE IF( NEP ) THEN
  1203. WRITE( NOUT, FMT = 9987 )
  1204. ELSE IF( SEP ) THEN
  1205. WRITE( NOUT, FMT = 9986 )
  1206. ELSE IF( SVD ) THEN
  1207. WRITE( NOUT, FMT = 9985 )
  1208. ELSE IF( DEV ) THEN
  1209. WRITE( NOUT, FMT = 9979 )
  1210. ELSE IF( DES ) THEN
  1211. WRITE( NOUT, FMT = 9978 )
  1212. ELSE IF( DVX ) THEN
  1213. WRITE( NOUT, FMT = 9977 )
  1214. ELSE IF( DSX ) THEN
  1215. WRITE( NOUT, FMT = 9976 )
  1216. ELSE IF( DGG ) THEN
  1217. WRITE( NOUT, FMT = 9975 )
  1218. ELSE IF( DGS ) THEN
  1219. WRITE( NOUT, FMT = 9964 )
  1220. ELSE IF( DGX ) THEN
  1221. WRITE( NOUT, FMT = 9965 )
  1222. ELSE IF( DGV ) THEN
  1223. WRITE( NOUT, FMT = 9963 )
  1224. ELSE IF( DXV ) THEN
  1225. WRITE( NOUT, FMT = 9962 )
  1226. ELSE IF( DSB ) THEN
  1227. WRITE( NOUT, FMT = 9974 )
  1228. ELSE IF( DBB ) THEN
  1229. WRITE( NOUT, FMT = 9967 )
  1230. ELSE IF( GLM ) THEN
  1231. WRITE( NOUT, FMT = 9971 )
  1232. ELSE IF( GQR ) THEN
  1233. WRITE( NOUT, FMT = 9970 )
  1234. ELSE IF( GSV ) THEN
  1235. WRITE( NOUT, FMT = 9969 )
  1236. ELSE IF( CSD ) THEN
  1237. WRITE( NOUT, FMT = 9960 )
  1238. ELSE IF( LSE ) THEN
  1239. WRITE( NOUT, FMT = 9968 )
  1240. ELSE IF( DBL ) THEN
  1241. *
  1242. * DGEBAL: Balancing
  1243. *
  1244. CALL DCHKBL( NIN, NOUT )
  1245. GO TO 10
  1246. ELSE IF( DBK ) THEN
  1247. *
  1248. * DGEBAK: Back transformation
  1249. *
  1250. CALL DCHKBK( NIN, NOUT )
  1251. GO TO 10
  1252. ELSE IF( DGL ) THEN
  1253. *
  1254. * DGGBAL: Balancing
  1255. *
  1256. CALL DCHKGL( NIN, NOUT )
  1257. GO TO 10
  1258. ELSE IF( DGK ) THEN
  1259. *
  1260. * DGGBAK: Back transformation
  1261. *
  1262. CALL DCHKGK( NIN, NOUT )
  1263. GO TO 10
  1264. ELSE IF( LSAMEN( 3, PATH, 'DEC' ) ) THEN
  1265. *
  1266. * DEC: Eigencondition estimation
  1267. *
  1268. READ( NIN, FMT = * )THRESH
  1269. CALL XLAENV( 1, 1 )
  1270. CALL XLAENV( 12, 11 )
  1271. CALL XLAENV( 13, 2 )
  1272. CALL XLAENV( 14, 0 )
  1273. CALL XLAENV( 15, 2 )
  1274. CALL XLAENV( 16, 2 )
  1275. TSTERR = .TRUE.
  1276. CALL DCHKEC( THRESH, TSTERR, NIN, NOUT )
  1277. GO TO 10
  1278. ELSE
  1279. WRITE( NOUT, FMT = 9992 )PATH
  1280. GO TO 10
  1281. END IF
  1282. CALL ILAVER( VERS_MAJOR, VERS_MINOR, VERS_PATCH )
  1283. WRITE( NOUT, FMT = 9972 ) VERS_MAJOR, VERS_MINOR, VERS_PATCH
  1284. WRITE( NOUT, FMT = 9984 )
  1285. *
  1286. * Read the number of values of M, P, and N.
  1287. *
  1288. READ( NIN, FMT = * )NN
  1289. IF( NN.LT.0 ) THEN
  1290. WRITE( NOUT, FMT = 9989 )' NN ', NN, 1
  1291. NN = 0
  1292. FATAL = .TRUE.
  1293. ELSE IF( NN.GT.MAXIN ) THEN
  1294. WRITE( NOUT, FMT = 9988 )' NN ', NN, MAXIN
  1295. NN = 0
  1296. FATAL = .TRUE.
  1297. END IF
  1298. *
  1299. * Read the values of M
  1300. *
  1301. IF( .NOT.( DGX .OR. DXV ) ) THEN
  1302. READ( NIN, FMT = * )( MVAL( I ), I = 1, NN )
  1303. IF( SVD ) THEN
  1304. VNAME = ' M '
  1305. ELSE
  1306. VNAME = ' N '
  1307. END IF
  1308. DO 20 I = 1, NN
  1309. IF( MVAL( I ).LT.0 ) THEN
  1310. WRITE( NOUT, FMT = 9989 )VNAME, MVAL( I ), 0
  1311. FATAL = .TRUE.
  1312. ELSE IF( MVAL( I ).GT.NMAX ) THEN
  1313. WRITE( NOUT, FMT = 9988 )VNAME, MVAL( I ), NMAX
  1314. FATAL = .TRUE.
  1315. END IF
  1316. 20 CONTINUE
  1317. WRITE( NOUT, FMT = 9983 )'M: ', ( MVAL( I ), I = 1, NN )
  1318. END IF
  1319. *
  1320. * Read the values of P
  1321. *
  1322. IF( GLM .OR. GQR .OR. GSV .OR. CSD .OR. LSE ) THEN
  1323. READ( NIN, FMT = * )( PVAL( I ), I = 1, NN )
  1324. DO 30 I = 1, NN
  1325. IF( PVAL( I ).LT.0 ) THEN
  1326. WRITE( NOUT, FMT = 9989 )' P ', PVAL( I ), 0
  1327. FATAL = .TRUE.
  1328. ELSE IF( PVAL( I ).GT.NMAX ) THEN
  1329. WRITE( NOUT, FMT = 9988 )' P ', PVAL( I ), NMAX
  1330. FATAL = .TRUE.
  1331. END IF
  1332. 30 CONTINUE
  1333. WRITE( NOUT, FMT = 9983 )'P: ', ( PVAL( I ), I = 1, NN )
  1334. END IF
  1335. *
  1336. * Read the values of N
  1337. *
  1338. IF( SVD .OR. DBB .OR. GLM .OR. GQR .OR. GSV .OR. CSD .OR.
  1339. $ LSE ) THEN
  1340. READ( NIN, FMT = * )( NVAL( I ), I = 1, NN )
  1341. DO 40 I = 1, NN
  1342. IF( NVAL( I ).LT.0 ) THEN
  1343. WRITE( NOUT, FMT = 9989 )' N ', NVAL( I ), 0
  1344. FATAL = .TRUE.
  1345. ELSE IF( NVAL( I ).GT.NMAX ) THEN
  1346. WRITE( NOUT, FMT = 9988 )' N ', NVAL( I ), NMAX
  1347. FATAL = .TRUE.
  1348. END IF
  1349. 40 CONTINUE
  1350. ELSE
  1351. DO 50 I = 1, NN
  1352. NVAL( I ) = MVAL( I )
  1353. 50 CONTINUE
  1354. END IF
  1355. IF( .NOT.( DGX .OR. DXV ) ) THEN
  1356. WRITE( NOUT, FMT = 9983 )'N: ', ( NVAL( I ), I = 1, NN )
  1357. ELSE
  1358. WRITE( NOUT, FMT = 9983 )'N: ', NN
  1359. END IF
  1360. *
  1361. * Read the number of values of K, followed by the values of K
  1362. *
  1363. IF( DSB .OR. DBB ) THEN
  1364. READ( NIN, FMT = * )NK
  1365. READ( NIN, FMT = * )( KVAL( I ), I = 1, NK )
  1366. DO 60 I = 1, NK
  1367. IF( KVAL( I ).LT.0 ) THEN
  1368. WRITE( NOUT, FMT = 9989 )' K ', KVAL( I ), 0
  1369. FATAL = .TRUE.
  1370. ELSE IF( KVAL( I ).GT.NMAX ) THEN
  1371. WRITE( NOUT, FMT = 9988 )' K ', KVAL( I ), NMAX
  1372. FATAL = .TRUE.
  1373. END IF
  1374. 60 CONTINUE
  1375. WRITE( NOUT, FMT = 9983 )'K: ', ( KVAL( I ), I = 1, NK )
  1376. END IF
  1377. *
  1378. IF( DEV .OR. DES .OR. DVX .OR. DSX ) THEN
  1379. *
  1380. * For the nonsymmetric QR driver routines, only one set of
  1381. * parameters is allowed.
  1382. *
  1383. READ( NIN, FMT = * )NBVAL( 1 ), NBMIN( 1 ), NXVAL( 1 ),
  1384. $ INMIN( 1 ), INWIN( 1 ), INIBL(1), ISHFTS(1), IACC22(1)
  1385. IF( NBVAL( 1 ).LT.1 ) THEN
  1386. WRITE( NOUT, FMT = 9989 )' NB ', NBVAL( 1 ), 1
  1387. FATAL = .TRUE.
  1388. ELSE IF( NBMIN( 1 ).LT.1 ) THEN
  1389. WRITE( NOUT, FMT = 9989 )'NBMIN ', NBMIN( 1 ), 1
  1390. FATAL = .TRUE.
  1391. ELSE IF( NXVAL( 1 ).LT.1 ) THEN
  1392. WRITE( NOUT, FMT = 9989 )' NX ', NXVAL( 1 ), 1
  1393. FATAL = .TRUE.
  1394. ELSE IF( INMIN( 1 ).LT.1 ) THEN
  1395. WRITE( NOUT, FMT = 9989 )' INMIN ', INMIN( 1 ), 1
  1396. FATAL = .TRUE.
  1397. ELSE IF( INWIN( 1 ).LT.1 ) THEN
  1398. WRITE( NOUT, FMT = 9989 )' INWIN ', INWIN( 1 ), 1
  1399. FATAL = .TRUE.
  1400. ELSE IF( INIBL( 1 ).LT.1 ) THEN
  1401. WRITE( NOUT, FMT = 9989 )' INIBL ', INIBL( 1 ), 1
  1402. FATAL = .TRUE.
  1403. ELSE IF( ISHFTS( 1 ).LT.1 ) THEN
  1404. WRITE( NOUT, FMT = 9989 )' ISHFTS ', ISHFTS( 1 ), 1
  1405. FATAL = .TRUE.
  1406. ELSE IF( IACC22( 1 ).LT.0 ) THEN
  1407. WRITE( NOUT, FMT = 9989 )' IACC22 ', IACC22( 1 ), 0
  1408. FATAL = .TRUE.
  1409. END IF
  1410. CALL XLAENV( 1, NBVAL( 1 ) )
  1411. CALL XLAENV( 2, NBMIN( 1 ) )
  1412. CALL XLAENV( 3, NXVAL( 1 ) )
  1413. CALL XLAENV(12, MAX( 11, INMIN( 1 ) ) )
  1414. CALL XLAENV(13, INWIN( 1 ) )
  1415. CALL XLAENV(14, INIBL( 1 ) )
  1416. CALL XLAENV(15, ISHFTS( 1 ) )
  1417. CALL XLAENV(16, IACC22( 1 ) )
  1418. WRITE( NOUT, FMT = 9983 )'NB: ', NBVAL( 1 )
  1419. WRITE( NOUT, FMT = 9983 )'NBMIN:', NBMIN( 1 )
  1420. WRITE( NOUT, FMT = 9983 )'NX: ', NXVAL( 1 )
  1421. WRITE( NOUT, FMT = 9983 )'INMIN: ', INMIN( 1 )
  1422. WRITE( NOUT, FMT = 9983 )'INWIN: ', INWIN( 1 )
  1423. WRITE( NOUT, FMT = 9983 )'INIBL: ', INIBL( 1 )
  1424. WRITE( NOUT, FMT = 9983 )'ISHFTS: ', ISHFTS( 1 )
  1425. WRITE( NOUT, FMT = 9983 )'IACC22: ', IACC22( 1 )
  1426. *
  1427. ELSEIF( DGS .OR. DGX .OR. DGV .OR. DXV ) THEN
  1428. *
  1429. * For the nonsymmetric generalized driver routines, only one set
  1430. * of parameters is allowed.
  1431. *
  1432. READ( NIN, FMT = * )NBVAL( 1 ), NBMIN( 1 ), NXVAL( 1 ),
  1433. $ NSVAL( 1 ), MXBVAL( 1 )
  1434. IF( NBVAL( 1 ).LT.1 ) THEN
  1435. WRITE( NOUT, FMT = 9989 )' NB ', NBVAL( 1 ), 1
  1436. FATAL = .TRUE.
  1437. ELSE IF( NBMIN( 1 ).LT.1 ) THEN
  1438. WRITE( NOUT, FMT = 9989 )'NBMIN ', NBMIN( 1 ), 1
  1439. FATAL = .TRUE.
  1440. ELSE IF( NXVAL( 1 ).LT.1 ) THEN
  1441. WRITE( NOUT, FMT = 9989 )' NX ', NXVAL( 1 ), 1
  1442. FATAL = .TRUE.
  1443. ELSE IF( NSVAL( 1 ).LT.2 ) THEN
  1444. WRITE( NOUT, FMT = 9989 )' NS ', NSVAL( 1 ), 2
  1445. FATAL = .TRUE.
  1446. ELSE IF( MXBVAL( 1 ).LT.1 ) THEN
  1447. WRITE( NOUT, FMT = 9989 )' MAXB ', MXBVAL( 1 ), 1
  1448. FATAL = .TRUE.
  1449. END IF
  1450. CALL XLAENV( 1, NBVAL( 1 ) )
  1451. CALL XLAENV( 2, NBMIN( 1 ) )
  1452. CALL XLAENV( 3, NXVAL( 1 ) )
  1453. CALL XLAENV( 4, NSVAL( 1 ) )
  1454. CALL XLAENV( 8, MXBVAL( 1 ) )
  1455. WRITE( NOUT, FMT = 9983 )'NB: ', NBVAL( 1 )
  1456. WRITE( NOUT, FMT = 9983 )'NBMIN:', NBMIN( 1 )
  1457. WRITE( NOUT, FMT = 9983 )'NX: ', NXVAL( 1 )
  1458. WRITE( NOUT, FMT = 9983 )'NS: ', NSVAL( 1 )
  1459. WRITE( NOUT, FMT = 9983 )'MAXB: ', MXBVAL( 1 )
  1460. *
  1461. ELSE IF( .NOT.DSB .AND. .NOT.GLM .AND. .NOT.GQR .AND. .NOT.
  1462. $ GSV .AND. .NOT.CSD .AND. .NOT.LSE ) THEN
  1463. *
  1464. * For the other paths, the number of parameters can be varied
  1465. * from the input file. Read the number of parameter values.
  1466. *
  1467. READ( NIN, FMT = * )NPARMS
  1468. IF( NPARMS.LT.1 ) THEN
  1469. WRITE( NOUT, FMT = 9989 )'NPARMS', NPARMS, 1
  1470. NPARMS = 0
  1471. FATAL = .TRUE.
  1472. ELSE IF( NPARMS.GT.MAXIN ) THEN
  1473. WRITE( NOUT, FMT = 9988 )'NPARMS', NPARMS, MAXIN
  1474. NPARMS = 0
  1475. FATAL = .TRUE.
  1476. END IF
  1477. *
  1478. * Read the values of NB
  1479. *
  1480. IF( .NOT.DBB ) THEN
  1481. READ( NIN, FMT = * )( NBVAL( I ), I = 1, NPARMS )
  1482. DO 70 I = 1, NPARMS
  1483. IF( NBVAL( I ).LT.0 ) THEN
  1484. WRITE( NOUT, FMT = 9989 )' NB ', NBVAL( I ), 0
  1485. FATAL = .TRUE.
  1486. ELSE IF( NBVAL( I ).GT.NMAX ) THEN
  1487. WRITE( NOUT, FMT = 9988 )' NB ', NBVAL( I ), NMAX
  1488. FATAL = .TRUE.
  1489. END IF
  1490. 70 CONTINUE
  1491. WRITE( NOUT, FMT = 9983 )'NB: ',
  1492. $ ( NBVAL( I ), I = 1, NPARMS )
  1493. END IF
  1494. *
  1495. * Read the values of NBMIN
  1496. *
  1497. IF( NEP .OR. SEP .OR. SVD .OR. DGG ) THEN
  1498. READ( NIN, FMT = * )( NBMIN( I ), I = 1, NPARMS )
  1499. DO 80 I = 1, NPARMS
  1500. IF( NBMIN( I ).LT.0 ) THEN
  1501. WRITE( NOUT, FMT = 9989 )'NBMIN ', NBMIN( I ), 0
  1502. FATAL = .TRUE.
  1503. ELSE IF( NBMIN( I ).GT.NMAX ) THEN
  1504. WRITE( NOUT, FMT = 9988 )'NBMIN ', NBMIN( I ), NMAX
  1505. FATAL = .TRUE.
  1506. END IF
  1507. 80 CONTINUE
  1508. WRITE( NOUT, FMT = 9983 )'NBMIN:',
  1509. $ ( NBMIN( I ), I = 1, NPARMS )
  1510. ELSE
  1511. DO 90 I = 1, NPARMS
  1512. NBMIN( I ) = 1
  1513. 90 CONTINUE
  1514. END IF
  1515. *
  1516. * Read the values of NX
  1517. *
  1518. IF( NEP .OR. SEP .OR. SVD ) THEN
  1519. READ( NIN, FMT = * )( NXVAL( I ), I = 1, NPARMS )
  1520. DO 100 I = 1, NPARMS
  1521. IF( NXVAL( I ).LT.0 ) THEN
  1522. WRITE( NOUT, FMT = 9989 )' NX ', NXVAL( I ), 0
  1523. FATAL = .TRUE.
  1524. ELSE IF( NXVAL( I ).GT.NMAX ) THEN
  1525. WRITE( NOUT, FMT = 9988 )' NX ', NXVAL( I ), NMAX
  1526. FATAL = .TRUE.
  1527. END IF
  1528. 100 CONTINUE
  1529. WRITE( NOUT, FMT = 9983 )'NX: ',
  1530. $ ( NXVAL( I ), I = 1, NPARMS )
  1531. ELSE
  1532. DO 110 I = 1, NPARMS
  1533. NXVAL( I ) = 1
  1534. 110 CONTINUE
  1535. END IF
  1536. *
  1537. * Read the values of NSHIFT (if DGG) or NRHS (if SVD
  1538. * or DBB).
  1539. *
  1540. IF( SVD .OR. DBB .OR. DGG ) THEN
  1541. READ( NIN, FMT = * )( NSVAL( I ), I = 1, NPARMS )
  1542. DO 120 I = 1, NPARMS
  1543. IF( NSVAL( I ).LT.0 ) THEN
  1544. WRITE( NOUT, FMT = 9989 )' NS ', NSVAL( I ), 0
  1545. FATAL = .TRUE.
  1546. ELSE IF( NSVAL( I ).GT.NMAX ) THEN
  1547. WRITE( NOUT, FMT = 9988 )' NS ', NSVAL( I ), NMAX
  1548. FATAL = .TRUE.
  1549. END IF
  1550. 120 CONTINUE
  1551. WRITE( NOUT, FMT = 9983 )'NS: ',
  1552. $ ( NSVAL( I ), I = 1, NPARMS )
  1553. ELSE
  1554. DO 130 I = 1, NPARMS
  1555. NSVAL( I ) = 1
  1556. 130 CONTINUE
  1557. END IF
  1558. *
  1559. * Read the values for MAXB.
  1560. *
  1561. IF( DGG ) THEN
  1562. READ( NIN, FMT = * )( MXBVAL( I ), I = 1, NPARMS )
  1563. DO 140 I = 1, NPARMS
  1564. IF( MXBVAL( I ).LT.0 ) THEN
  1565. WRITE( NOUT, FMT = 9989 )' MAXB ', MXBVAL( I ), 0
  1566. FATAL = .TRUE.
  1567. ELSE IF( MXBVAL( I ).GT.NMAX ) THEN
  1568. WRITE( NOUT, FMT = 9988 )' MAXB ', MXBVAL( I ), NMAX
  1569. FATAL = .TRUE.
  1570. END IF
  1571. 140 CONTINUE
  1572. WRITE( NOUT, FMT = 9983 )'MAXB: ',
  1573. $ ( MXBVAL( I ), I = 1, NPARMS )
  1574. ELSE
  1575. DO 150 I = 1, NPARMS
  1576. MXBVAL( I ) = 1
  1577. 150 CONTINUE
  1578. END IF
  1579. *
  1580. * Read the values for INMIN.
  1581. *
  1582. IF( NEP ) THEN
  1583. READ( NIN, FMT = * )( INMIN( I ), I = 1, NPARMS )
  1584. DO 540 I = 1, NPARMS
  1585. IF( INMIN( I ).LT.0 ) THEN
  1586. WRITE( NOUT, FMT = 9989 )' INMIN ', INMIN( I ), 0
  1587. FATAL = .TRUE.
  1588. END IF
  1589. 540 CONTINUE
  1590. WRITE( NOUT, FMT = 9983 )'INMIN: ',
  1591. $ ( INMIN( I ), I = 1, NPARMS )
  1592. ELSE
  1593. DO 550 I = 1, NPARMS
  1594. INMIN( I ) = 1
  1595. 550 CONTINUE
  1596. END IF
  1597. *
  1598. * Read the values for INWIN.
  1599. *
  1600. IF( NEP ) THEN
  1601. READ( NIN, FMT = * )( INWIN( I ), I = 1, NPARMS )
  1602. DO 560 I = 1, NPARMS
  1603. IF( INWIN( I ).LT.0 ) THEN
  1604. WRITE( NOUT, FMT = 9989 )' INWIN ', INWIN( I ), 0
  1605. FATAL = .TRUE.
  1606. END IF
  1607. 560 CONTINUE
  1608. WRITE( NOUT, FMT = 9983 )'INWIN: ',
  1609. $ ( INWIN( I ), I = 1, NPARMS )
  1610. ELSE
  1611. DO 570 I = 1, NPARMS
  1612. INWIN( I ) = 1
  1613. 570 CONTINUE
  1614. END IF
  1615. *
  1616. * Read the values for INIBL.
  1617. *
  1618. IF( NEP ) THEN
  1619. READ( NIN, FMT = * )( INIBL( I ), I = 1, NPARMS )
  1620. DO 580 I = 1, NPARMS
  1621. IF( INIBL( I ).LT.0 ) THEN
  1622. WRITE( NOUT, FMT = 9989 )' INIBL ', INIBL( I ), 0
  1623. FATAL = .TRUE.
  1624. END IF
  1625. 580 CONTINUE
  1626. WRITE( NOUT, FMT = 9983 )'INIBL: ',
  1627. $ ( INIBL( I ), I = 1, NPARMS )
  1628. ELSE
  1629. DO 590 I = 1, NPARMS
  1630. INIBL( I ) = 1
  1631. 590 CONTINUE
  1632. END IF
  1633. *
  1634. * Read the values for ISHFTS.
  1635. *
  1636. IF( NEP ) THEN
  1637. READ( NIN, FMT = * )( ISHFTS( I ), I = 1, NPARMS )
  1638. DO 600 I = 1, NPARMS
  1639. IF( ISHFTS( I ).LT.0 ) THEN
  1640. WRITE( NOUT, FMT = 9989 )' ISHFTS ', ISHFTS( I ), 0
  1641. FATAL = .TRUE.
  1642. END IF
  1643. 600 CONTINUE
  1644. WRITE( NOUT, FMT = 9983 )'ISHFTS: ',
  1645. $ ( ISHFTS( I ), I = 1, NPARMS )
  1646. ELSE
  1647. DO 610 I = 1, NPARMS
  1648. ISHFTS( I ) = 1
  1649. 610 CONTINUE
  1650. END IF
  1651. *
  1652. * Read the values for IACC22.
  1653. *
  1654. IF( NEP .OR. DGG ) THEN
  1655. READ( NIN, FMT = * )( IACC22( I ), I = 1, NPARMS )
  1656. DO 620 I = 1, NPARMS
  1657. IF( IACC22( I ).LT.0 ) THEN
  1658. WRITE( NOUT, FMT = 9989 )' IACC22 ', IACC22( I ), 0
  1659. FATAL = .TRUE.
  1660. END IF
  1661. 620 CONTINUE
  1662. WRITE( NOUT, FMT = 9983 )'IACC22: ',
  1663. $ ( IACC22( I ), I = 1, NPARMS )
  1664. ELSE
  1665. DO 630 I = 1, NPARMS
  1666. IACC22( I ) = 1
  1667. 630 CONTINUE
  1668. END IF
  1669. *
  1670. * Read the values for NBCOL.
  1671. *
  1672. IF( DGG ) THEN
  1673. READ( NIN, FMT = * )( NBCOL( I ), I = 1, NPARMS )
  1674. DO 160 I = 1, NPARMS
  1675. IF( NBCOL( I ).LT.0 ) THEN
  1676. WRITE( NOUT, FMT = 9989 )'NBCOL ', NBCOL( I ), 0
  1677. FATAL = .TRUE.
  1678. ELSE IF( NBCOL( I ).GT.NMAX ) THEN
  1679. WRITE( NOUT, FMT = 9988 )'NBCOL ', NBCOL( I ), NMAX
  1680. FATAL = .TRUE.
  1681. END IF
  1682. 160 CONTINUE
  1683. WRITE( NOUT, FMT = 9983 )'NBCOL:',
  1684. $ ( NBCOL( I ), I = 1, NPARMS )
  1685. ELSE
  1686. DO 170 I = 1, NPARMS
  1687. NBCOL( I ) = 1
  1688. 170 CONTINUE
  1689. END IF
  1690. END IF
  1691. *
  1692. * Calculate and print the machine dependent constants.
  1693. *
  1694. WRITE( NOUT, FMT = * )
  1695. EPS = DLAMCH( 'Underflow threshold' )
  1696. WRITE( NOUT, FMT = 9981 )'underflow', EPS
  1697. EPS = DLAMCH( 'Overflow threshold' )
  1698. WRITE( NOUT, FMT = 9981 )'overflow ', EPS
  1699. EPS = DLAMCH( 'Epsilon' )
  1700. WRITE( NOUT, FMT = 9981 )'precision', EPS
  1701. *
  1702. * Read the threshold value for the test ratios.
  1703. *
  1704. READ( NIN, FMT = * )THRESH
  1705. WRITE( NOUT, FMT = 9982 )THRESH
  1706. IF( SEP .OR. SVD .OR. DGG ) THEN
  1707. *
  1708. * Read the flag that indicates whether to test LAPACK routines.
  1709. *
  1710. READ( NIN, FMT = * )TSTCHK
  1711. *
  1712. * Read the flag that indicates whether to test driver routines.
  1713. *
  1714. READ( NIN, FMT = * )TSTDRV
  1715. END IF
  1716. *
  1717. * Read the flag that indicates whether to test the error exits.
  1718. *
  1719. READ( NIN, FMT = * )TSTERR
  1720. *
  1721. * Read the code describing how to set the random number seed.
  1722. *
  1723. READ( NIN, FMT = * )NEWSD
  1724. *
  1725. * If NEWSD = 2, read another line with 4 integers for the seed.
  1726. *
  1727. IF( NEWSD.EQ.2 )
  1728. $ READ( NIN, FMT = * )( IOLDSD( I ), I = 1, 4 )
  1729. *
  1730. DO 180 I = 1, 4
  1731. ISEED( I ) = IOLDSD( I )
  1732. 180 CONTINUE
  1733. *
  1734. IF( FATAL ) THEN
  1735. WRITE( NOUT, FMT = 9999 )
  1736. STOP
  1737. END IF
  1738. *
  1739. * Read the input lines indicating the test path and its parameters.
  1740. * The first three characters indicate the test path, and the number
  1741. * of test matrix types must be the first nonblank item in columns
  1742. * 4-80.
  1743. *
  1744. 190 CONTINUE
  1745. *
  1746. IF( .NOT.( DGX .OR. DXV ) ) THEN
  1747. *
  1748. 200 CONTINUE
  1749. READ( NIN, FMT = '(A80)', END = 380 )LINE
  1750. C3 = LINE( 1: 3 )
  1751. LENP = LEN( LINE )
  1752. I = 3
  1753. ITMP = 0
  1754. I1 = 0
  1755. 210 CONTINUE
  1756. I = I + 1
  1757. IF( I.GT.LENP ) THEN
  1758. IF( I1.GT.0 ) THEN
  1759. GO TO 240
  1760. ELSE
  1761. NTYPES = MAXT
  1762. GO TO 240
  1763. END IF
  1764. END IF
  1765. IF( LINE( I: I ).NE.' ' .AND. LINE( I: I ).NE.',' ) THEN
  1766. I1 = I
  1767. C1 = LINE( I1: I1 )
  1768. *
  1769. * Check that a valid integer was read
  1770. *
  1771. DO 220 K = 1, 10
  1772. IF( C1.EQ.INTSTR( K: K ) ) THEN
  1773. IC = K - 1
  1774. GO TO 230
  1775. END IF
  1776. 220 CONTINUE
  1777. WRITE( NOUT, FMT = 9991 )I, LINE
  1778. GO TO 200
  1779. 230 CONTINUE
  1780. ITMP = 10*ITMP + IC
  1781. GO TO 210
  1782. ELSE IF( I1.GT.0 ) THEN
  1783. GO TO 240
  1784. ELSE
  1785. GO TO 210
  1786. END IF
  1787. 240 CONTINUE
  1788. NTYPES = ITMP
  1789. *
  1790. * Skip the tests if NTYPES is <= 0.
  1791. *
  1792. IF( .NOT.( DEV .OR. DES .OR. DVX .OR. DSX .OR. DGV .OR.
  1793. $ DGS ) .AND. NTYPES.LE.0 ) THEN
  1794. WRITE( NOUT, FMT = 9990 )C3
  1795. GO TO 200
  1796. END IF
  1797. *
  1798. ELSE
  1799. IF( DXV )
  1800. $ C3 = 'DXV'
  1801. IF( DGX )
  1802. $ C3 = 'DGX'
  1803. END IF
  1804. *
  1805. * Reset the random number seed.
  1806. *
  1807. IF( NEWSD.EQ.0 ) THEN
  1808. DO 250 K = 1, 4
  1809. ISEED( K ) = IOLDSD( K )
  1810. 250 CONTINUE
  1811. END IF
  1812. *
  1813. IF( LSAMEN( 3, C3, 'DHS' ) .OR. LSAMEN( 3, C3, 'NEP' ) ) THEN
  1814. *
  1815. * -------------------------------------
  1816. * NEP: Nonsymmetric Eigenvalue Problem
  1817. * -------------------------------------
  1818. * Vary the parameters
  1819. * NB = block size
  1820. * NBMIN = minimum block size
  1821. * NX = crossover point
  1822. * NS = number of shifts
  1823. * MAXB = minimum submatrix size
  1824. *
  1825. MAXTYP = 21
  1826. NTYPES = MIN( MAXTYP, NTYPES )
  1827. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  1828. CALL XLAENV( 1, 1 )
  1829. IF( TSTERR )
  1830. $ CALL DERRHS( 'DHSEQR', NOUT )
  1831. DO 270 I = 1, NPARMS
  1832. CALL XLAENV( 1, NBVAL( I ) )
  1833. CALL XLAENV( 2, NBMIN( I ) )
  1834. CALL XLAENV( 3, NXVAL( I ) )
  1835. CALL XLAENV(12, MAX( 11, INMIN( I ) ) )
  1836. CALL XLAENV(13, INWIN( I ) )
  1837. CALL XLAENV(14, INIBL( I ) )
  1838. CALL XLAENV(15, ISHFTS( I ) )
  1839. CALL XLAENV(16, IACC22( I ) )
  1840. *
  1841. IF( NEWSD.EQ.0 ) THEN
  1842. DO 260 K = 1, 4
  1843. ISEED( K ) = IOLDSD( K )
  1844. 260 CONTINUE
  1845. END IF
  1846. WRITE( NOUT, FMT = 9961 )C3, NBVAL( I ), NBMIN( I ),
  1847. $ NXVAL( I ), MAX( 11, INMIN(I)),
  1848. $ INWIN( I ), INIBL( I ), ISHFTS( I ), IACC22( I )
  1849. CALL DCHKHS( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH, NOUT,
  1850. $ A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  1851. $ A( 1, 4 ), A( 1, 5 ), NMAX, A( 1, 6 ),
  1852. $ A( 1, 7 ), D( 1, 1 ), D( 1, 2 ), D( 1, 3 ),
  1853. $ D( 1, 4 ), D( 1, 5 ), D( 1, 6 ), A( 1, 8 ),
  1854. $ A( 1, 9 ), A( 1, 10 ), A( 1, 11 ), A( 1, 12 ),
  1855. $ D( 1, 7 ), WORK, LWORK, IWORK, LOGWRK, RESULT,
  1856. $ INFO )
  1857. IF( INFO.NE.0 )
  1858. $ WRITE( NOUT, FMT = 9980 )'DCHKHS', INFO
  1859. 270 CONTINUE
  1860. *
  1861. ELSE IF( LSAMEN( 3, C3, 'DST' ) .OR. LSAMEN( 3, C3, 'SEP' )
  1862. $ .OR. LSAMEN( 3, C3, 'SE2' ) ) THEN
  1863. *
  1864. * ----------------------------------
  1865. * SEP: Symmetric Eigenvalue Problem
  1866. * ----------------------------------
  1867. * Vary the parameters
  1868. * NB = block size
  1869. * NBMIN = minimum block size
  1870. * NX = crossover point
  1871. *
  1872. MAXTYP = 21
  1873. NTYPES = MIN( MAXTYP, NTYPES )
  1874. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  1875. CALL XLAENV( 1, 1 )
  1876. CALL XLAENV( 9, 25 )
  1877. IF( TSTERR ) THEN
  1878. #if defined(_OPENMP)
  1879. N_THREADS = OMP_GET_MAX_THREADS()
  1880. CALL OMP_SET_NUM_THREADS(1)
  1881. #endif
  1882. CALL DERRST( 'DST', NOUT )
  1883. #if defined(_OPENMP)
  1884. CALL OMP_SET_NUM_THREADS(N_THREADS)
  1885. #endif
  1886. END IF
  1887. DO 290 I = 1, NPARMS
  1888. CALL XLAENV( 1, NBVAL( I ) )
  1889. CALL XLAENV( 2, NBMIN( I ) )
  1890. CALL XLAENV( 3, NXVAL( I ) )
  1891. *
  1892. IF( NEWSD.EQ.0 ) THEN
  1893. DO 280 K = 1, 4
  1894. ISEED( K ) = IOLDSD( K )
  1895. 280 CONTINUE
  1896. END IF
  1897. WRITE( NOUT, FMT = 9997 )C3, NBVAL( I ), NBMIN( I ),
  1898. $ NXVAL( I )
  1899. IF( TSTCHK ) THEN
  1900. IF( LSAMEN( 3, C3, 'SE2' ) ) THEN
  1901. CALL DCHKST2STG( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH,
  1902. $ NOUT, A( 1, 1 ), NMAX, A( 1, 2 ), D( 1, 1 ),
  1903. $ D( 1, 2 ), D( 1, 3 ), D( 1, 4 ), D( 1, 5 ),
  1904. $ D( 1, 6 ), D( 1, 7 ), D( 1, 8 ), D( 1, 9 ),
  1905. $ D( 1, 10 ), D( 1, 11 ), A( 1, 3 ), NMAX,
  1906. $ A( 1, 4 ), A( 1, 5 ), D( 1, 12 ), A( 1, 6 ),
  1907. $ WORK, LWORK, IWORK, LIWORK, RESULT, INFO )
  1908. ELSE
  1909. CALL DCHKST( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH,
  1910. $ NOUT, A( 1, 1 ), NMAX, A( 1, 2 ), D( 1, 1 ),
  1911. $ D( 1, 2 ), D( 1, 3 ), D( 1, 4 ), D( 1, 5 ),
  1912. $ D( 1, 6 ), D( 1, 7 ), D( 1, 8 ), D( 1, 9 ),
  1913. $ D( 1, 10 ), D( 1, 11 ), A( 1, 3 ), NMAX,
  1914. $ A( 1, 4 ), A( 1, 5 ), D( 1, 12 ), A( 1, 6 ),
  1915. $ WORK, LWORK, IWORK, LIWORK, RESULT, INFO )
  1916. ENDIF
  1917. IF( INFO.NE.0 )
  1918. $ WRITE( NOUT, FMT = 9980 )'DCHKST', INFO
  1919. END IF
  1920. IF( TSTDRV ) THEN
  1921. IF( LSAMEN( 3, C3, 'SE2' ) ) THEN
  1922. CALL DDRVST2STG( NN, NVAL, 18, DOTYPE, ISEED, THRESH,
  1923. $ NOUT, A( 1, 1 ), NMAX, D( 1, 3 ), D( 1, 4 ),
  1924. $ D( 1, 5 ), D( 1, 6 ), D( 1, 8 ), D( 1, 9 ),
  1925. $ D( 1, 10 ), D( 1, 11 ), A( 1, 2 ), NMAX,
  1926. $ A( 1, 3 ), D( 1, 12 ), A( 1, 4 ), WORK,
  1927. $ LWORK, IWORK, LIWORK, RESULT, INFO )
  1928. ELSE
  1929. CALL DDRVST( NN, NVAL, 18, DOTYPE, ISEED, THRESH, NOUT,
  1930. $ A( 1, 1 ), NMAX, D( 1, 3 ), D( 1, 4 ),
  1931. $ D( 1, 5 ), D( 1, 6 ), D( 1, 8 ), D( 1, 9 ),
  1932. $ D( 1, 10 ), D( 1, 11 ), A( 1, 2 ), NMAX,
  1933. $ A( 1, 3 ), D( 1, 12 ), A( 1, 4 ), WORK,
  1934. $ LWORK, IWORK, LIWORK, RESULT, INFO )
  1935. ENDIF
  1936. IF( INFO.NE.0 )
  1937. $ WRITE( NOUT, FMT = 9980 )'DDRVST', INFO
  1938. END IF
  1939. 290 CONTINUE
  1940. *
  1941. ELSE IF( LSAMEN( 3, C3, 'DSG' ) ) THEN
  1942. *
  1943. * ----------------------------------------------
  1944. * DSG: Symmetric Generalized Eigenvalue Problem
  1945. * ----------------------------------------------
  1946. * Vary the parameters
  1947. * NB = block size
  1948. * NBMIN = minimum block size
  1949. * NX = crossover point
  1950. *
  1951. MAXTYP = 21
  1952. NTYPES = MIN( MAXTYP, NTYPES )
  1953. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  1954. CALL XLAENV( 9, 25 )
  1955. DO 310 I = 1, NPARMS
  1956. CALL XLAENV( 1, NBVAL( I ) )
  1957. CALL XLAENV( 2, NBMIN( I ) )
  1958. CALL XLAENV( 3, NXVAL( I ) )
  1959. *
  1960. IF( NEWSD.EQ.0 ) THEN
  1961. DO 300 K = 1, 4
  1962. ISEED( K ) = IOLDSD( K )
  1963. 300 CONTINUE
  1964. END IF
  1965. WRITE( NOUT, FMT = 9997 )C3, NBVAL( I ), NBMIN( I ),
  1966. $ NXVAL( I )
  1967. IF( TSTCHK ) THEN
  1968. * CALL DDRVSG( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH,
  1969. * $ NOUT, A( 1, 1 ), NMAX, A( 1, 2 ), NMAX,
  1970. * $ D( 1, 3 ), A( 1, 3 ), NMAX, A( 1, 4 ),
  1971. * $ A( 1, 5 ), A( 1, 6 ), A( 1, 7 ), WORK,
  1972. * $ LWORK, IWORK, LIWORK, RESULT, INFO )
  1973. CALL DDRVSG2STG( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH,
  1974. $ NOUT, A( 1, 1 ), NMAX, A( 1, 2 ), NMAX,
  1975. $ D( 1, 3 ), D( 1, 3 ), A( 1, 3 ), NMAX,
  1976. $ A( 1, 4 ), A( 1, 5 ), A( 1, 6 ),
  1977. $ A( 1, 7 ), WORK, LWORK, IWORK, LIWORK,
  1978. $ RESULT, INFO )
  1979. IF( INFO.NE.0 )
  1980. $ WRITE( NOUT, FMT = 9980 )'DDRVSG', INFO
  1981. END IF
  1982. 310 CONTINUE
  1983. *
  1984. ELSE IF( LSAMEN( 3, C3, 'DBD' ) .OR. LSAMEN( 3, C3, 'SVD' ) ) THEN
  1985. *
  1986. * ----------------------------------
  1987. * SVD: Singular Value Decomposition
  1988. * ----------------------------------
  1989. * Vary the parameters
  1990. * NB = block size
  1991. * NBMIN = minimum block size
  1992. * NX = crossover point
  1993. * NRHS = number of right hand sides
  1994. *
  1995. MAXTYP = 16
  1996. NTYPES = MIN( MAXTYP, NTYPES )
  1997. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  1998. CALL XLAENV( 1, 1 )
  1999. CALL XLAENV( 9, 25 )
  2000. *
  2001. * Test the error exits
  2002. *
  2003. IF( TSTERR .AND. TSTCHK )
  2004. $ CALL DERRBD( 'DBD', NOUT )
  2005. IF( TSTERR .AND. TSTDRV )
  2006. $ CALL DERRED( 'DBD', NOUT )
  2007. *
  2008. DO 330 I = 1, NPARMS
  2009. NRHS = NSVAL( I )
  2010. CALL XLAENV( 1, NBVAL( I ) )
  2011. CALL XLAENV( 2, NBMIN( I ) )
  2012. CALL XLAENV( 3, NXVAL( I ) )
  2013. IF( NEWSD.EQ.0 ) THEN
  2014. DO 320 K = 1, 4
  2015. ISEED( K ) = IOLDSD( K )
  2016. 320 CONTINUE
  2017. END IF
  2018. WRITE( NOUT, FMT = 9995 )C3, NBVAL( I ), NBMIN( I ),
  2019. $ NXVAL( I ), NRHS
  2020. IF( TSTCHK ) THEN
  2021. CALL DCHKBD( NN, MVAL, NVAL, MAXTYP, DOTYPE, NRHS, ISEED,
  2022. $ THRESH, A( 1, 1 ), NMAX, D( 1, 1 ),
  2023. $ D( 1, 2 ), D( 1, 3 ), D( 1, 4 ), A( 1, 2 ),
  2024. $ NMAX, A( 1, 3 ), A( 1, 4 ), A( 1, 5 ), NMAX,
  2025. $ A( 1, 6 ), NMAX, A( 1, 7 ), A( 1, 8 ), WORK,
  2026. $ LWORK, IWORK, NOUT, INFO )
  2027. IF( INFO.NE.0 )
  2028. $ WRITE( NOUT, FMT = 9980 )'DCHKBD', INFO
  2029. END IF
  2030. IF( TSTDRV )
  2031. $ CALL DDRVBD( NN, MVAL, NVAL, MAXTYP, DOTYPE, ISEED,
  2032. $ THRESH, A( 1, 1 ), NMAX, A( 1, 2 ), NMAX,
  2033. $ A( 1, 3 ), NMAX, A( 1, 4 ), A( 1, 5 ),
  2034. $ A( 1, 6 ), D( 1, 1 ), D( 1, 2 ), D( 1, 3 ),
  2035. $ WORK, LWORK, IWORK, NOUT, INFO )
  2036. 330 CONTINUE
  2037. *
  2038. ELSE IF( LSAMEN( 3, C3, 'DEV' ) ) THEN
  2039. *
  2040. * --------------------------------------------
  2041. * DEV: Nonsymmetric Eigenvalue Problem Driver
  2042. * DGEEV (eigenvalues and eigenvectors)
  2043. * --------------------------------------------
  2044. *
  2045. MAXTYP = 21
  2046. NTYPES = MIN( MAXTYP, NTYPES )
  2047. IF( NTYPES.LE.0 ) THEN
  2048. WRITE( NOUT, FMT = 9990 )C3
  2049. ELSE
  2050. IF( TSTERR )
  2051. $ CALL DERRED( C3, NOUT )
  2052. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2053. CALL DDRVEV( NN, NVAL, NTYPES, DOTYPE, ISEED, THRESH, NOUT,
  2054. $ A( 1, 1 ), NMAX, A( 1, 2 ), D( 1, 1 ),
  2055. $ D( 1, 2 ), D( 1, 3 ), D( 1, 4 ), A( 1, 3 ),
  2056. $ NMAX, A( 1, 4 ), NMAX, A( 1, 5 ), NMAX, RESULT,
  2057. $ WORK, LWORK, IWORK, INFO )
  2058. IF( INFO.NE.0 )
  2059. $ WRITE( NOUT, FMT = 9980 )'DGEEV', INFO
  2060. END IF
  2061. WRITE( NOUT, FMT = 9973 )
  2062. GO TO 10
  2063. *
  2064. ELSE IF( LSAMEN( 3, C3, 'DES' ) ) THEN
  2065. *
  2066. * --------------------------------------------
  2067. * DES: Nonsymmetric Eigenvalue Problem Driver
  2068. * DGEES (Schur form)
  2069. * --------------------------------------------
  2070. *
  2071. MAXTYP = 21
  2072. NTYPES = MIN( MAXTYP, NTYPES )
  2073. IF( NTYPES.LE.0 ) THEN
  2074. WRITE( NOUT, FMT = 9990 )C3
  2075. ELSE
  2076. IF( TSTERR )
  2077. $ CALL DERRED( C3, NOUT )
  2078. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2079. CALL DDRVES( NN, NVAL, NTYPES, DOTYPE, ISEED, THRESH, NOUT,
  2080. $ A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  2081. $ D( 1, 1 ), D( 1, 2 ), D( 1, 3 ), D( 1, 4 ),
  2082. $ A( 1, 4 ), NMAX, RESULT, WORK, LWORK, IWORK,
  2083. $ LOGWRK, INFO )
  2084. IF( INFO.NE.0 )
  2085. $ WRITE( NOUT, FMT = 9980 )'DGEES', INFO
  2086. END IF
  2087. WRITE( NOUT, FMT = 9973 )
  2088. GO TO 10
  2089. *
  2090. ELSE IF( LSAMEN( 3, C3, 'DVX' ) ) THEN
  2091. *
  2092. * --------------------------------------------------------------
  2093. * DVX: Nonsymmetric Eigenvalue Problem Expert Driver
  2094. * DGEEVX (eigenvalues, eigenvectors and condition numbers)
  2095. * --------------------------------------------------------------
  2096. *
  2097. MAXTYP = 21
  2098. NTYPES = MIN( MAXTYP, NTYPES )
  2099. IF( NTYPES.LT.0 ) THEN
  2100. WRITE( NOUT, FMT = 9990 )C3
  2101. ELSE
  2102. IF( TSTERR )
  2103. $ CALL DERRED( C3, NOUT )
  2104. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2105. CALL DDRVVX( NN, NVAL, NTYPES, DOTYPE, ISEED, THRESH, NIN,
  2106. $ NOUT, A( 1, 1 ), NMAX, A( 1, 2 ), D( 1, 1 ),
  2107. $ D( 1, 2 ), D( 1, 3 ), D( 1, 4 ), A( 1, 3 ),
  2108. $ NMAX, A( 1, 4 ), NMAX, A( 1, 5 ), NMAX,
  2109. $ D( 1, 5 ), D( 1, 6 ), D( 1, 7 ), D( 1, 8 ),
  2110. $ D( 1, 9 ), D( 1, 10 ), D( 1, 11 ), D( 1, 12 ),
  2111. $ RESULT, WORK, LWORK, IWORK, INFO )
  2112. IF( INFO.NE.0 )
  2113. $ WRITE( NOUT, FMT = 9980 )'DGEEVX', INFO
  2114. END IF
  2115. WRITE( NOUT, FMT = 9973 )
  2116. GO TO 10
  2117. *
  2118. ELSE IF( LSAMEN( 3, C3, 'DSX' ) ) THEN
  2119. *
  2120. * ---------------------------------------------------
  2121. * DSX: Nonsymmetric Eigenvalue Problem Expert Driver
  2122. * DGEESX (Schur form and condition numbers)
  2123. * ---------------------------------------------------
  2124. *
  2125. MAXTYP = 21
  2126. NTYPES = MIN( MAXTYP, NTYPES )
  2127. IF( NTYPES.LT.0 ) THEN
  2128. WRITE( NOUT, FMT = 9990 )C3
  2129. ELSE
  2130. IF( TSTERR )
  2131. $ CALL DERRED( C3, NOUT )
  2132. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2133. CALL DDRVSX( NN, NVAL, NTYPES, DOTYPE, ISEED, THRESH, NIN,
  2134. $ NOUT, A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  2135. $ D( 1, 1 ), D( 1, 2 ), D( 1, 3 ), D( 1, 4 ),
  2136. $ D( 1, 5 ), D( 1, 6 ), A( 1, 4 ), NMAX,
  2137. $ A( 1, 5 ), RESULT, WORK, LWORK, IWORK, LOGWRK,
  2138. $ INFO )
  2139. IF( INFO.NE.0 )
  2140. $ WRITE( NOUT, FMT = 9980 )'DGEESX', INFO
  2141. END IF
  2142. WRITE( NOUT, FMT = 9973 )
  2143. GO TO 10
  2144. *
  2145. ELSE IF( LSAMEN( 3, C3, 'DGG' ) ) THEN
  2146. *
  2147. * -------------------------------------------------
  2148. * DGG: Generalized Nonsymmetric Eigenvalue Problem
  2149. * -------------------------------------------------
  2150. * Vary the parameters
  2151. * NB = block size
  2152. * NBMIN = minimum block size
  2153. * NS = number of shifts
  2154. * MAXB = minimum submatrix size
  2155. * IACC22: structured matrix multiply
  2156. * NBCOL = minimum column dimension for blocks
  2157. *
  2158. MAXTYP = 26
  2159. NTYPES = MIN( MAXTYP, NTYPES )
  2160. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2161. CALL XLAENV(1,1)
  2162. IF( TSTCHK .AND. TSTERR )
  2163. $ CALL DERRGG( C3, NOUT )
  2164. DO 350 I = 1, NPARMS
  2165. CALL XLAENV( 1, NBVAL( I ) )
  2166. CALL XLAENV( 2, NBMIN( I ) )
  2167. CALL XLAENV( 4, NSVAL( I ) )
  2168. CALL XLAENV( 8, MXBVAL( I ) )
  2169. CALL XLAENV( 16, IACC22( I ) )
  2170. CALL XLAENV( 5, NBCOL( I ) )
  2171. *
  2172. IF( NEWSD.EQ.0 ) THEN
  2173. DO 340 K = 1, 4
  2174. ISEED( K ) = IOLDSD( K )
  2175. 340 CONTINUE
  2176. END IF
  2177. WRITE( NOUT, FMT = 9996 )C3, NBVAL( I ), NBMIN( I ),
  2178. $ NSVAL( I ), MXBVAL( I ), IACC22( I ), NBCOL( I )
  2179. TSTDIF = .FALSE.
  2180. THRSHN = 10.D0
  2181. IF( TSTCHK ) THEN
  2182. CALL DCHKGG( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH,
  2183. $ TSTDIF, THRSHN, NOUT, A( 1, 1 ), NMAX,
  2184. $ A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), A( 1, 5 ),
  2185. $ A( 1, 6 ), A( 1, 7 ), A( 1, 8 ), A( 1, 9 ),
  2186. $ NMAX, A( 1, 10 ), A( 1, 11 ), A( 1, 12 ),
  2187. $ D( 1, 1 ), D( 1, 2 ), D( 1, 3 ), D( 1, 4 ),
  2188. $ D( 1, 5 ), D( 1, 6 ), A( 1, 13 ),
  2189. $ A( 1, 14 ), WORK, LWORK, LOGWRK, RESULT,
  2190. $ INFO )
  2191. IF( INFO.NE.0 )
  2192. $ WRITE( NOUT, FMT = 9980 )'DCHKGG', INFO
  2193. END IF
  2194. 350 CONTINUE
  2195. *
  2196. ELSE IF( LSAMEN( 3, C3, 'DGS' ) ) THEN
  2197. *
  2198. * -------------------------------------------------
  2199. * DGS: Generalized Nonsymmetric Eigenvalue Problem
  2200. * DGGES (Schur form)
  2201. * -------------------------------------------------
  2202. *
  2203. MAXTYP = 26
  2204. NTYPES = MIN( MAXTYP, NTYPES )
  2205. IF( NTYPES.LE.0 ) THEN
  2206. WRITE( NOUT, FMT = 9990 )C3
  2207. ELSE
  2208. IF( TSTERR )
  2209. $ CALL DERRGG( C3, NOUT )
  2210. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2211. CALL DDRGES( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH, NOUT,
  2212. $ A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  2213. $ A( 1, 4 ), A( 1, 7 ), NMAX, A( 1, 8 ),
  2214. $ D( 1, 1 ), D( 1, 2 ), D( 1, 3 ), WORK, LWORK,
  2215. $ RESULT, LOGWRK, INFO )
  2216. IF( INFO.NE.0 )
  2217. $ WRITE( NOUT, FMT = 9980 )'DDRGES', INFO
  2218. *
  2219. * Blocked version
  2220. *
  2221. CALL XLAENV(16, 2)
  2222. CALL DDRGES3( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH, NOUT,
  2223. $ A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  2224. $ A( 1, 4 ), A( 1, 7 ), NMAX, A( 1, 8 ),
  2225. $ D( 1, 1 ), D( 1, 2 ), D( 1, 3 ), WORK, LWORK,
  2226. $ RESULT, LOGWRK, INFO )
  2227. IF( INFO.NE.0 )
  2228. $ WRITE( NOUT, FMT = 9980 )'DDRGES3', INFO
  2229. END IF
  2230. WRITE( NOUT, FMT = 9973 )
  2231. GO TO 10
  2232. *
  2233. ELSE IF( DGX ) THEN
  2234. *
  2235. * -------------------------------------------------
  2236. * DGX: Generalized Nonsymmetric Eigenvalue Problem
  2237. * DGGESX (Schur form and condition numbers)
  2238. * -------------------------------------------------
  2239. *
  2240. MAXTYP = 5
  2241. NTYPES = MAXTYP
  2242. IF( NN.LT.0 ) THEN
  2243. WRITE( NOUT, FMT = 9990 )C3
  2244. ELSE
  2245. IF( TSTERR )
  2246. $ CALL DERRGG( C3, NOUT )
  2247. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2248. CALL XLAENV( 5, 2 )
  2249. CALL DDRGSX( NN, NCMAX, THRESH, NIN, NOUT, A( 1, 1 ), NMAX,
  2250. $ A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), A( 1, 5 ),
  2251. $ A( 1, 6 ), D( 1, 1 ), D( 1, 2 ), D( 1, 3 ),
  2252. $ C( 1, 1 ), NCMAX*NCMAX, A( 1, 12 ), WORK,
  2253. $ LWORK, IWORK, LIWORK, LOGWRK, INFO )
  2254. IF( INFO.NE.0 )
  2255. $ WRITE( NOUT, FMT = 9980 )'DDRGSX', INFO
  2256. END IF
  2257. WRITE( NOUT, FMT = 9973 )
  2258. GO TO 10
  2259. *
  2260. ELSE IF( LSAMEN( 3, C3, 'DGV' ) ) THEN
  2261. *
  2262. * -------------------------------------------------
  2263. * DGV: Generalized Nonsymmetric Eigenvalue Problem
  2264. * DGGEV (Eigenvalue/vector form)
  2265. * -------------------------------------------------
  2266. *
  2267. MAXTYP = 26
  2268. NTYPES = MIN( MAXTYP, NTYPES )
  2269. IF( NTYPES.LE.0 ) THEN
  2270. WRITE( NOUT, FMT = 9990 )C3
  2271. ELSE
  2272. IF( TSTERR )
  2273. $ CALL DERRGG( C3, NOUT )
  2274. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2275. CALL DDRGEV( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH, NOUT,
  2276. $ A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  2277. $ A( 1, 4 ), A( 1, 7 ), NMAX, A( 1, 8 ),
  2278. $ A( 1, 9 ), NMAX, D( 1, 1 ), D( 1, 2 ),
  2279. $ D( 1, 3 ), D( 1, 4 ), D( 1, 5 ), D( 1, 6 ),
  2280. $ WORK, LWORK, RESULT, INFO )
  2281. IF( INFO.NE.0 )
  2282. $ WRITE( NOUT, FMT = 9980 )'DDRGEV', INFO
  2283. *
  2284. * Blocked version
  2285. *
  2286. CALL DDRGEV3( NN, NVAL, MAXTYP, DOTYPE, ISEED, THRESH, NOUT,
  2287. $ A( 1, 1 ), NMAX, A( 1, 2 ), A( 1, 3 ),
  2288. $ A( 1, 4 ), A( 1, 7 ), NMAX, A( 1, 8 ),
  2289. $ A( 1, 9 ), NMAX, D( 1, 1 ), D( 1, 2 ),
  2290. $ D( 1, 3 ), D( 1, 4 ), D( 1, 5 ), D( 1, 6 ),
  2291. $ WORK, LWORK, RESULT, INFO )
  2292. IF( INFO.NE.0 )
  2293. $ WRITE( NOUT, FMT = 9980 )'DDRGEV3', INFO
  2294. END IF
  2295. WRITE( NOUT, FMT = 9973 )
  2296. GO TO 10
  2297. *
  2298. ELSE IF( DXV ) THEN
  2299. *
  2300. * -------------------------------------------------
  2301. * DXV: Generalized Nonsymmetric Eigenvalue Problem
  2302. * DGGEVX (eigenvalue/vector with condition numbers)
  2303. * -------------------------------------------------
  2304. *
  2305. MAXTYP = 2
  2306. NTYPES = MAXTYP
  2307. IF( NN.LT.0 ) THEN
  2308. WRITE( NOUT, FMT = 9990 )C3
  2309. ELSE
  2310. IF( TSTERR )
  2311. $ CALL DERRGG( C3, NOUT )
  2312. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2313. CALL DDRGVX( NN, THRESH, NIN, NOUT, A( 1, 1 ), NMAX,
  2314. $ A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), D( 1, 1 ),
  2315. $ D( 1, 2 ), D( 1, 3 ), A( 1, 5 ), A( 1, 6 ),
  2316. $ IWORK( 1 ), IWORK( 2 ), D( 1, 4 ), D( 1, 5 ),
  2317. $ D( 1, 6 ), D( 1, 7 ), D( 1, 8 ), D( 1, 9 ),
  2318. $ WORK, LWORK, IWORK( 3 ), LIWORK-2, RESULT,
  2319. $ LOGWRK, INFO )
  2320. *
  2321. IF( INFO.NE.0 )
  2322. $ WRITE( NOUT, FMT = 9980 )'DDRGVX', INFO
  2323. END IF
  2324. WRITE( NOUT, FMT = 9973 )
  2325. GO TO 10
  2326. *
  2327. ELSE IF( LSAMEN( 3, C3, 'DSB' ) ) THEN
  2328. *
  2329. * ------------------------------
  2330. * DSB: Symmetric Band Reduction
  2331. * ------------------------------
  2332. *
  2333. MAXTYP = 15
  2334. NTYPES = MIN( MAXTYP, NTYPES )
  2335. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2336. IF( TSTERR )
  2337. $ CALL DERRST( 'DSB', NOUT )
  2338. * CALL DCHKSB( NN, NVAL, NK, KVAL, MAXTYP, DOTYPE, ISEED, THRESH,
  2339. * $ NOUT, A( 1, 1 ), NMAX, D( 1, 1 ), D( 1, 2 ),
  2340. * $ A( 1, 2 ), NMAX, WORK, LWORK, RESULT, INFO )
  2341. CALL DCHKSB2STG( NN, NVAL, NK, KVAL, MAXTYP, DOTYPE, ISEED,
  2342. $ THRESH, NOUT, A( 1, 1 ), NMAX, D( 1, 1 ),
  2343. $ D( 1, 2 ), D( 1, 3 ), D( 1, 4 ), D( 1, 5 ),
  2344. $ A( 1, 2 ), NMAX, WORK, LWORK, RESULT, INFO )
  2345. IF( INFO.NE.0 )
  2346. $ WRITE( NOUT, FMT = 9980 )'DCHKSB', INFO
  2347. *
  2348. ELSE IF( LSAMEN( 3, C3, 'DBB' ) ) THEN
  2349. *
  2350. * ------------------------------
  2351. * DBB: General Band Reduction
  2352. * ------------------------------
  2353. *
  2354. MAXTYP = 15
  2355. NTYPES = MIN( MAXTYP, NTYPES )
  2356. CALL ALAREQ( C3, NTYPES, DOTYPE, MAXTYP, NIN, NOUT )
  2357. DO 370 I = 1, NPARMS
  2358. NRHS = NSVAL( I )
  2359. *
  2360. IF( NEWSD.EQ.0 ) THEN
  2361. DO 360 K = 1, 4
  2362. ISEED( K ) = IOLDSD( K )
  2363. 360 CONTINUE
  2364. END IF
  2365. WRITE( NOUT, FMT = 9966 )C3, NRHS
  2366. CALL DCHKBB( NN, MVAL, NVAL, NK, KVAL, MAXTYP, DOTYPE, NRHS,
  2367. $ ISEED, THRESH, NOUT, A( 1, 1 ), NMAX,
  2368. $ A( 1, 2 ), 2*NMAX, D( 1, 1 ), D( 1, 2 ),
  2369. $ A( 1, 4 ), NMAX, A( 1, 5 ), NMAX, A( 1, 6 ),
  2370. $ NMAX, A( 1, 7 ), WORK, LWORK, RESULT, INFO )
  2371. IF( INFO.NE.0 )
  2372. $ WRITE( NOUT, FMT = 9980 )'DCHKBB', INFO
  2373. 370 CONTINUE
  2374. *
  2375. ELSE IF( LSAMEN( 3, C3, 'GLM' ) ) THEN
  2376. *
  2377. * -----------------------------------------
  2378. * GLM: Generalized Linear Regression Model
  2379. * -----------------------------------------
  2380. *
  2381. CALL XLAENV( 1, 1 )
  2382. IF( TSTERR )
  2383. $ CALL DERRGG( 'GLM', NOUT )
  2384. CALL DCKGLM( NN, MVAL, PVAL, NVAL, NTYPES, ISEED, THRESH, NMAX,
  2385. $ A( 1, 1 ), A( 1, 2 ), B( 1, 1 ), B( 1, 2 ), X,
  2386. $ WORK, D( 1, 1 ), NIN, NOUT, INFO )
  2387. IF( INFO.NE.0 )
  2388. $ WRITE( NOUT, FMT = 9980 )'DCKGLM', INFO
  2389. *
  2390. ELSE IF( LSAMEN( 3, C3, 'GQR' ) ) THEN
  2391. *
  2392. * ------------------------------------------
  2393. * GQR: Generalized QR and RQ factorizations
  2394. * ------------------------------------------
  2395. *
  2396. CALL XLAENV( 1, 1 )
  2397. IF( TSTERR )
  2398. $ CALL DERRGG( 'GQR', NOUT )
  2399. CALL DCKGQR( NN, MVAL, NN, PVAL, NN, NVAL, NTYPES, ISEED,
  2400. $ THRESH, NMAX, A( 1, 1 ), A( 1, 2 ), A( 1, 3 ),
  2401. $ A( 1, 4 ), TAUA, B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
  2402. $ B( 1, 4 ), B( 1, 5 ), TAUB, WORK, D( 1, 1 ), NIN,
  2403. $ NOUT, INFO )
  2404. IF( INFO.NE.0 )
  2405. $ WRITE( NOUT, FMT = 9980 )'DCKGQR', INFO
  2406. *
  2407. ELSE IF( LSAMEN( 3, C3, 'GSV' ) ) THEN
  2408. *
  2409. * ----------------------------------------------
  2410. * GSV: Generalized Singular Value Decomposition
  2411. * ----------------------------------------------
  2412. *
  2413. CALL XLAENV(1,1)
  2414. IF( TSTERR )
  2415. $ CALL DERRGG( 'GSV', NOUT )
  2416. CALL DCKGSV( NN, MVAL, PVAL, NVAL, NTYPES, ISEED, THRESH, NMAX,
  2417. $ A( 1, 1 ), A( 1, 2 ), B( 1, 1 ), B( 1, 2 ),
  2418. $ A( 1, 3 ), B( 1, 3 ), A( 1, 4 ), TAUA, TAUB,
  2419. $ B( 1, 4 ), IWORK, WORK, D( 1, 1 ), NIN, NOUT,
  2420. $ INFO )
  2421. IF( INFO.NE.0 )
  2422. $ WRITE( NOUT, FMT = 9980 )'DCKGSV', INFO
  2423. *
  2424. ELSE IF( LSAMEN( 3, C3, 'CSD' ) ) THEN
  2425. *
  2426. * ----------------------------------------------
  2427. * CSD: CS Decomposition
  2428. * ----------------------------------------------
  2429. *
  2430. CALL XLAENV(1,1)
  2431. IF( TSTERR )
  2432. $ CALL DERRGG( 'CSD', NOUT )
  2433. CALL DCKCSD( NN, MVAL, PVAL, NVAL, NTYPES, ISEED, THRESH, NMAX,
  2434. $ A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), A( 1, 4 ),
  2435. $ A( 1, 5 ), A( 1, 6 ), A( 1, 7 ), IWORK, WORK,
  2436. $ D( 1, 1 ), NIN, NOUT, INFO )
  2437. IF( INFO.NE.0 )
  2438. $ WRITE( NOUT, FMT = 9980 )'DCKCSD', INFO
  2439. *
  2440. ELSE IF( LSAMEN( 3, C3, 'LSE' ) ) THEN
  2441. *
  2442. * --------------------------------------
  2443. * LSE: Constrained Linear Least Squares
  2444. * --------------------------------------
  2445. *
  2446. CALL XLAENV( 1, 1 )
  2447. IF( TSTERR )
  2448. $ CALL DERRGG( 'LSE', NOUT )
  2449. CALL DCKLSE( NN, MVAL, PVAL, NVAL, NTYPES, ISEED, THRESH, NMAX,
  2450. $ A( 1, 1 ), A( 1, 2 ), B( 1, 1 ), B( 1, 2 ), X,
  2451. $ WORK, D( 1, 1 ), NIN, NOUT, INFO )
  2452. IF( INFO.NE.0 )
  2453. $ WRITE( NOUT, FMT = 9980 )'DCKLSE', INFO
  2454. *
  2455. ELSE
  2456. WRITE( NOUT, FMT = * )
  2457. WRITE( NOUT, FMT = * )
  2458. WRITE( NOUT, FMT = 9992 )C3
  2459. END IF
  2460. IF( .NOT.( DGX .OR. DXV ) )
  2461. $ GO TO 190
  2462. 380 CONTINUE
  2463. WRITE( NOUT, FMT = 9994 )
  2464. S2 = DSECND( )
  2465. WRITE( NOUT, FMT = 9993 )S2 - S1
  2466. *
  2467. DEALLOCATE (A, STAT = AllocateStatus)
  2468. DEALLOCATE (B, STAT = AllocateStatus)
  2469. DEALLOCATE (C, STAT = AllocateStatus)
  2470. DEALLOCATE (WORK, STAT = AllocateStatus)
  2471. *
  2472. 9999 FORMAT( / ' Execution not attempted due to input errors' )
  2473. 9997 FORMAT( / / 1X, A3, ': NB =', I4, ', NBMIN =', I4, ', NX =', I4 )
  2474. 9996 FORMAT( / / 1X, A3, ': NB =', I4, ', NBMIN =', I4, ', NS =', I4,
  2475. $ ', MAXB =', I4, ', IACC22 =', I4, ', NBCOL =', I4 )
  2476. 9995 FORMAT( / / 1X, A3, ': NB =', I4, ', NBMIN =', I4, ', NX =', I4,
  2477. $ ', NRHS =', I4 )
  2478. 9994 FORMAT( / / ' End of tests' )
  2479. 9993 FORMAT( ' Total time used = ', F12.2, ' seconds', / )
  2480. 9992 FORMAT( 1X, A3, ': Unrecognized path name' )
  2481. 9991 FORMAT( / / ' *** Invalid integer value in column ', I2,
  2482. $ ' of input', ' line:', / A79 )
  2483. 9990 FORMAT( / / 1X, A3, ' routines were not tested' )
  2484. 9989 FORMAT( ' Invalid input value: ', A, '=', I6, '; must be >=',
  2485. $ I6 )
  2486. 9988 FORMAT( ' Invalid input value: ', A, '=', I6, '; must be <=',
  2487. $ I6 )
  2488. 9987 FORMAT( ' Tests of the Nonsymmetric Eigenvalue Problem routines' )
  2489. 9986 FORMAT( ' Tests of the Symmetric Eigenvalue Problem routines' )
  2490. 9985 FORMAT( ' Tests of the Singular Value Decomposition routines' )
  2491. 9984 FORMAT( / ' The following parameter values will be used:' )
  2492. 9983 FORMAT( 4X, A, 10I6, / 10X, 10I6 )
  2493. 9982 FORMAT( / ' Routines pass computational tests if test ratio is ',
  2494. $ 'less than', F8.2, / )
  2495. 9981 FORMAT( ' Relative machine ', A, ' is taken to be', D16.6 )
  2496. 9980 FORMAT( ' *** Error code from ', A, ' = ', I4 )
  2497. 9979 FORMAT( / ' Tests of the Nonsymmetric Eigenvalue Problem Driver',
  2498. $ / ' DGEEV (eigenvalues and eigevectors)' )
  2499. 9978 FORMAT( / ' Tests of the Nonsymmetric Eigenvalue Problem Driver',
  2500. $ / ' DGEES (Schur form)' )
  2501. 9977 FORMAT( / ' Tests of the Nonsymmetric Eigenvalue Problem Expert',
  2502. $ ' Driver', / ' DGEEVX (eigenvalues, eigenvectors and',
  2503. $ ' condition numbers)' )
  2504. 9976 FORMAT( / ' Tests of the Nonsymmetric Eigenvalue Problem Expert',
  2505. $ ' Driver', / ' DGEESX (Schur form and condition',
  2506. $ ' numbers)' )
  2507. 9975 FORMAT( / ' Tests of the Generalized Nonsymmetric Eigenvalue ',
  2508. $ 'Problem routines' )
  2509. 9974 FORMAT( ' Tests of DSBTRD', / ' (reduction of a symmetric band ',
  2510. $ 'matrix to tridiagonal form)' )
  2511. 9973 FORMAT( / 1X, 71( '-' ) )
  2512. 9972 FORMAT( / ' LAPACK VERSION ', I1, '.', I1, '.', I1 )
  2513. 9971 FORMAT( / ' Tests of the Generalized Linear Regression Model ',
  2514. $ 'routines' )
  2515. 9970 FORMAT( / ' Tests of the Generalized QR and RQ routines' )
  2516. 9969 FORMAT( / ' Tests of the Generalized Singular Value',
  2517. $ ' Decomposition routines' )
  2518. 9968 FORMAT( / ' Tests of the Linear Least Squares routines' )
  2519. 9967 FORMAT( ' Tests of DGBBRD', / ' (reduction of a general band ',
  2520. $ 'matrix to real bidiagonal form)' )
  2521. 9966 FORMAT( / / 1X, A3, ': NRHS =', I4 )
  2522. 9965 FORMAT( / ' Tests of the Generalized Nonsymmetric Eigenvalue ',
  2523. $ 'Problem Expert Driver DGGESX' )
  2524. 9964 FORMAT( / ' Tests of the Generalized Nonsymmetric Eigenvalue ',
  2525. $ 'Problem Driver DGGES' )
  2526. 9963 FORMAT( / ' Tests of the Generalized Nonsymmetric Eigenvalue ',
  2527. $ 'Problem Driver DGGEV' )
  2528. 9962 FORMAT( / ' Tests of the Generalized Nonsymmetric Eigenvalue ',
  2529. $ 'Problem Expert Driver DGGEVX' )
  2530. 9961 FORMAT( / / 1X, A3, ': NB =', I4, ', NBMIN =', I4, ', NX =', I4,
  2531. $ ', INMIN=', I4,
  2532. $ ', INWIN =', I4, ', INIBL =', I4, ', ISHFTS =', I4,
  2533. $ ', IACC22 =', I4)
  2534. 9960 FORMAT( / ' Tests of the CS Decomposition routines' )
  2535. *
  2536. * End of DCHKEE
  2537. *
  2538. END