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