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.

ddrgev3.f 35 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972
  1. *> \brief \b DDRGEV3
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. * Definition:
  9. * ===========
  10. *
  11. * SUBROUTINE DDRGEV3( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
  12. * NOUNIT, A, LDA, B, S, T, Q, LDQ, Z, QE, LDQE,
  13. * ALPHAR, ALPHAI, BETA, ALPHR1, ALPHI1, BETA1,
  14. * WORK, LWORK, RESULT, INFO )
  15. *
  16. * .. Scalar Arguments ..
  17. * INTEGER INFO, LDA, LDQ, LDQE, LWORK, NOUNIT, NSIZES,
  18. * $ NTYPES
  19. * DOUBLE PRECISION THRESH
  20. * ..
  21. * .. Array Arguments ..
  22. * LOGICAL DOTYPE( * )
  23. * INTEGER ISEED( 4 ), NN( * )
  24. * DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
  25. * $ ALPHI1( * ), ALPHR1( * ), B( LDA, * ),
  26. * $ BETA( * ), BETA1( * ), Q( LDQ, * ),
  27. * $ QE( LDQE, * ), RESULT( * ), S( LDA, * ),
  28. * $ T( LDA, * ), WORK( * ), Z( LDQ, * )
  29. * ..
  30. *
  31. *
  32. *> \par Purpose:
  33. * =============
  34. *>
  35. *> \verbatim
  36. *>
  37. *> DDRGEV3 checks the nonsymmetric generalized eigenvalue problem driver
  38. *> routine DGGEV3.
  39. *>
  40. *> DGGEV3 computes for a pair of n-by-n nonsymmetric matrices (A,B) the
  41. *> generalized eigenvalues and, optionally, the left and right
  42. *> eigenvectors.
  43. *>
  44. *> A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
  45. *> or a ratio alpha/beta = w, such that A - w*B is singular. It is
  46. *> usually represented as the pair (alpha,beta), as there is reasonable
  47. *> interpretation for beta=0, and even for both being zero.
  48. *>
  49. *> A right generalized eigenvector corresponding to a generalized
  50. *> eigenvalue w for a pair of matrices (A,B) is a vector r such that
  51. *> (A - wB) * r = 0. A left generalized eigenvector is a vector l such
  52. *> that l**H * (A - wB) = 0, where l**H is the conjugate-transpose of l.
  53. *>
  54. *> When DDRGEV3 is called, a number of matrix "sizes" ("n's") and a
  55. *> number of matrix "types" are specified. For each size ("n")
  56. *> and each type of matrix, a pair of matrices (A, B) will be generated
  57. *> and used for testing. For each matrix pair, the following tests
  58. *> will be performed and compared with the threshold THRESH.
  59. *>
  60. *> Results from DGGEV3:
  61. *>
  62. *> (1) max over all left eigenvalue/-vector pairs (alpha/beta,l) of
  63. *>
  64. *> | VL**H * (beta A - alpha B) |/( ulp max(|beta A|, |alpha B|) )
  65. *>
  66. *> where VL**H is the conjugate-transpose of VL.
  67. *>
  68. *> (2) | |VL(i)| - 1 | / ulp and whether largest component real
  69. *>
  70. *> VL(i) denotes the i-th column of VL.
  71. *>
  72. *> (3) max over all left eigenvalue/-vector pairs (alpha/beta,r) of
  73. *>
  74. *> | (beta A - alpha B) * VR | / ( ulp max(|beta A|, |alpha B|) )
  75. *>
  76. *> (4) | |VR(i)| - 1 | / ulp and whether largest component real
  77. *>
  78. *> VR(i) denotes the i-th column of VR.
  79. *>
  80. *> (5) W(full) = W(partial)
  81. *> W(full) denotes the eigenvalues computed when both l and r
  82. *> are also computed, and W(partial) denotes the eigenvalues
  83. *> computed when only W, only W and r, or only W and l are
  84. *> computed.
  85. *>
  86. *> (6) VL(full) = VL(partial)
  87. *> VL(full) denotes the left eigenvectors computed when both l
  88. *> and r are computed, and VL(partial) denotes the result
  89. *> when only l is computed.
  90. *>
  91. *> (7) VR(full) = VR(partial)
  92. *> VR(full) denotes the right eigenvectors computed when both l
  93. *> and r are also computed, and VR(partial) denotes the result
  94. *> when only l is computed.
  95. *>
  96. *>
  97. *> Test Matrices
  98. *> ---- --------
  99. *>
  100. *> The sizes of the test matrices are specified by an array
  101. *> NN(1:NSIZES); the value of each element NN(j) specifies one size.
  102. *> The "types" are specified by a logical array DOTYPE( 1:NTYPES ); if
  103. *> DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
  104. *> Currently, the list of possible types is:
  105. *>
  106. *> (1) ( 0, 0 ) (a pair of zero matrices)
  107. *>
  108. *> (2) ( I, 0 ) (an identity and a zero matrix)
  109. *>
  110. *> (3) ( 0, I ) (an identity and a zero matrix)
  111. *>
  112. *> (4) ( I, I ) (a pair of identity matrices)
  113. *>
  114. *> t t
  115. *> (5) ( J , J ) (a pair of transposed Jordan blocks)
  116. *>
  117. *> t ( I 0 )
  118. *> (6) ( X, Y ) where X = ( J 0 ) and Y = ( t )
  119. *> ( 0 I ) ( 0 J )
  120. *> and I is a k x k identity and J a (k+1)x(k+1)
  121. *> Jordan block; k=(N-1)/2
  122. *>
  123. *> (7) ( D, I ) where D is diag( 0, 1,..., N-1 ) (a diagonal
  124. *> matrix with those diagonal entries.)
  125. *> (8) ( I, D )
  126. *>
  127. *> (9) ( big*D, small*I ) where "big" is near overflow and small=1/big
  128. *>
  129. *> (10) ( small*D, big*I )
  130. *>
  131. *> (11) ( big*I, small*D )
  132. *>
  133. *> (12) ( small*I, big*D )
  134. *>
  135. *> (13) ( big*D, big*I )
  136. *>
  137. *> (14) ( small*D, small*I )
  138. *>
  139. *> (15) ( D1, D2 ) where D1 is diag( 0, 0, 1, ..., N-3, 0 ) and
  140. *> D2 is diag( 0, N-3, N-4,..., 1, 0, 0 )
  141. *> t t
  142. *> (16) Q ( J , J ) Z where Q and Z are random orthogonal matrices.
  143. *>
  144. *> (17) Q ( T1, T2 ) Z where T1 and T2 are upper triangular matrices
  145. *> with random O(1) entries above the diagonal
  146. *> and diagonal entries diag(T1) =
  147. *> ( 0, 0, 1, ..., N-3, 0 ) and diag(T2) =
  148. *> ( 0, N-3, N-4,..., 1, 0, 0 )
  149. *>
  150. *> (18) Q ( T1, T2 ) Z diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
  151. *> diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
  152. *> s = machine precision.
  153. *>
  154. *> (19) Q ( T1, T2 ) Z diag(T1)=( 0,0,1,1, 1-d, ..., 1-(N-5)*d=s, 0 )
  155. *> diag(T2) = ( 0, 1, 0, 1, ..., 1, 0 )
  156. *>
  157. *> N-5
  158. *> (20) Q ( T1, T2 ) Z diag(T1)=( 0, 0, 1, 1, a, ..., a =s, 0 )
  159. *> diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
  160. *>
  161. *> (21) Q ( T1, T2 ) Z diag(T1)=( 0, 0, 1, r1, r2, ..., r(N-4), 0 )
  162. *> diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
  163. *> where r1,..., r(N-4) are random.
  164. *>
  165. *> (22) Q ( big*T1, small*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
  166. *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
  167. *>
  168. *> (23) Q ( small*T1, big*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
  169. *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
  170. *>
  171. *> (24) Q ( small*T1, small*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
  172. *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
  173. *>
  174. *> (25) Q ( big*T1, big*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
  175. *> diag(T2) = ( 0, 1, ..., 1, 0, 0 )
  176. *>
  177. *> (26) Q ( T1, T2 ) Z where T1 and T2 are random upper-triangular
  178. *> matrices.
  179. *>
  180. *> \endverbatim
  181. *
  182. * Arguments:
  183. * ==========
  184. *
  185. *> \param[in] NSIZES
  186. *> \verbatim
  187. *> NSIZES is INTEGER
  188. *> The number of sizes of matrices to use. If it is zero,
  189. *> DDRGEV3 does nothing. NSIZES >= 0.
  190. *> \endverbatim
  191. *>
  192. *> \param[in] NN
  193. *> \verbatim
  194. *> NN is INTEGER array, dimension (NSIZES)
  195. *> An array containing the sizes to be used for the matrices.
  196. *> Zero values will be skipped. NN >= 0.
  197. *> \endverbatim
  198. *>
  199. *> \param[in] NTYPES
  200. *> \verbatim
  201. *> NTYPES is INTEGER
  202. *> The number of elements in DOTYPE. If it is zero, DDRGEV3
  203. *> does nothing. It must be at least zero. If it is MAXTYP+1
  204. *> and NSIZES is 1, then an additional type, MAXTYP+1 is
  205. *> defined, which is to use whatever matrix is in A. This
  206. *> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
  207. *> DOTYPE(MAXTYP+1) is .TRUE. .
  208. *> \endverbatim
  209. *>
  210. *> \param[in] DOTYPE
  211. *> \verbatim
  212. *> DOTYPE is LOGICAL array, dimension (NTYPES)
  213. *> If DOTYPE(j) is .TRUE., then for each size in NN a
  214. *> matrix of that size and of type j will be generated.
  215. *> If NTYPES is smaller than the maximum number of types
  216. *> defined (PARAMETER MAXTYP), then types NTYPES+1 through
  217. *> MAXTYP will not be generated. If NTYPES is larger
  218. *> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
  219. *> will be ignored.
  220. *> \endverbatim
  221. *>
  222. *> \param[in,out] ISEED
  223. *> \verbatim
  224. *> ISEED is INTEGER array, dimension (4)
  225. *> On entry ISEED specifies the seed of the random number
  226. *> generator. The array elements should be between 0 and 4095;
  227. *> if not they will be reduced mod 4096. Also, ISEED(4) must
  228. *> be odd. The random number generator uses a linear
  229. *> congruential sequence limited to small integers, and so
  230. *> should produce machine independent random numbers. The
  231. *> values of ISEED are changed on exit, and can be used in the
  232. *> next call to DDRGEV3 to continue the same random number
  233. *> sequence.
  234. *> \endverbatim
  235. *>
  236. *> \param[in] THRESH
  237. *> \verbatim
  238. *> THRESH is DOUBLE PRECISION
  239. *> A test will count as "failed" if the "error", computed as
  240. *> described above, exceeds THRESH. Note that the error is
  241. *> scaled to be O(1), so THRESH should be a reasonably small
  242. *> multiple of 1, e.g., 10 or 100. In particular, it should
  243. *> not depend on the precision (single vs. double) or the size
  244. *> of the matrix. It must be at least zero.
  245. *> \endverbatim
  246. *>
  247. *> \param[in] NOUNIT
  248. *> \verbatim
  249. *> NOUNIT is INTEGER
  250. *> The FORTRAN unit number for printing out error messages
  251. *> (e.g., if a routine returns IERR not equal to 0.)
  252. *> \endverbatim
  253. *>
  254. *> \param[in,out] A
  255. *> \verbatim
  256. *> A is DOUBLE PRECISION array,
  257. *> dimension(LDA, max(NN))
  258. *> Used to hold the original A matrix. Used as input only
  259. *> if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
  260. *> DOTYPE(MAXTYP+1)=.TRUE.
  261. *> \endverbatim
  262. *>
  263. *> \param[in] LDA
  264. *> \verbatim
  265. *> LDA is INTEGER
  266. *> The leading dimension of A, B, S, and T.
  267. *> It must be at least 1 and at least max( NN ).
  268. *> \endverbatim
  269. *>
  270. *> \param[in,out] B
  271. *> \verbatim
  272. *> B is DOUBLE PRECISION array,
  273. *> dimension(LDA, max(NN))
  274. *> Used to hold the original B matrix. Used as input only
  275. *> if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
  276. *> DOTYPE(MAXTYP+1)=.TRUE.
  277. *> \endverbatim
  278. *>
  279. *> \param[out] S
  280. *> \verbatim
  281. *> S is DOUBLE PRECISION array,
  282. *> dimension (LDA, max(NN))
  283. *> The Schur form matrix computed from A by DGGEV3. On exit, S
  284. *> contains the Schur form matrix corresponding to the matrix
  285. *> in A.
  286. *> \endverbatim
  287. *>
  288. *> \param[out] T
  289. *> \verbatim
  290. *> T is DOUBLE PRECISION array,
  291. *> dimension (LDA, max(NN))
  292. *> The upper triangular matrix computed from B by DGGEV3.
  293. *> \endverbatim
  294. *>
  295. *> \param[out] Q
  296. *> \verbatim
  297. *> Q is DOUBLE PRECISION array,
  298. *> dimension (LDQ, max(NN))
  299. *> The (left) eigenvectors matrix computed by DGGEV3.
  300. *> \endverbatim
  301. *>
  302. *> \param[in] LDQ
  303. *> \verbatim
  304. *> LDQ is INTEGER
  305. *> The leading dimension of Q and Z. It must
  306. *> be at least 1 and at least max( NN ).
  307. *> \endverbatim
  308. *>
  309. *> \param[out] Z
  310. *> \verbatim
  311. *> Z is DOUBLE PRECISION array, dimension( LDQ, max(NN) )
  312. *> The (right) orthogonal matrix computed by DGGEV3.
  313. *> \endverbatim
  314. *>
  315. *> \param[out] QE
  316. *> \verbatim
  317. *> QE is DOUBLE PRECISION array, dimension( LDQ, max(NN) )
  318. *> QE holds the computed right or left eigenvectors.
  319. *> \endverbatim
  320. *>
  321. *> \param[in] LDQE
  322. *> \verbatim
  323. *> LDQE is INTEGER
  324. *> The leading dimension of QE. LDQE >= max(1,max(NN)).
  325. *> \endverbatim
  326. *>
  327. *> \param[out] ALPHAR
  328. *> \verbatim
  329. *> ALPHAR is DOUBLE PRECISION array, dimension (max(NN))
  330. *> \endverbatim
  331. *>
  332. *> \param[out] ALPHAI
  333. *> \verbatim
  334. *> ALPHAI is DOUBLE PRECISION array, dimension (max(NN))
  335. *> \endverbatim
  336. *>
  337. *> \param[out] BETA
  338. *> \verbatim
  339. *> BETA is DOUBLE PRECISION array, dimension (max(NN))
  340. *>
  341. *> The generalized eigenvalues of (A,B) computed by DGGEV3.
  342. *> ( ALPHAR(k)+ALPHAI(k)*i ) / BETA(k) is the k-th
  343. *> generalized eigenvalue of A and B.
  344. *> \endverbatim
  345. *>
  346. *> \param[out] ALPHR1
  347. *> \verbatim
  348. *> ALPHR1 is DOUBLE PRECISION array, dimension (max(NN))
  349. *> \endverbatim
  350. *>
  351. *> \param[out] ALPHI1
  352. *> \verbatim
  353. *> ALPHI1 is DOUBLE PRECISION array, dimension (max(NN))
  354. *> \endverbatim
  355. *>
  356. *> \param[out] BETA1
  357. *> \verbatim
  358. *> BETA1 is DOUBLE PRECISION array, dimension (max(NN))
  359. *>
  360. *> Like ALPHAR, ALPHAI, BETA, these arrays contain the
  361. *> eigenvalues of A and B, but those computed when DGGEV3 only
  362. *> computes a partial eigendecomposition, i.e. not the
  363. *> eigenvalues and left and right eigenvectors.
  364. *> \endverbatim
  365. *>
  366. *> \param[out] WORK
  367. *> \verbatim
  368. *> WORK is DOUBLE PRECISION array, dimension (LWORK)
  369. *> \endverbatim
  370. *>
  371. *> \param[in] LWORK
  372. *> \verbatim
  373. *> LWORK is INTEGER
  374. *> The number of entries in WORK. LWORK >= MAX( 8*N, N*(N+1) ).
  375. *> \endverbatim
  376. *>
  377. *> \param[out] RESULT
  378. *> \verbatim
  379. *> RESULT is DOUBLE PRECISION array, dimension (2)
  380. *> The values computed by the tests described above.
  381. *> The values are currently limited to 1/ulp, to avoid overflow.
  382. *> \endverbatim
  383. *>
  384. *> \param[out] INFO
  385. *> \verbatim
  386. *> INFO is INTEGER
  387. *> = 0: successful exit
  388. *> < 0: if INFO = -i, the i-th argument had an illegal value.
  389. *> > 0: A routine returned an error code. INFO is the
  390. *> absolute value of the INFO value returned.
  391. *> \endverbatim
  392. *
  393. * Authors:
  394. * ========
  395. *
  396. *> \author Univ. of Tennessee
  397. *> \author Univ. of California Berkeley
  398. *> \author Univ. of Colorado Denver
  399. *> \author NAG Ltd.
  400. *
  401. *> \ingroup double_eig
  402. *
  403. * =====================================================================
  404. SUBROUTINE DDRGEV3( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
  405. $ NOUNIT, A, LDA, B, S, T, Q, LDQ, Z, QE, LDQE,
  406. $ ALPHAR, ALPHAI, BETA, ALPHR1, ALPHI1, BETA1,
  407. $ WORK, LWORK, RESULT, INFO )
  408. *
  409. * -- LAPACK test routine --
  410. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  411. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  412. *
  413. * .. Scalar Arguments ..
  414. INTEGER INFO, LDA, LDQ, LDQE, LWORK, NOUNIT, NSIZES,
  415. $ NTYPES
  416. DOUBLE PRECISION THRESH
  417. * ..
  418. * .. Array Arguments ..
  419. LOGICAL DOTYPE( * )
  420. INTEGER ISEED( 4 ), NN( * )
  421. DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
  422. $ ALPHI1( * ), ALPHR1( * ), B( LDA, * ),
  423. $ BETA( * ), BETA1( * ), Q( LDQ, * ),
  424. $ QE( LDQE, * ), RESULT( * ), S( LDA, * ),
  425. $ T( LDA, * ), WORK( * ), Z( LDQ, * )
  426. * ..
  427. *
  428. * =====================================================================
  429. *
  430. * .. Parameters ..
  431. DOUBLE PRECISION ZERO, ONE
  432. PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
  433. INTEGER MAXTYP
  434. PARAMETER ( MAXTYP = 27 )
  435. * ..
  436. * .. Local Scalars ..
  437. LOGICAL BADNN
  438. INTEGER I, IADD, IERR, IN, J, JC, JR, JSIZE, JTYPE,
  439. $ MAXWRK, MINWRK, MTYPES, N, N1, NERRS, NMATS,
  440. $ NMAX, NTESTT
  441. DOUBLE PRECISION SAFMAX, SAFMIN, ULP, ULPINV
  442. * ..
  443. * .. Local Arrays ..
  444. INTEGER IASIGN( MAXTYP ), IBSIGN( MAXTYP ),
  445. $ IOLDSD( 4 ), KADD( 6 ), KAMAGN( MAXTYP ),
  446. $ KATYPE( MAXTYP ), KAZERO( MAXTYP ),
  447. $ KBMAGN( MAXTYP ), KBTYPE( MAXTYP ),
  448. $ KBZERO( MAXTYP ), KCLASS( MAXTYP ),
  449. $ KTRIAN( MAXTYP ), KZ1( 6 ), KZ2( 6 )
  450. DOUBLE PRECISION RMAGN( 0: 3 )
  451. * ..
  452. * .. External Functions ..
  453. INTEGER ILAENV
  454. DOUBLE PRECISION DLAMCH, DLARND
  455. EXTERNAL ILAENV, DLAMCH, DLARND
  456. * ..
  457. * .. External Subroutines ..
  458. EXTERNAL ALASVM, DGET52, DGGEV3, DLABAD, DLACPY, DLARFG,
  459. $ DLASET, DLATM4, DORM2R, XERBLA
  460. * ..
  461. * .. Intrinsic Functions ..
  462. INTRINSIC ABS, DBLE, MAX, MIN, SIGN
  463. * ..
  464. * .. Data statements ..
  465. DATA KCLASS / 15*1, 10*2, 1*3, 1*4 /
  466. DATA KZ1 / 0, 1, 2, 1, 3, 3 /
  467. DATA KZ2 / 0, 0, 1, 2, 1, 1 /
  468. DATA KADD / 0, 0, 0, 0, 3, 2 /
  469. DATA KATYPE / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4,
  470. $ 4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0, 0 /
  471. DATA KBTYPE / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4,
  472. $ 1, 1, -4, 2, -4, 8*8, 0, 0 /
  473. DATA KAZERO / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3,
  474. $ 4*5, 4*3, 1, 1 /
  475. DATA KBZERO / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
  476. $ 4*6, 4*4, 1, 1 /
  477. DATA KAMAGN / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
  478. $ 2, 1, 3 /
  479. DATA KBMAGN / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
  480. $ 2, 1, 3 /
  481. DATA KTRIAN / 16*0, 11*1 /
  482. DATA IASIGN / 6*0, 2, 0, 2*2, 2*0, 3*2, 0, 2, 3*0,
  483. $ 5*2, 2*0 /
  484. DATA IBSIGN / 7*0, 2, 2*0, 2*2, 2*0, 2, 0, 2, 10*0 /
  485. * ..
  486. * .. Executable Statements ..
  487. *
  488. * Check for errors
  489. *
  490. INFO = 0
  491. *
  492. BADNN = .FALSE.
  493. NMAX = 1
  494. DO 10 J = 1, NSIZES
  495. NMAX = MAX( NMAX, NN( J ) )
  496. IF( NN( J ).LT.0 )
  497. $ BADNN = .TRUE.
  498. 10 CONTINUE
  499. *
  500. IF( NSIZES.LT.0 ) THEN
  501. INFO = -1
  502. ELSE IF( BADNN ) THEN
  503. INFO = -2
  504. ELSE IF( NTYPES.LT.0 ) THEN
  505. INFO = -3
  506. ELSE IF( THRESH.LT.ZERO ) THEN
  507. INFO = -6
  508. ELSE IF( LDA.LE.1 .OR. LDA.LT.NMAX ) THEN
  509. INFO = -9
  510. ELSE IF( LDQ.LE.1 .OR. LDQ.LT.NMAX ) THEN
  511. INFO = -14
  512. ELSE IF( LDQE.LE.1 .OR. LDQE.LT.NMAX ) THEN
  513. INFO = -17
  514. END IF
  515. *
  516. * Compute workspace
  517. * (Note: Comments in the code beginning "Workspace:" describe the
  518. * minimal amount of workspace needed at that point in the code,
  519. * as well as the preferred amount for good performance.
  520. * NB refers to the optimal block size for the immediately
  521. * following subroutine, as returned by ILAENV.
  522. *
  523. MINWRK = 1
  524. IF( INFO.EQ.0 .AND. LWORK.GE.1 ) THEN
  525. MINWRK = MAX( 1, 8*NMAX, NMAX*( NMAX+1 ) )
  526. MAXWRK = 7*NMAX + NMAX*ILAENV( 1, 'DGEQRF', ' ', NMAX, 1, NMAX,
  527. $ 0 )
  528. MAXWRK = MAX( MAXWRK, NMAX*( NMAX+1 ) )
  529. WORK( 1 ) = MAXWRK
  530. END IF
  531. *
  532. IF( LWORK.LT.MINWRK )
  533. $ INFO = -25
  534. *
  535. IF( INFO.NE.0 ) THEN
  536. CALL XERBLA( 'DDRGEV3', -INFO )
  537. RETURN
  538. END IF
  539. *
  540. * Quick return if possible
  541. *
  542. IF( NSIZES.EQ.0 .OR. NTYPES.EQ.0 )
  543. $ RETURN
  544. *
  545. SAFMIN = DLAMCH( 'Safe minimum' )
  546. ULP = DLAMCH( 'Epsilon' )*DLAMCH( 'Base' )
  547. SAFMIN = SAFMIN / ULP
  548. SAFMAX = ONE / SAFMIN
  549. CALL DLABAD( SAFMIN, SAFMAX )
  550. ULPINV = ONE / ULP
  551. *
  552. * The values RMAGN(2:3) depend on N, see below.
  553. *
  554. RMAGN( 0 ) = ZERO
  555. RMAGN( 1 ) = ONE
  556. *
  557. * Loop over sizes, types
  558. *
  559. NTESTT = 0
  560. NERRS = 0
  561. NMATS = 0
  562. *
  563. DO 220 JSIZE = 1, NSIZES
  564. N = NN( JSIZE )
  565. N1 = MAX( 1, N )
  566. RMAGN( 2 ) = SAFMAX*ULP / DBLE( N1 )
  567. RMAGN( 3 ) = SAFMIN*ULPINV*N1
  568. *
  569. IF( NSIZES.NE.1 ) THEN
  570. MTYPES = MIN( MAXTYP, NTYPES )
  571. ELSE
  572. MTYPES = MIN( MAXTYP+1, NTYPES )
  573. END IF
  574. *
  575. DO 210 JTYPE = 1, MTYPES
  576. IF( .NOT.DOTYPE( JTYPE ) )
  577. $ GO TO 210
  578. NMATS = NMATS + 1
  579. *
  580. * Save ISEED in case of an error.
  581. *
  582. DO 20 J = 1, 4
  583. IOLDSD( J ) = ISEED( J )
  584. 20 CONTINUE
  585. *
  586. * Generate test matrices A and B
  587. *
  588. * Description of control parameters:
  589. *
  590. * KZLASS: =1 means w/o rotation, =2 means w/ rotation,
  591. * =3 means random, =4 means random generalized
  592. * upper Hessenberg.
  593. * KATYPE: the "type" to be passed to DLATM4 for computing A.
  594. * KAZERO: the pattern of zeros on the diagonal for A:
  595. * =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ),
  596. * =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ),
  597. * =6: ( 0, 1, 0, xxx, 0 ). (xxx means a string of
  598. * non-zero entries.)
  599. * KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1),
  600. * =2: large, =3: small.
  601. * IASIGN: 1 if the diagonal elements of A are to be
  602. * multiplied by a random magnitude 1 number, =2 if
  603. * randomly chosen diagonal blocks are to be rotated
  604. * to form 2x2 blocks.
  605. * KBTYPE, KBZERO, KBMAGN, IBSIGN: the same, but for B.
  606. * KTRIAN: =0: don't fill in the upper triangle, =1: do.
  607. * KZ1, KZ2, KADD: used to implement KAZERO and KBZERO.
  608. * RMAGN: used to implement KAMAGN and KBMAGN.
  609. *
  610. IF( MTYPES.GT.MAXTYP )
  611. $ GO TO 100
  612. IERR = 0
  613. IF( KCLASS( JTYPE ).LT.3 ) THEN
  614. *
  615. * Generate A (w/o rotation)
  616. *
  617. IF( ABS( KATYPE( JTYPE ) ).EQ.3 ) THEN
  618. IN = 2*( ( N-1 ) / 2 ) + 1
  619. IF( IN.NE.N )
  620. $ CALL DLASET( 'Full', N, N, ZERO, ZERO, A, LDA )
  621. ELSE
  622. IN = N
  623. END IF
  624. CALL DLATM4( KATYPE( JTYPE ), IN, KZ1( KAZERO( JTYPE ) ),
  625. $ KZ2( KAZERO( JTYPE ) ), IASIGN( JTYPE ),
  626. $ RMAGN( KAMAGN( JTYPE ) ), ULP,
  627. $ RMAGN( KTRIAN( JTYPE )*KAMAGN( JTYPE ) ), 2,
  628. $ ISEED, A, LDA )
  629. IADD = KADD( KAZERO( JTYPE ) )
  630. IF( IADD.GT.0 .AND. IADD.LE.N )
  631. $ A( IADD, IADD ) = ONE
  632. *
  633. * Generate B (w/o rotation)
  634. *
  635. IF( ABS( KBTYPE( JTYPE ) ).EQ.3 ) THEN
  636. IN = 2*( ( N-1 ) / 2 ) + 1
  637. IF( IN.NE.N )
  638. $ CALL DLASET( 'Full', N, N, ZERO, ZERO, B, LDA )
  639. ELSE
  640. IN = N
  641. END IF
  642. CALL DLATM4( KBTYPE( JTYPE ), IN, KZ1( KBZERO( JTYPE ) ),
  643. $ KZ2( KBZERO( JTYPE ) ), IBSIGN( JTYPE ),
  644. $ RMAGN( KBMAGN( JTYPE ) ), ONE,
  645. $ RMAGN( KTRIAN( JTYPE )*KBMAGN( JTYPE ) ), 2,
  646. $ ISEED, B, LDA )
  647. IADD = KADD( KBZERO( JTYPE ) )
  648. IF( IADD.NE.0 .AND. IADD.LE.N )
  649. $ B( IADD, IADD ) = ONE
  650. *
  651. IF( KCLASS( JTYPE ).EQ.2 .AND. N.GT.0 ) THEN
  652. *
  653. * Include rotations
  654. *
  655. * Generate Q, Z as Householder transformations times
  656. * a diagonal matrix.
  657. *
  658. DO 40 JC = 1, N - 1
  659. DO 30 JR = JC, N
  660. Q( JR, JC ) = DLARND( 3, ISEED )
  661. Z( JR, JC ) = DLARND( 3, ISEED )
  662. 30 CONTINUE
  663. CALL DLARFG( N+1-JC, Q( JC, JC ), Q( JC+1, JC ), 1,
  664. $ WORK( JC ) )
  665. WORK( 2*N+JC ) = SIGN( ONE, Q( JC, JC ) )
  666. Q( JC, JC ) = ONE
  667. CALL DLARFG( N+1-JC, Z( JC, JC ), Z( JC+1, JC ), 1,
  668. $ WORK( N+JC ) )
  669. WORK( 3*N+JC ) = SIGN( ONE, Z( JC, JC ) )
  670. Z( JC, JC ) = ONE
  671. 40 CONTINUE
  672. Q( N, N ) = ONE
  673. WORK( N ) = ZERO
  674. WORK( 3*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
  675. Z( N, N ) = ONE
  676. WORK( 2*N ) = ZERO
  677. WORK( 4*N ) = SIGN( ONE, DLARND( 2, ISEED ) )
  678. *
  679. * Apply the diagonal matrices
  680. *
  681. DO 60 JC = 1, N
  682. DO 50 JR = 1, N
  683. A( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
  684. $ A( JR, JC )
  685. B( JR, JC ) = WORK( 2*N+JR )*WORK( 3*N+JC )*
  686. $ B( JR, JC )
  687. 50 CONTINUE
  688. 60 CONTINUE
  689. CALL DORM2R( 'L', 'N', N, N, N-1, Q, LDQ, WORK, A,
  690. $ LDA, WORK( 2*N+1 ), IERR )
  691. IF( IERR.NE.0 )
  692. $ GO TO 90
  693. CALL DORM2R( 'R', 'T', N, N, N-1, Z, LDQ, WORK( N+1 ),
  694. $ A, LDA, WORK( 2*N+1 ), IERR )
  695. IF( IERR.NE.0 )
  696. $ GO TO 90
  697. CALL DORM2R( 'L', 'N', N, N, N-1, Q, LDQ, WORK, B,
  698. $ LDA, WORK( 2*N+1 ), IERR )
  699. IF( IERR.NE.0 )
  700. $ GO TO 90
  701. CALL DORM2R( 'R', 'T', N, N, N-1, Z, LDQ, WORK( N+1 ),
  702. $ B, LDA, WORK( 2*N+1 ), IERR )
  703. IF( IERR.NE.0 )
  704. $ GO TO 90
  705. END IF
  706. ELSE IF (KCLASS( JTYPE ).EQ.3) THEN
  707. *
  708. * Random matrices
  709. *
  710. DO 80 JC = 1, N
  711. DO 70 JR = 1, N
  712. A( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )*
  713. $ DLARND( 2, ISEED )
  714. B( JR, JC ) = RMAGN( KBMAGN( JTYPE ) )*
  715. $ DLARND( 2, ISEED )
  716. 70 CONTINUE
  717. 80 CONTINUE
  718. ELSE
  719. *
  720. * Random upper Hessenberg pencil with singular B
  721. *
  722. DO 81 JC = 1, N
  723. DO 71 JR = 1, min( JC + 1, N)
  724. A( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )*
  725. $ DLARND( 2, ISEED )
  726. 71 CONTINUE
  727. DO 72 JR = JC + 2, N
  728. A( JR, JC ) = ZERO
  729. 72 CONTINUE
  730. 81 CONTINUE
  731. DO 82 JC = 1, N
  732. DO 73 JR = 1, JC
  733. B( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )*
  734. $ DLARND( 2, ISEED )
  735. 73 CONTINUE
  736. DO 74 JR = JC + 1, N
  737. B( JR, JC ) = ZERO
  738. 74 CONTINUE
  739. 82 CONTINUE
  740. DO 83 JC = 1, N, 4
  741. B( JC, JC ) = ZERO
  742. 83 CONTINUE
  743. END IF
  744. *
  745. 90 CONTINUE
  746. *
  747. IF( IERR.NE.0 ) THEN
  748. WRITE( NOUNIT, FMT = 9999 )'Generator', IERR, N, JTYPE,
  749. $ IOLDSD
  750. INFO = ABS( IERR )
  751. RETURN
  752. END IF
  753. *
  754. 100 CONTINUE
  755. *
  756. DO 110 I = 1, 7
  757. RESULT( I ) = -ONE
  758. 110 CONTINUE
  759. *
  760. * Call XLAENV to set the parameters used in DLAQZ0
  761. *
  762. CALL XLAENV( 12, 10 )
  763. CALL XLAENV( 13, 12 )
  764. CALL XLAENV( 14, 13 )
  765. CALL XLAENV( 15, 2 )
  766. CALL XLAENV( 17, 10 )
  767. *
  768. * Call DGGEV3 to compute eigenvalues and eigenvectors.
  769. *
  770. CALL DLACPY( ' ', N, N, A, LDA, S, LDA )
  771. CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
  772. CALL DGGEV3( 'V', 'V', N, S, LDA, T, LDA, ALPHAR, ALPHAI,
  773. $ BETA, Q, LDQ, Z, LDQ, WORK, LWORK, IERR )
  774. IF( IERR.NE.0 .AND. IERR.NE.N+1 ) THEN
  775. RESULT( 1 ) = ULPINV
  776. WRITE( NOUNIT, FMT = 9999 )'DGGEV31', IERR, N, JTYPE,
  777. $ IOLDSD
  778. INFO = ABS( IERR )
  779. GO TO 190
  780. END IF
  781. *
  782. * Do the tests (1) and (2)
  783. *
  784. CALL DGET52( .TRUE., N, A, LDA, B, LDA, Q, LDQ, ALPHAR,
  785. $ ALPHAI, BETA, WORK, RESULT( 1 ) )
  786. IF( RESULT( 2 ).GT.THRESH ) THEN
  787. WRITE( NOUNIT, FMT = 9998 )'Left', 'DGGEV31',
  788. $ RESULT( 2 ), N, JTYPE, IOLDSD
  789. END IF
  790. *
  791. * Do the tests (3) and (4)
  792. *
  793. CALL DGET52( .FALSE., N, A, LDA, B, LDA, Z, LDQ, ALPHAR,
  794. $ ALPHAI, BETA, WORK, RESULT( 3 ) )
  795. IF( RESULT( 4 ).GT.THRESH ) THEN
  796. WRITE( NOUNIT, FMT = 9998 )'Right', 'DGGEV31',
  797. $ RESULT( 4 ), N, JTYPE, IOLDSD
  798. END IF
  799. *
  800. * Do the test (5)
  801. *
  802. CALL DLACPY( ' ', N, N, A, LDA, S, LDA )
  803. CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
  804. CALL DGGEV3( 'N', 'N', N, S, LDA, T, LDA, ALPHR1, ALPHI1,
  805. $ BETA1, Q, LDQ, Z, LDQ, WORK, LWORK, IERR )
  806. IF( IERR.NE.0 .AND. IERR.NE.N+1 ) THEN
  807. RESULT( 1 ) = ULPINV
  808. WRITE( NOUNIT, FMT = 9999 )'DGGEV32', IERR, N, JTYPE,
  809. $ IOLDSD
  810. INFO = ABS( IERR )
  811. GO TO 190
  812. END IF
  813. *
  814. DO 120 J = 1, N
  815. IF( ALPHAR( J ).NE.ALPHR1( J ) .OR. ALPHAI( J ).NE.
  816. $ ALPHI1( J ) .OR. BETA( J ).NE.BETA1( J ) )RESULT( 5 )
  817. $ = ULPINV
  818. 120 CONTINUE
  819. *
  820. * Do the test (6): Compute eigenvalues and left eigenvectors,
  821. * and test them
  822. *
  823. CALL DLACPY( ' ', N, N, A, LDA, S, LDA )
  824. CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
  825. CALL DGGEV3( 'V', 'N', N, S, LDA, T, LDA, ALPHR1, ALPHI1,
  826. $ BETA1, QE, LDQE, Z, LDQ, WORK, LWORK, IERR )
  827. IF( IERR.NE.0 .AND. IERR.NE.N+1 ) THEN
  828. RESULT( 1 ) = ULPINV
  829. WRITE( NOUNIT, FMT = 9999 )'DGGEV33', IERR, N, JTYPE,
  830. $ IOLDSD
  831. INFO = ABS( IERR )
  832. GO TO 190
  833. END IF
  834. *
  835. DO 130 J = 1, N
  836. IF( ALPHAR( J ).NE.ALPHR1( J ) .OR. ALPHAI( J ).NE.
  837. $ ALPHI1( J ) .OR. BETA( J ).NE.BETA1( J ) )RESULT( 6 )
  838. $ = ULPINV
  839. 130 CONTINUE
  840. *
  841. DO 150 J = 1, N
  842. DO 140 JC = 1, N
  843. IF( Q( J, JC ).NE.QE( J, JC ) )
  844. $ RESULT( 6 ) = ULPINV
  845. 140 CONTINUE
  846. 150 CONTINUE
  847. *
  848. * DO the test (7): Compute eigenvalues and right eigenvectors,
  849. * and test them
  850. *
  851. CALL DLACPY( ' ', N, N, A, LDA, S, LDA )
  852. CALL DLACPY( ' ', N, N, B, LDA, T, LDA )
  853. CALL DGGEV3( 'N', 'V', N, S, LDA, T, LDA, ALPHR1, ALPHI1,
  854. $ BETA1, Q, LDQ, QE, LDQE, WORK, LWORK, IERR )
  855. IF( IERR.NE.0 .AND. IERR.NE.N+1 ) THEN
  856. RESULT( 1 ) = ULPINV
  857. WRITE( NOUNIT, FMT = 9999 )'DGGEV34', IERR, N, JTYPE,
  858. $ IOLDSD
  859. INFO = ABS( IERR )
  860. GO TO 190
  861. END IF
  862. *
  863. DO 160 J = 1, N
  864. IF( ALPHAR( J ).NE.ALPHR1( J ) .OR. ALPHAI( J ).NE.
  865. $ ALPHI1( J ) .OR. BETA( J ).NE.BETA1( J ) )RESULT( 7 )
  866. $ = ULPINV
  867. 160 CONTINUE
  868. *
  869. DO 180 J = 1, N
  870. DO 170 JC = 1, N
  871. IF( Z( J, JC ).NE.QE( J, JC ) )
  872. $ RESULT( 7 ) = ULPINV
  873. 170 CONTINUE
  874. 180 CONTINUE
  875. *
  876. * End of Loop -- Check for RESULT(j) > THRESH
  877. *
  878. 190 CONTINUE
  879. *
  880. NTESTT = NTESTT + 7
  881. *
  882. * Print out tests which fail.
  883. *
  884. DO 200 JR = 1, 7
  885. IF( RESULT( JR ).GE.THRESH ) THEN
  886. *
  887. * If this is the first test to fail,
  888. * print a header to the data file.
  889. *
  890. IF( NERRS.EQ.0 ) THEN
  891. WRITE( NOUNIT, FMT = 9997 )'DGV'
  892. *
  893. * Matrix types
  894. *
  895. WRITE( NOUNIT, FMT = 9996 )
  896. WRITE( NOUNIT, FMT = 9995 )
  897. WRITE( NOUNIT, FMT = 9994 )'Orthogonal'
  898. *
  899. * Tests performed
  900. *
  901. WRITE( NOUNIT, FMT = 9993 )
  902. *
  903. END IF
  904. NERRS = NERRS + 1
  905. IF( RESULT( JR ).LT.10000.0D0 ) THEN
  906. WRITE( NOUNIT, FMT = 9992 )N, JTYPE, IOLDSD, JR,
  907. $ RESULT( JR )
  908. ELSE
  909. WRITE( NOUNIT, FMT = 9991 )N, JTYPE, IOLDSD, JR,
  910. $ RESULT( JR )
  911. END IF
  912. END IF
  913. 200 CONTINUE
  914. *
  915. 210 CONTINUE
  916. 220 CONTINUE
  917. *
  918. * Summary
  919. *
  920. CALL ALASVM( 'DGV', NOUNIT, NERRS, NTESTT, 0 )
  921. *
  922. WORK( 1 ) = MAXWRK
  923. *
  924. RETURN
  925. *
  926. 9999 FORMAT( ' DDRGEV3: ', A, ' returned INFO=', I6, '.', / 3X, 'N=',
  927. $ I6, ', JTYPE=', I6, ', ISEED=(', 4( I4, ',' ), I5, ')' )
  928. *
  929. 9998 FORMAT( ' DDRGEV3: ', A, ' Eigenvectors from ', A,
  930. $ ' incorrectly normalized.', / ' Bits of error=', 0P, G10.3,
  931. $ ',', 3X, 'N=', I4, ', JTYPE=', I3, ', ISEED=(',
  932. $ 4( I4, ',' ), I5, ')' )
  933. *
  934. 9997 FORMAT( / 1X, A3, ' -- Real Generalized eigenvalue problem driver'
  935. $ )
  936. *
  937. 9996 FORMAT( ' Matrix types (see DDRGEV3 for details): ' )
  938. *
  939. 9995 FORMAT( ' Special Matrices:', 23X,
  940. $ '(J''=transposed Jordan block)',
  941. $ / ' 1=(0,0) 2=(I,0) 3=(0,I) 4=(I,I) 5=(J'',J'') ',
  942. $ '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices: ( ',
  943. $ 'D=diag(0,1,2,...) )', / ' 7=(D,I) 9=(large*D, small*I',
  944. $ ') 11=(large*I, small*D) 13=(large*D, large*I)', /
  945. $ ' 8=(I,D) 10=(small*D, large*I) 12=(small*I, large*D) ',
  946. $ ' 14=(small*D, small*I)', / ' 15=(D, reversed D)' )
  947. 9994 FORMAT( ' Matrices Rotated by Random ', A, ' Matrices U, V:',
  948. $ / ' 16=Transposed Jordan Blocks 19=geometric ',
  949. $ 'alpha, beta=0,1', / ' 17=arithm. alpha&beta ',
  950. $ ' 20=arithmetic alpha, beta=0,1', / ' 18=clustered ',
  951. $ 'alpha, beta=0,1 21=random alpha, beta=0,1',
  952. $ / ' Large & Small Matrices:', / ' 22=(large, small) ',
  953. $ '23=(small,large) 24=(small,small) 25=(large,large)',
  954. $ / ' 26=random O(1) matrices.' )
  955. *
  956. 9993 FORMAT( / ' Tests performed: ',
  957. $ / ' 1 = max | ( b A - a B )''*l | / const.,',
  958. $ / ' 2 = | |VR(i)| - 1 | / ulp,',
  959. $ / ' 3 = max | ( b A - a B )*r | / const.',
  960. $ / ' 4 = | |VL(i)| - 1 | / ulp,',
  961. $ / ' 5 = 0 if W same no matter if r or l computed,',
  962. $ / ' 6 = 0 if l same no matter if l computed,',
  963. $ / ' 7 = 0 if r same no matter if r computed,', / 1X )
  964. 9992 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
  965. $ 4( I4, ',' ), ' result ', I2, ' is', 0P, F8.2 )
  966. 9991 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
  967. $ 4( I4, ',' ), ' result ', I2, ' is', 1P, D10.3 )
  968. *
  969. * End of DDRGEV3
  970. *
  971. END