LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
ddrvvx.f
Go to the documentation of this file.
1 *> \brief \b DDRVVX
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 DDRVVX( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12 * NIUNIT, NOUNIT, A, LDA, H, WR, WI, WR1, WI1,
13 * VL, LDVL, VR, LDVR, LRE, LDLRE, RCONDV, RCNDV1,
14 * RCDVIN, RCONDE, RCNDE1, RCDEIN, SCALE, SCALE1,
15 * RESULT, WORK, NWORK, IWORK, INFO )
16 *
17 * .. Scalar Arguments ..
18 * INTEGER INFO, LDA, LDLRE, LDVL, LDVR, NIUNIT, NOUNIT,
19 * $ NSIZES, NTYPES, NWORK
20 * DOUBLE PRECISION THRESH
21 * ..
22 * .. Array Arguments ..
23 * LOGICAL DOTYPE( * )
24 * INTEGER ISEED( 4 ), IWORK( * ), NN( * )
25 * DOUBLE PRECISION A( LDA, * ), H( LDA, * ), LRE( LDLRE, * ),
26 * $ RCDEIN( * ), RCDVIN( * ), RCNDE1( * ),
27 * $ RCNDV1( * ), RCONDE( * ), RCONDV( * ),
28 * $ RESULT( 11 ), SCALE( * ), SCALE1( * ),
29 * $ VL( LDVL, * ), VR( LDVR, * ), WI( * ),
30 * $ WI1( * ), WORK( * ), WR( * ), WR1( * )
31 * ..
32 *
33 *
34 *> \par Purpose:
35 * =============
36 *>
37 *> \verbatim
38 *>
39 *> DDRVVX checks the nonsymmetric eigenvalue problem expert driver
40 *> DGEEVX.
41 *>
42 *> DDRVVX uses both test matrices generated randomly depending on
43 *> data supplied in the calling sequence, as well as on data
44 *> read from an input file and including precomputed condition
45 *> numbers to which it compares the ones it computes.
46 *>
47 *> When DDRVVX is called, a number of matrix "sizes" ("n's") and a
48 *> number of matrix "types" are specified in the calling sequence.
49 *> For each size ("n") and each type of matrix, one matrix will be
50 *> generated and used to test the nonsymmetric eigenroutines. For
51 *> each matrix, 9 tests will be performed:
52 *>
53 *> (1) | A * VR - VR * W | / ( n |A| ulp )
54 *>
55 *> Here VR is the matrix of unit right eigenvectors.
56 *> W is a block diagonal matrix, with a 1x1 block for each
57 *> real eigenvalue and a 2x2 block for each complex conjugate
58 *> pair. If eigenvalues j and j+1 are a complex conjugate pair,
59 *> so WR(j) = WR(j+1) = wr and WI(j) = - WI(j+1) = wi, then the
60 *> 2 x 2 block corresponding to the pair will be:
61 *>
62 *> ( wr wi )
63 *> ( -wi wr )
64 *>
65 *> Such a block multiplying an n x 2 matrix ( ur ui ) on the
66 *> right will be the same as multiplying ur + i*ui by wr + i*wi.
67 *>
68 *> (2) | A**H * VL - VL * W**H | / ( n |A| ulp )
69 *>
70 *> Here VL is the matrix of unit left eigenvectors, A**H is the
71 *> conjugate transpose of A, and W is as above.
72 *>
73 *> (3) | |VR(i)| - 1 | / ulp and largest component real
74 *>
75 *> VR(i) denotes the i-th column of VR.
76 *>
77 *> (4) | |VL(i)| - 1 | / ulp and largest component real
78 *>
79 *> VL(i) denotes the i-th column of VL.
80 *>
81 *> (5) W(full) = W(partial)
82 *>
83 *> W(full) denotes the eigenvalues computed when VR, VL, RCONDV
84 *> and RCONDE are also computed, and W(partial) denotes the
85 *> eigenvalues computed when only some of VR, VL, RCONDV, and
86 *> RCONDE are computed.
87 *>
88 *> (6) VR(full) = VR(partial)
89 *>
90 *> VR(full) denotes the right eigenvectors computed when VL, RCONDV
91 *> and RCONDE are computed, and VR(partial) denotes the result
92 *> when only some of VL and RCONDV are computed.
93 *>
94 *> (7) VL(full) = VL(partial)
95 *>
96 *> VL(full) denotes the left eigenvectors computed when VR, RCONDV
97 *> and RCONDE are computed, and VL(partial) denotes the result
98 *> when only some of VR and RCONDV are computed.
99 *>
100 *> (8) 0 if SCALE, ILO, IHI, ABNRM (full) =
101 *> SCALE, ILO, IHI, ABNRM (partial)
102 *> 1/ulp otherwise
103 *>
104 *> SCALE, ILO, IHI and ABNRM describe how the matrix is balanced.
105 *> (full) is when VR, VL, RCONDE and RCONDV are also computed, and
106 *> (partial) is when some are not computed.
107 *>
108 *> (9) RCONDV(full) = RCONDV(partial)
109 *>
110 *> RCONDV(full) denotes the reciprocal condition numbers of the
111 *> right eigenvectors computed when VR, VL and RCONDE are also
112 *> computed. RCONDV(partial) denotes the reciprocal condition
113 *> numbers when only some of VR, VL and RCONDE are computed.
114 *>
115 *> The "sizes" are specified by an array NN(1:NSIZES); the value of
116 *> each element NN(j) specifies one size.
117 *> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
118 *> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
119 *> Currently, the list of possible types is:
120 *>
121 *> (1) The zero matrix.
122 *> (2) The identity matrix.
123 *> (3) A (transposed) Jordan block, with 1's on the diagonal.
124 *>
125 *> (4) A diagonal matrix with evenly spaced entries
126 *> 1, ..., ULP and random signs.
127 *> (ULP = (first number larger than 1) - 1 )
128 *> (5) A diagonal matrix with geometrically spaced entries
129 *> 1, ..., ULP and random signs.
130 *> (6) A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
131 *> and random signs.
132 *>
133 *> (7) Same as (4), but multiplied by a constant near
134 *> the overflow threshold
135 *> (8) Same as (4), but multiplied by a constant near
136 *> the underflow threshold
137 *>
138 *> (9) A matrix of the form U' T U, where U is orthogonal and
139 *> T has evenly spaced entries 1, ..., ULP with random signs
140 *> on the diagonal and random O(1) entries in the upper
141 *> triangle.
142 *>
143 *> (10) A matrix of the form U' T U, where U is orthogonal and
144 *> T has geometrically spaced entries 1, ..., ULP with random
145 *> signs on the diagonal and random O(1) entries in the upper
146 *> triangle.
147 *>
148 *> (11) A matrix of the form U' T U, where U is orthogonal and
149 *> T has "clustered" entries 1, ULP,..., ULP with random
150 *> signs on the diagonal and random O(1) entries in the upper
151 *> triangle.
152 *>
153 *> (12) A matrix of the form U' T U, where U is orthogonal and
154 *> T has real or complex conjugate paired eigenvalues randomly
155 *> chosen from ( ULP, 1 ) and random O(1) entries in the upper
156 *> triangle.
157 *>
158 *> (13) A matrix of the form X' T X, where X has condition
159 *> SQRT( ULP ) and T has evenly spaced entries 1, ..., ULP
160 *> with random signs on the diagonal and random O(1) entries
161 *> in the upper triangle.
162 *>
163 *> (14) A matrix of the form X' T X, where X has condition
164 *> SQRT( ULP ) and T has geometrically spaced entries
165 *> 1, ..., ULP with random signs on the diagonal and random
166 *> O(1) entries in the upper triangle.
167 *>
168 *> (15) A matrix of the form X' T X, where X has condition
169 *> SQRT( ULP ) and T has "clustered" entries 1, ULP,..., ULP
170 *> with random signs on the diagonal and random O(1) entries
171 *> in the upper triangle.
172 *>
173 *> (16) A matrix of the form X' T X, where X has condition
174 *> SQRT( ULP ) and T has real or complex conjugate paired
175 *> eigenvalues randomly chosen from ( ULP, 1 ) and random
176 *> O(1) entries in the upper triangle.
177 *>
178 *> (17) Same as (16), but multiplied by a constant
179 *> near the overflow threshold
180 *> (18) Same as (16), but multiplied by a constant
181 *> near the underflow threshold
182 *>
183 *> (19) Nonsymmetric matrix with random entries chosen from (-1,1).
184 *> If N is at least 4, all entries in first two rows and last
185 *> row, and first column and last two columns are zero.
186 *> (20) Same as (19), but multiplied by a constant
187 *> near the overflow threshold
188 *> (21) Same as (19), but multiplied by a constant
189 *> near the underflow threshold
190 *>
191 *> In addition, an input file will be read from logical unit number
192 *> NIUNIT. The file contains matrices along with precomputed
193 *> eigenvalues and reciprocal condition numbers for the eigenvalues
194 *> and right eigenvectors. For these matrices, in addition to tests
195 *> (1) to (9) we will compute the following two tests:
196 *>
197 *> (10) |RCONDV - RCDVIN| / cond(RCONDV)
198 *>
199 *> RCONDV is the reciprocal right eigenvector condition number
200 *> computed by DGEEVX and RCDVIN (the precomputed true value)
201 *> is supplied as input. cond(RCONDV) is the condition number of
202 *> RCONDV, and takes errors in computing RCONDV into account, so
203 *> that the resulting quantity should be O(ULP). cond(RCONDV) is
204 *> essentially given by norm(A)/RCONDE.
205 *>
206 *> (11) |RCONDE - RCDEIN| / cond(RCONDE)
207 *>
208 *> RCONDE is the reciprocal eigenvalue condition number
209 *> computed by DGEEVX and RCDEIN (the precomputed true value)
210 *> is supplied as input. cond(RCONDE) is the condition number
211 *> of RCONDE, and takes errors in computing RCONDE into account,
212 *> so that the resulting quantity should be O(ULP). cond(RCONDE)
213 *> is essentially given by norm(A)/RCONDV.
214 *> \endverbatim
215 *
216 * Arguments:
217 * ==========
218 *
219 *> \param[in] NSIZES
220 *> \verbatim
221 *> NSIZES is INTEGER
222 *> The number of sizes of matrices to use. NSIZES must be at
223 *> least zero. If it is zero, no randomly generated matrices
224 *> are tested, but any test matrices read from NIUNIT will be
225 *> tested.
226 *> \endverbatim
227 *>
228 *> \param[in] NN
229 *> \verbatim
230 *> NN is INTEGER array, dimension (NSIZES)
231 *> An array containing the sizes to be used for the matrices.
232 *> Zero values will be skipped. The values must be at least
233 *> zero.
234 *> \endverbatim
235 *>
236 *> \param[in] NTYPES
237 *> \verbatim
238 *> NTYPES is INTEGER
239 *> The number of elements in DOTYPE. NTYPES must be at least
240 *> zero. If it is zero, no randomly generated test matrices
241 *> are tested, but and test matrices read from NIUNIT will be
242 *> tested. If it is MAXTYP+1 and NSIZES is 1, then an
243 *> additional type, MAXTYP+1 is defined, which is to use
244 *> whatever matrix is in A. This is only useful if
245 *> DOTYPE(1:MAXTYP) is .FALSE. and DOTYPE(MAXTYP+1) is .TRUE. .
246 *> \endverbatim
247 *>
248 *> \param[in] DOTYPE
249 *> \verbatim
250 *> DOTYPE is LOGICAL array, dimension (NTYPES)
251 *> If DOTYPE(j) is .TRUE., then for each size in NN a
252 *> matrix of that size and of type j will be generated.
253 *> If NTYPES is smaller than the maximum number of types
254 *> defined (PARAMETER MAXTYP), then types NTYPES+1 through
255 *> MAXTYP will not be generated. If NTYPES is larger
256 *> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
257 *> will be ignored.
258 *> \endverbatim
259 *>
260 *> \param[in,out] ISEED
261 *> \verbatim
262 *> ISEED is INTEGER array, dimension (4)
263 *> On entry ISEED specifies the seed of the random number
264 *> generator. The array elements should be between 0 and 4095;
265 *> if not they will be reduced mod 4096. Also, ISEED(4) must
266 *> be odd. The random number generator uses a linear
267 *> congruential sequence limited to small integers, and so
268 *> should produce machine independent random numbers. The
269 *> values of ISEED are changed on exit, and can be used in the
270 *> next call to DDRVVX to continue the same random number
271 *> sequence.
272 *> \endverbatim
273 *>
274 *> \param[in] THRESH
275 *> \verbatim
276 *> THRESH is DOUBLE PRECISION
277 *> A test will count as "failed" if the "error", computed as
278 *> described above, exceeds THRESH. Note that the error
279 *> is scaled to be O(1), so THRESH should be a reasonably
280 *> small multiple of 1, e.g., 10 or 100. In particular,
281 *> it should not depend on the precision (single vs. double)
282 *> or the size of the matrix. It must be at least zero.
283 *> \endverbatim
284 *>
285 *> \param[in] NIUNIT
286 *> \verbatim
287 *> NIUNIT is INTEGER
288 *> The FORTRAN unit number for reading in the data file of
289 *> problems to solve.
290 *> \endverbatim
291 *>
292 *> \param[in] NOUNIT
293 *> \verbatim
294 *> NOUNIT is INTEGER
295 *> The FORTRAN unit number for printing out error messages
296 *> (e.g., if a routine returns INFO not equal to 0.)
297 *> \endverbatim
298 *>
299 *> \param[out] A
300 *> \verbatim
301 *> A is DOUBLE PRECISION array, dimension
302 *> (LDA, max(NN,12))
303 *> Used to hold the matrix whose eigenvalues are to be
304 *> computed. On exit, A contains the last matrix actually used.
305 *> \endverbatim
306 *>
307 *> \param[in] LDA
308 *> \verbatim
309 *> LDA is INTEGER
310 *> The leading dimension of the arrays A and H.
311 *> LDA >= max(NN,12), since 12 is the dimension of the largest
312 *> matrix in the precomputed input file.
313 *> \endverbatim
314 *>
315 *> \param[out] H
316 *> \verbatim
317 *> H is DOUBLE PRECISION array, dimension
318 *> (LDA, max(NN,12))
319 *> Another copy of the test matrix A, modified by DGEEVX.
320 *> \endverbatim
321 *>
322 *> \param[out] WR
323 *> \verbatim
324 *> WR is DOUBLE PRECISION array, dimension (max(NN))
325 *> \endverbatim
326 *>
327 *> \param[out] WI
328 *> \verbatim
329 *> WI is DOUBLE PRECISION array, dimension (max(NN))
330 *>
331 *> The real and imaginary parts of the eigenvalues of A.
332 *> On exit, WR + WI*i are the eigenvalues of the matrix in A.
333 *> \endverbatim
334 *>
335 *> \param[out] WR1
336 *> \verbatim
337 *> WR1 is DOUBLE PRECISION array, dimension (max(NN,12))
338 *> \endverbatim
339 *>
340 *> \param[out] WI1
341 *> \verbatim
342 *> WI1 is DOUBLE PRECISION array, dimension (max(NN,12))
343 *>
344 *> Like WR, WI, these arrays contain the eigenvalues of A,
345 *> but those computed when DGEEVX only computes a partial
346 *> eigendecomposition, i.e. not the eigenvalues and left
347 *> and right eigenvectors.
348 *> \endverbatim
349 *>
350 *> \param[out] VL
351 *> \verbatim
352 *> VL is DOUBLE PRECISION array, dimension
353 *> (LDVL, max(NN,12))
354 *> VL holds the computed left eigenvectors.
355 *> \endverbatim
356 *>
357 *> \param[in] LDVL
358 *> \verbatim
359 *> LDVL is INTEGER
360 *> Leading dimension of VL. Must be at least max(1,max(NN,12)).
361 *> \endverbatim
362 *>
363 *> \param[out] VR
364 *> \verbatim
365 *> VR is DOUBLE PRECISION array, dimension
366 *> (LDVR, max(NN,12))
367 *> VR holds the computed right eigenvectors.
368 *> \endverbatim
369 *>
370 *> \param[in] LDVR
371 *> \verbatim
372 *> LDVR is INTEGER
373 *> Leading dimension of VR. Must be at least max(1,max(NN,12)).
374 *> \endverbatim
375 *>
376 *> \param[out] LRE
377 *> \verbatim
378 *> LRE is DOUBLE PRECISION array, dimension
379 *> (LDLRE, max(NN,12))
380 *> LRE holds the computed right or left eigenvectors.
381 *> \endverbatim
382 *>
383 *> \param[in] LDLRE
384 *> \verbatim
385 *> LDLRE is INTEGER
386 *> Leading dimension of LRE. Must be at least max(1,max(NN,12))
387 *> \endverbatim
388 *>
389 *> \param[out] RCONDV
390 *> \verbatim
391 *> RCONDV is DOUBLE PRECISION array, dimension (N)
392 *> RCONDV holds the computed reciprocal condition numbers
393 *> for eigenvectors.
394 *> \endverbatim
395 *>
396 *> \param[out] RCNDV1
397 *> \verbatim
398 *> RCNDV1 is DOUBLE PRECISION array, dimension (N)
399 *> RCNDV1 holds more computed reciprocal condition numbers
400 *> for eigenvectors.
401 *> \endverbatim
402 *>
403 *> \param[out] RCDVIN
404 *> \verbatim
405 *> RCDVIN is DOUBLE PRECISION array, dimension (N)
406 *> When COMP = .TRUE. RCDVIN holds the precomputed reciprocal
407 *> condition numbers for eigenvectors to be compared with
408 *> RCONDV.
409 *> \endverbatim
410 *>
411 *> \param[out] RCONDE
412 *> \verbatim
413 *> RCONDE is DOUBLE PRECISION array, dimension (N)
414 *> RCONDE holds the computed reciprocal condition numbers
415 *> for eigenvalues.
416 *> \endverbatim
417 *>
418 *> \param[out] RCNDE1
419 *> \verbatim
420 *> RCNDE1 is DOUBLE PRECISION array, dimension (N)
421 *> RCNDE1 holds more computed reciprocal condition numbers
422 *> for eigenvalues.
423 *> \endverbatim
424 *>
425 *> \param[out] RCDEIN
426 *> \verbatim
427 *> RCDEIN is DOUBLE PRECISION array, dimension (N)
428 *> When COMP = .TRUE. RCDEIN holds the precomputed reciprocal
429 *> condition numbers for eigenvalues to be compared with
430 *> RCONDE.
431 *> \endverbatim
432 *>
433 *> \param[out] SCALE
434 *> \verbatim
435 *> SCALE is DOUBLE PRECISION array, dimension (N)
436 *> Holds information describing balancing of matrix.
437 *> \endverbatim
438 *>
439 *> \param[out] SCALE1
440 *> \verbatim
441 *> SCALE1 is DOUBLE PRECISION array, dimension (N)
442 *> Holds information describing balancing of matrix.
443 *> \endverbatim
444 *>
445 *> \param[out] RESULT
446 *> \verbatim
447 *> RESULT is DOUBLE PRECISION array, dimension (11)
448 *> The values computed by the seven tests described above.
449 *> The values are currently limited to 1/ulp, to avoid overflow.
450 *> \endverbatim
451 *>
452 *> \param[out] WORK
453 *> \verbatim
454 *> WORK is DOUBLE PRECISION array, dimension (NWORK)
455 *> \endverbatim
456 *>
457 *> \param[in] NWORK
458 *> \verbatim
459 *> NWORK is INTEGER
460 *> The number of entries in WORK. This must be at least
461 *> max(6*12+2*12**2,6*NN(j)+2*NN(j)**2) =
462 *> max( 360 ,6*NN(j)+2*NN(j)**2) for all j.
463 *> \endverbatim
464 *>
465 *> \param[out] IWORK
466 *> \verbatim
467 *> IWORK is INTEGER array, dimension (2*max(NN,12))
468 *> \endverbatim
469 *>
470 *> \param[out] INFO
471 *> \verbatim
472 *> INFO is INTEGER
473 *> If 0, then successful exit.
474 *> If <0, then input parameter -INFO is incorrect.
475 *> If >0, DLATMR, SLATMS, SLATME or DGET23 returned an error
476 *> code, and INFO is its absolute value.
477 *>
478 *>-----------------------------------------------------------------------
479 *>
480 *> Some Local Variables and Parameters:
481 *> ---- ----- --------- --- ----------
482 *>
483 *> ZERO, ONE Real 0 and 1.
484 *> MAXTYP The number of types defined.
485 *> NMAX Largest value in NN or 12.
486 *> NERRS The number of tests which have exceeded THRESH
487 *> COND, CONDS,
488 *> IMODE Values to be passed to the matrix generators.
489 *> ANORM Norm of A; passed to matrix generators.
490 *>
491 *> OVFL, UNFL Overflow and underflow thresholds.
492 *> ULP, ULPINV Finest relative precision and its inverse.
493 *> RTULP, RTULPI Square roots of the previous 4 values.
494 *>
495 *> The following four arrays decode JTYPE:
496 *> KTYPE(j) The general type (1-10) for type "j".
497 *> KMODE(j) The MODE value to be passed to the matrix
498 *> generator for type "j".
499 *> KMAGN(j) The order of magnitude ( O(1),
500 *> O(overflow^(1/2) ), O(underflow^(1/2) )
501 *> KCONDS(j) Selectw whether CONDS is to be 1 or
502 *> 1/sqrt(ulp). (0 means irrelevant.)
503 *> \endverbatim
504 *
505 * Authors:
506 * ========
507 *
508 *> \author Univ. of Tennessee
509 *> \author Univ. of California Berkeley
510 *> \author Univ. of Colorado Denver
511 *> \author NAG Ltd.
512 *
513 *> \date June 2016
514 *
515 *> \ingroup double_eig
516 *
517 * =====================================================================
518  SUBROUTINE ddrvvx( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
519  $ niunit, nounit, a, lda, h, wr, wi, wr1, wi1,
520  $ vl, ldvl, vr, ldvr, lre, ldlre, rcondv, rcndv1,
521  $ rcdvin, rconde, rcnde1, rcdein, scale, scale1,
522  $ result, work, nwork, iwork, info )
523 *
524 * -- LAPACK test routine (version 3.6.1) --
525 * -- LAPACK is a software package provided by Univ. of Tennessee, --
526 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
527 * June 2016
528 *
529 * .. Scalar Arguments ..
530  INTEGER INFO, LDA, LDLRE, LDVL, LDVR, NIUNIT, NOUNIT,
531  $ nsizes, ntypes, nwork
532  DOUBLE PRECISION THRESH
533 * ..
534 * .. Array Arguments ..
535  LOGICAL DOTYPE( * )
536  INTEGER ISEED( 4 ), IWORK( * ), NN( * )
537  DOUBLE PRECISION A( lda, * ), H( lda, * ), LRE( ldlre, * ),
538  $ rcdein( * ), rcdvin( * ), rcnde1( * ),
539  $ rcndv1( * ), rconde( * ), rcondv( * ),
540  $ result( 11 ), scale( * ), scale1( * ),
541  $ vl( ldvl, * ), vr( ldvr, * ), wi( * ),
542  $ wi1( * ), work( * ), wr( * ), wr1( * )
543 * ..
544 *
545 * =====================================================================
546 *
547 * .. Parameters ..
548  DOUBLE PRECISION ZERO, ONE
549  parameter ( zero = 0.0d0, one = 1.0d0 )
550  INTEGER MAXTYP
551  parameter ( maxtyp = 21 )
552 * ..
553 * .. Local Scalars ..
554  LOGICAL BADNN
555  CHARACTER BALANC
556  CHARACTER*3 PATH
557  INTEGER I, IBAL, IINFO, IMODE, ITYPE, IWK, J, JCOL,
558  $ jsize, jtype, mtypes, n, nerrs, nfail, nmax,
559  $ nnwork, ntest, ntestf, ntestt
560  DOUBLE PRECISION ANORM, COND, CONDS, OVFL, RTULP, RTULPI, ULP,
561  $ ulpinv, unfl
562 * ..
563 * .. Local Arrays ..
564  CHARACTER ADUMMA( 1 ), BAL( 4 )
565  INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KCONDS( maxtyp ),
566  $ kmagn( maxtyp ), kmode( maxtyp ),
567  $ ktype( maxtyp )
568 * ..
569 * .. External Functions ..
570  DOUBLE PRECISION DLAMCH
571  EXTERNAL dlamch
572 * ..
573 * .. External Subroutines ..
574  EXTERNAL dget23, dlabad, dlaset, dlasum, dlatme, dlatmr,
575  $ dlatms, xerbla
576 * ..
577 * .. Intrinsic Functions ..
578  INTRINSIC abs, max, min, sqrt
579 * ..
580 * .. Data statements ..
581  DATA ktype / 1, 2, 3, 5*4, 4*6, 6*6, 3*9 /
582  DATA kmagn / 3*1, 1, 1, 1, 2, 3, 4*1, 1, 1, 1, 1, 2,
583  $ 3, 1, 2, 3 /
584  DATA kmode / 3*0, 4, 3, 1, 4, 4, 4, 3, 1, 5, 4, 3,
585  $ 1, 5, 5, 5, 4, 3, 1 /
586  DATA kconds / 3*0, 5*0, 4*1, 6*2, 3*0 /
587  DATA bal / 'N', 'P', 'S', 'B' /
588 * ..
589 * .. Executable Statements ..
590 *
591  path( 1: 1 ) = 'Double precision'
592  path( 2: 3 ) = 'VX'
593 *
594 * Check for errors
595 *
596  ntestt = 0
597  ntestf = 0
598  info = 0
599 *
600 * Important constants
601 *
602  badnn = .false.
603 *
604 * 12 is the largest dimension in the input file of precomputed
605 * problems
606 *
607  nmax = 12
608  DO 10 j = 1, nsizes
609  nmax = max( nmax, nn( j ) )
610  IF( nn( j ).LT.0 )
611  $ badnn = .true.
612  10 CONTINUE
613 *
614 * Check for errors
615 *
616  IF( nsizes.LT.0 ) THEN
617  info = -1
618  ELSE IF( badnn ) THEN
619  info = -2
620  ELSE IF( ntypes.LT.0 ) THEN
621  info = -3
622  ELSE IF( thresh.LT.zero ) THEN
623  info = -6
624  ELSE IF( lda.LT.1 .OR. lda.LT.nmax ) THEN
625  info = -10
626  ELSE IF( ldvl.LT.1 .OR. ldvl.LT.nmax ) THEN
627  info = -17
628  ELSE IF( ldvr.LT.1 .OR. ldvr.LT.nmax ) THEN
629  info = -19
630  ELSE IF( ldlre.LT.1 .OR. ldlre.LT.nmax ) THEN
631  info = -21
632  ELSE IF( 6*nmax+2*nmax**2.GT.nwork ) THEN
633  info = -32
634  END IF
635 *
636  IF( info.NE.0 ) THEN
637  CALL xerbla( 'DDRVVX', -info )
638  RETURN
639  END IF
640 *
641 * If nothing to do check on NIUNIT
642 *
643  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
644  $ GO TO 160
645 *
646 * More Important constants
647 *
648  unfl = dlamch( 'Safe minimum' )
649  ovfl = one / unfl
650  CALL dlabad( unfl, ovfl )
651  ulp = dlamch( 'Precision' )
652  ulpinv = one / ulp
653  rtulp = sqrt( ulp )
654  rtulpi = one / rtulp
655 *
656 * Loop over sizes, types
657 *
658  nerrs = 0
659 *
660  DO 150 jsize = 1, nsizes
661  n = nn( jsize )
662  IF( nsizes.NE.1 ) THEN
663  mtypes = min( maxtyp, ntypes )
664  ELSE
665  mtypes = min( maxtyp+1, ntypes )
666  END IF
667 *
668  DO 140 jtype = 1, mtypes
669  IF( .NOT.dotype( jtype ) )
670  $ GO TO 140
671 *
672 * Save ISEED in case of an error.
673 *
674  DO 20 j = 1, 4
675  ioldsd( j ) = iseed( j )
676  20 CONTINUE
677 *
678 * Compute "A"
679 *
680 * Control parameters:
681 *
682 * KMAGN KCONDS KMODE KTYPE
683 * =1 O(1) 1 clustered 1 zero
684 * =2 large large clustered 2 identity
685 * =3 small exponential Jordan
686 * =4 arithmetic diagonal, (w/ eigenvalues)
687 * =5 random log symmetric, w/ eigenvalues
688 * =6 random general, w/ eigenvalues
689 * =7 random diagonal
690 * =8 random symmetric
691 * =9 random general
692 * =10 random triangular
693 *
694  IF( mtypes.GT.maxtyp )
695  $ GO TO 90
696 *
697  itype = ktype( jtype )
698  imode = kmode( jtype )
699 *
700 * Compute norm
701 *
702  GO TO ( 30, 40, 50 )kmagn( jtype )
703 *
704  30 CONTINUE
705  anorm = one
706  GO TO 60
707 *
708  40 CONTINUE
709  anorm = ovfl*ulp
710  GO TO 60
711 *
712  50 CONTINUE
713  anorm = unfl*ulpinv
714  GO TO 60
715 *
716  60 CONTINUE
717 *
718  CALL dlaset( 'Full', lda, n, zero, zero, a, lda )
719  iinfo = 0
720  cond = ulpinv
721 *
722 * Special Matrices -- Identity & Jordan block
723 *
724 * Zero
725 *
726  IF( itype.EQ.1 ) THEN
727  iinfo = 0
728 *
729  ELSE IF( itype.EQ.2 ) THEN
730 *
731 * Identity
732 *
733  DO 70 jcol = 1, n
734  a( jcol, jcol ) = anorm
735  70 CONTINUE
736 *
737  ELSE IF( itype.EQ.3 ) THEN
738 *
739 * Jordan Block
740 *
741  DO 80 jcol = 1, n
742  a( jcol, jcol ) = anorm
743  IF( jcol.GT.1 )
744  $ a( jcol, jcol-1 ) = one
745  80 CONTINUE
746 *
747  ELSE IF( itype.EQ.4 ) THEN
748 *
749 * Diagonal Matrix, [Eigen]values Specified
750 *
751  CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
752  $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
753  $ iinfo )
754 *
755  ELSE IF( itype.EQ.5 ) THEN
756 *
757 * Symmetric, eigenvalues specified
758 *
759  CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
760  $ anorm, n, n, 'N', a, lda, work( n+1 ),
761  $ iinfo )
762 *
763  ELSE IF( itype.EQ.6 ) THEN
764 *
765 * General, eigenvalues specified
766 *
767  IF( kconds( jtype ).EQ.1 ) THEN
768  conds = one
769  ELSE IF( kconds( jtype ).EQ.2 ) THEN
770  conds = rtulpi
771  ELSE
772  conds = zero
773  END IF
774 *
775  adumma( 1 ) = ' '
776  CALL dlatme( n, 'S', iseed, work, imode, cond, one,
777  $ adumma, 'T', 'T', 'T', work( n+1 ), 4,
778  $ conds, n, n, anorm, a, lda, work( 2*n+1 ),
779  $ iinfo )
780 *
781  ELSE IF( itype.EQ.7 ) THEN
782 *
783 * Diagonal, random eigenvalues
784 *
785  CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
786  $ 'T', 'N', work( n+1 ), 1, one,
787  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
788  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
789 *
790  ELSE IF( itype.EQ.8 ) THEN
791 *
792 * Symmetric, random eigenvalues
793 *
794  CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
795  $ 'T', 'N', work( n+1 ), 1, one,
796  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
797  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
798 *
799  ELSE IF( itype.EQ.9 ) THEN
800 *
801 * General, random eigenvalues
802 *
803  CALL dlatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
804  $ 'T', 'N', work( n+1 ), 1, one,
805  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
806  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
807  IF( n.GE.4 ) THEN
808  CALL dlaset( 'Full', 2, n, zero, zero, a, lda )
809  CALL dlaset( 'Full', n-3, 1, zero, zero, a( 3, 1 ),
810  $ lda )
811  CALL dlaset( 'Full', n-3, 2, zero, zero, a( 3, n-1 ),
812  $ lda )
813  CALL dlaset( 'Full', 1, n, zero, zero, a( n, 1 ),
814  $ lda )
815  END IF
816 *
817  ELSE IF( itype.EQ.10 ) THEN
818 *
819 * Triangular, random eigenvalues
820 *
821  CALL dlatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
822  $ 'T', 'N', work( n+1 ), 1, one,
823  $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
824  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
825 *
826  ELSE
827 *
828  iinfo = 1
829  END IF
830 *
831  IF( iinfo.NE.0 ) THEN
832  WRITE( nounit, fmt = 9992 )'Generator', iinfo, n, jtype,
833  $ ioldsd
834  info = abs( iinfo )
835  RETURN
836  END IF
837 *
838  90 CONTINUE
839 *
840 * Test for minimal and generous workspace
841 *
842  DO 130 iwk = 1, 3
843  IF( iwk.EQ.1 ) THEN
844  nnwork = 3*n
845  ELSE IF( iwk.EQ.2 ) THEN
846  nnwork = 6*n + n**2
847  ELSE
848  nnwork = 6*n + 2*n**2
849  END IF
850  nnwork = max( nnwork, 1 )
851 *
852 * Test for all balancing options
853 *
854  DO 120 ibal = 1, 4
855  balanc = bal( ibal )
856 *
857 * Perform tests
858 *
859  CALL dget23( .false., balanc, jtype, thresh, ioldsd,
860  $ nounit, n, a, lda, h, wr, wi, wr1, wi1,
861  $ vl, ldvl, vr, ldvr, lre, ldlre, rcondv,
862  $ rcndv1, rcdvin, rconde, rcnde1, rcdein,
863  $ scale, scale1, result, work, nnwork,
864  $ iwork, info )
865 *
866 * Check for RESULT(j) > THRESH
867 *
868  ntest = 0
869  nfail = 0
870  DO 100 j = 1, 9
871  IF( result( j ).GE.zero )
872  $ ntest = ntest + 1
873  IF( result( j ).GE.thresh )
874  $ nfail = nfail + 1
875  100 CONTINUE
876 *
877  IF( nfail.GT.0 )
878  $ ntestf = ntestf + 1
879  IF( ntestf.EQ.1 ) THEN
880  WRITE( nounit, fmt = 9999 )path
881  WRITE( nounit, fmt = 9998 )
882  WRITE( nounit, fmt = 9997 )
883  WRITE( nounit, fmt = 9996 )
884  WRITE( nounit, fmt = 9995 )thresh
885  ntestf = 2
886  END IF
887 *
888  DO 110 j = 1, 9
889  IF( result( j ).GE.thresh ) THEN
890  WRITE( nounit, fmt = 9994 )balanc, n, iwk,
891  $ ioldsd, jtype, j, result( j )
892  END IF
893  110 CONTINUE
894 *
895  nerrs = nerrs + nfail
896  ntestt = ntestt + ntest
897 *
898  120 CONTINUE
899  130 CONTINUE
900  140 CONTINUE
901  150 CONTINUE
902 *
903  160 CONTINUE
904 *
905 * Read in data from file to check accuracy of condition estimation.
906 * Assume input eigenvalues are sorted lexicographically (increasing
907 * by real part, then decreasing by imaginary part)
908 *
909  jtype = 0
910  170 CONTINUE
911  READ( niunit, fmt = *, end = 220 )n
912 *
913 * Read input data until N=0
914 *
915  IF( n.EQ.0 )
916  $ GO TO 220
917  jtype = jtype + 1
918  iseed( 1 ) = jtype
919  DO 180 i = 1, n
920  READ( niunit, fmt = * )( a( i, j ), j = 1, n )
921  180 CONTINUE
922  DO 190 i = 1, n
923  READ( niunit, fmt = * )wr1( i ), wi1( i ), rcdein( i ),
924  $ rcdvin( i )
925  190 CONTINUE
926  CALL dget23( .true., 'N', 22, thresh, iseed, nounit, n, a, lda, h,
927  $ wr, wi, wr1, wi1, vl, ldvl, vr, ldvr, lre, ldlre,
928  $ rcondv, rcndv1, rcdvin, rconde, rcnde1, rcdein,
929  $ scale, scale1, result, work, 6*n+2*n**2, iwork,
930  $ info )
931 *
932 * Check for RESULT(j) > THRESH
933 *
934  ntest = 0
935  nfail = 0
936  DO 200 j = 1, 11
937  IF( result( j ).GE.zero )
938  $ ntest = ntest + 1
939  IF( result( j ).GE.thresh )
940  $ nfail = nfail + 1
941  200 CONTINUE
942 *
943  IF( nfail.GT.0 )
944  $ ntestf = ntestf + 1
945  IF( ntestf.EQ.1 ) THEN
946  WRITE( nounit, fmt = 9999 )path
947  WRITE( nounit, fmt = 9998 )
948  WRITE( nounit, fmt = 9997 )
949  WRITE( nounit, fmt = 9996 )
950  WRITE( nounit, fmt = 9995 )thresh
951  ntestf = 2
952  END IF
953 *
954  DO 210 j = 1, 11
955  IF( result( j ).GE.thresh ) THEN
956  WRITE( nounit, fmt = 9993 )n, jtype, j, result( j )
957  END IF
958  210 CONTINUE
959 *
960  nerrs = nerrs + nfail
961  ntestt = ntestt + ntest
962  GO TO 170
963  220 CONTINUE
964 *
965 * Summary
966 *
967  CALL dlasum( path, nounit, nerrs, ntestt )
968 *
969  9999 FORMAT( / 1x, a3, ' -- Real Eigenvalue-Eigenvector Decomposition',
970  $ ' Expert Driver', /
971  $ ' Matrix types (see DDRVVX for details): ' )
972 *
973  9998 FORMAT( / ' Special Matrices:', / ' 1=Zero matrix. ',
974  $ ' ', ' 5=Diagonal: geometr. spaced entries.',
975  $ / ' 2=Identity matrix. ', ' 6=Diagona',
976  $ 'l: clustered entries.', / ' 3=Transposed Jordan block. ',
977  $ ' ', ' 7=Diagonal: large, evenly spaced.', / ' ',
978  $ '4=Diagonal: evenly spaced entries. ', ' 8=Diagonal: s',
979  $ 'mall, evenly spaced.' )
980  9997 FORMAT( ' Dense, Non-Symmetric Matrices:', / ' 9=Well-cond., ev',
981  $ 'enly spaced eigenvals.', ' 14=Ill-cond., geomet. spaced e',
982  $ 'igenals.', / ' 10=Well-cond., geom. spaced eigenvals. ',
983  $ ' 15=Ill-conditioned, clustered e.vals.', / ' 11=Well-cond',
984  $ 'itioned, clustered e.vals. ', ' 16=Ill-cond., random comp',
985  $ 'lex ', / ' 12=Well-cond., random complex ', ' ',
986  $ ' 17=Ill-cond., large rand. complx ', / ' 13=Ill-condi',
987  $ 'tioned, evenly spaced. ', ' 18=Ill-cond., small rand.',
988  $ ' complx ' )
989  9996 FORMAT( ' 19=Matrix with random O(1) entries. ', ' 21=Matrix ',
990  $ 'with small random entries.', / ' 20=Matrix with large ran',
991  $ 'dom entries. ', ' 22=Matrix read from input file', / )
992  9995 FORMAT( ' Tests performed with test threshold =', f8.2,
993  $ / / ' 1 = | A VR - VR W | / ( n |A| ulp ) ',
994  $ / ' 2 = | transpose(A) VL - VL W | / ( n |A| ulp ) ',
995  $ / ' 3 = | |VR(i)| - 1 | / ulp ',
996  $ / ' 4 = | |VL(i)| - 1 | / ulp ',
997  $ / ' 5 = 0 if W same no matter if VR or VL computed,',
998  $ ' 1/ulp otherwise', /
999  $ ' 6 = 0 if VR same no matter what else computed,',
1000  $ ' 1/ulp otherwise', /
1001  $ ' 7 = 0 if VL same no matter what else computed,',
1002  $ ' 1/ulp otherwise', /
1003  $ ' 8 = 0 if RCONDV same no matter what else computed,',
1004  $ ' 1/ulp otherwise', /
1005  $ ' 9 = 0 if SCALE, ILO, IHI, ABNRM same no matter what else',
1006  $ ' computed, 1/ulp otherwise',
1007  $ / ' 10 = | RCONDV - RCONDV(precomputed) | / cond(RCONDV),',
1008  $ / ' 11 = | RCONDE - RCONDE(precomputed) | / cond(RCONDE),' )
1009  9994 FORMAT( ' BALANC=''', a1, ''',N=', i4, ',IWK=', i1, ', seed=',
1010  $ 4( i4, ',' ), ' type ', i2, ', test(', i2, ')=', g10.3 )
1011  9993 FORMAT( ' N=', i5, ', input example =', i3, ', test(', i2, ')=',
1012  $ g10.3 )
1013  9992 FORMAT( ' DDRVVX: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
1014  $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1015 *
1016  RETURN
1017 *
1018 * End of DDRVVX
1019 *
1020  END
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine dlatme(N, DIST, ISEED, D, MODE, COND, DMAX, EI, RSIGN, UPPER, SIM, DS, MODES, CONDS, KL, KU, ANORM, A, LDA, WORK, INFO)
DLATME
Definition: dlatme.f:334
subroutine dlatmr(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, RSIGN, GRADE, DL, MODEL, CONDL, DR, MODER, CONDR, PIVTNG, IPIVOT, KL, KU, SPARSE, ANORM, PACK, A, LDA, IWORK, INFO)
DLATMR
Definition: dlatmr.f:473
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine dget23(COMP, BALANC, JTYPE, THRESH, ISEED, NOUNIT, N, A, LDA, H, WR, WI, WR1, WI1, VL, LDVL, VR, LDVR, LRE, LDLRE, RCONDV, RCNDV1, RCDVIN, RCONDE, RCNDE1, RCDEIN, SCALE, SCALE1, RESULT, WORK, LWORK, IWORK, INFO)
DGET23
Definition: dget23.f:380
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:45
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine ddrvvx(NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH, NIUNIT, NOUNIT, A, LDA, H, WR, WI, WR1, WI1, VL, LDVL, VR, LDVR, LRE, LDLRE, RCONDV, RCNDV1, RCDVIN, RCONDE, RCNDE1, RCDEIN, SCALE, SCALE1, RESULT, WORK, NWORK, IWORK, INFO)
DDRVVX
Definition: ddrvvx.f:523