LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
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*> \ingroup double_eig
514*
515* =====================================================================
516 SUBROUTINE ddrvvx( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
517 $ NIUNIT, NOUNIT, A, LDA, H, WR, WI, WR1, WI1,
518 $ VL, LDVL, VR, LDVR, LRE, LDLRE, RCONDV, RCNDV1,
519 $ RCDVIN, RCONDE, RCNDE1, RCDEIN, SCALE, SCALE1,
520 $ RESULT, WORK, NWORK, IWORK, INFO )
521*
522* -- LAPACK test routine --
523* -- LAPACK is a software package provided by Univ. of Tennessee, --
524* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
525*
526* .. Scalar Arguments ..
527 INTEGER INFO, LDA, LDLRE, LDVL, LDVR, NIUNIT, NOUNIT,
528 $ NSIZES, NTYPES, NWORK
529 DOUBLE PRECISION THRESH
530* ..
531* .. Array Arguments ..
532 LOGICAL DOTYPE( * )
533 INTEGER ISEED( 4 ), IWORK( * ), NN( * )
534 DOUBLE PRECISION A( LDA, * ), H( LDA, * ), LRE( LDLRE, * ),
535 $ RCDEIN( * ), RCDVIN( * ), RCNDE1( * ),
536 $ rcndv1( * ), rconde( * ), rcondv( * ),
537 $ result( 11 ), scale( * ), scale1( * ),
538 $ vl( ldvl, * ), vr( ldvr, * ), wi( * ),
539 $ wi1( * ), work( * ), wr( * ), wr1( * )
540* ..
541*
542* =====================================================================
543*
544* .. Parameters ..
545 DOUBLE PRECISION ZERO, ONE
546 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
547 INTEGER MAXTYP
548 PARAMETER ( MAXTYP = 21 )
549* ..
550* .. Local Scalars ..
551 LOGICAL BADNN
552 CHARACTER BALANC
553 CHARACTER*3 PATH
554 INTEGER I, IBAL, IINFO, IMODE, ITYPE, IWK, J, JCOL,
555 $ jsize, jtype, mtypes, n, nerrs, nfail, nmax,
556 $ nnwork, ntest, ntestf, ntestt
557 DOUBLE PRECISION ANORM, COND, CONDS, OVFL, RTULP, RTULPI, ULP,
558 $ ulpinv, unfl
559* ..
560* .. Local Arrays ..
561 CHARACTER ADUMMA( 1 ), BAL( 4 )
562 INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KCONDS( MAXTYP ),
563 $ KMAGN( MAXTYP ), KMODE( MAXTYP ),
564 $ KTYPE( MAXTYP )
565* ..
566* .. External Functions ..
567 DOUBLE PRECISION DLAMCH
568 EXTERNAL DLAMCH
569* ..
570* .. External Subroutines ..
571 EXTERNAL dget23, dlabad, dlaset, dlasum, dlatme, dlatmr,
572 $ dlatms, xerbla
573* ..
574* .. Intrinsic Functions ..
575 INTRINSIC abs, max, min, sqrt
576* ..
577* .. Data statements ..
578 DATA ktype / 1, 2, 3, 5*4, 4*6, 6*6, 3*9 /
579 DATA kmagn / 3*1, 1, 1, 1, 2, 3, 4*1, 1, 1, 1, 1, 2,
580 $ 3, 1, 2, 3 /
581 DATA kmode / 3*0, 4, 3, 1, 4, 4, 4, 3, 1, 5, 4, 3,
582 $ 1, 5, 5, 5, 4, 3, 1 /
583 DATA kconds / 3*0, 5*0, 4*1, 6*2, 3*0 /
584 DATA bal / 'N', 'P', 'S', 'B' /
585* ..
586* .. Executable Statements ..
587*
588 path( 1: 1 ) = 'Double precision'
589 path( 2: 3 ) = 'VX'
590*
591* Check for errors
592*
593 ntestt = 0
594 ntestf = 0
595 info = 0
596*
597* Important constants
598*
599 badnn = .false.
600*
601* 12 is the largest dimension in the input file of precomputed
602* problems
603*
604 nmax = 12
605 DO 10 j = 1, nsizes
606 nmax = max( nmax, nn( j ) )
607 IF( nn( j ).LT.0 )
608 $ badnn = .true.
609 10 CONTINUE
610*
611* Check for errors
612*
613 IF( nsizes.LT.0 ) THEN
614 info = -1
615 ELSE IF( badnn ) THEN
616 info = -2
617 ELSE IF( ntypes.LT.0 ) THEN
618 info = -3
619 ELSE IF( thresh.LT.zero ) THEN
620 info = -6
621 ELSE IF( lda.LT.1 .OR. lda.LT.nmax ) THEN
622 info = -10
623 ELSE IF( ldvl.LT.1 .OR. ldvl.LT.nmax ) THEN
624 info = -17
625 ELSE IF( ldvr.LT.1 .OR. ldvr.LT.nmax ) THEN
626 info = -19
627 ELSE IF( ldlre.LT.1 .OR. ldlre.LT.nmax ) THEN
628 info = -21
629 ELSE IF( 6*nmax+2*nmax**2.GT.nwork ) THEN
630 info = -32
631 END IF
632*
633 IF( info.NE.0 ) THEN
634 CALL xerbla( 'DDRVVX', -info )
635 RETURN
636 END IF
637*
638* If nothing to do check on NIUNIT
639*
640 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
641 $ GO TO 160
642*
643* More Important constants
644*
645 unfl = dlamch( 'Safe minimum' )
646 ovfl = one / unfl
647 CALL dlabad( unfl, ovfl )
648 ulp = dlamch( 'Precision' )
649 ulpinv = one / ulp
650 rtulp = sqrt( ulp )
651 rtulpi = one / rtulp
652*
653* Loop over sizes, types
654*
655 nerrs = 0
656*
657 DO 150 jsize = 1, nsizes
658 n = nn( jsize )
659 IF( nsizes.NE.1 ) THEN
660 mtypes = min( maxtyp, ntypes )
661 ELSE
662 mtypes = min( maxtyp+1, ntypes )
663 END IF
664*
665 DO 140 jtype = 1, mtypes
666 IF( .NOT.dotype( jtype ) )
667 $ GO TO 140
668*
669* Save ISEED in case of an error.
670*
671 DO 20 j = 1, 4
672 ioldsd( j ) = iseed( j )
673 20 CONTINUE
674*
675* Compute "A"
676*
677* Control parameters:
678*
679* KMAGN KCONDS KMODE KTYPE
680* =1 O(1) 1 clustered 1 zero
681* =2 large large clustered 2 identity
682* =3 small exponential Jordan
683* =4 arithmetic diagonal, (w/ eigenvalues)
684* =5 random log symmetric, w/ eigenvalues
685* =6 random general, w/ eigenvalues
686* =7 random diagonal
687* =8 random symmetric
688* =9 random general
689* =10 random triangular
690*
691 IF( mtypes.GT.maxtyp )
692 $ GO TO 90
693*
694 itype = ktype( jtype )
695 imode = kmode( jtype )
696*
697* Compute norm
698*
699 GO TO ( 30, 40, 50 )kmagn( jtype )
700*
701 30 CONTINUE
702 anorm = one
703 GO TO 60
704*
705 40 CONTINUE
706 anorm = ovfl*ulp
707 GO TO 60
708*
709 50 CONTINUE
710 anorm = unfl*ulpinv
711 GO TO 60
712*
713 60 CONTINUE
714*
715 CALL dlaset( 'Full', lda, n, zero, zero, a, lda )
716 iinfo = 0
717 cond = ulpinv
718*
719* Special Matrices -- Identity & Jordan block
720*
721* Zero
722*
723 IF( itype.EQ.1 ) THEN
724 iinfo = 0
725*
726 ELSE IF( itype.EQ.2 ) THEN
727*
728* Identity
729*
730 DO 70 jcol = 1, n
731 a( jcol, jcol ) = anorm
732 70 CONTINUE
733*
734 ELSE IF( itype.EQ.3 ) THEN
735*
736* Jordan Block
737*
738 DO 80 jcol = 1, n
739 a( jcol, jcol ) = anorm
740 IF( jcol.GT.1 )
741 $ a( jcol, jcol-1 ) = one
742 80 CONTINUE
743*
744 ELSE IF( itype.EQ.4 ) THEN
745*
746* Diagonal Matrix, [Eigen]values Specified
747*
748 CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
749 $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
750 $ iinfo )
751*
752 ELSE IF( itype.EQ.5 ) THEN
753*
754* Symmetric, eigenvalues specified
755*
756 CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
757 $ anorm, n, n, 'N', a, lda, work( n+1 ),
758 $ iinfo )
759*
760 ELSE IF( itype.EQ.6 ) THEN
761*
762* General, eigenvalues specified
763*
764 IF( kconds( jtype ).EQ.1 ) THEN
765 conds = one
766 ELSE IF( kconds( jtype ).EQ.2 ) THEN
767 conds = rtulpi
768 ELSE
769 conds = zero
770 END IF
771*
772 adumma( 1 ) = ' '
773 CALL dlatme( n, 'S', iseed, work, imode, cond, one,
774 $ adumma, 'T', 'T', 'T', work( n+1 ), 4,
775 $ conds, n, n, anorm, a, lda, work( 2*n+1 ),
776 $ iinfo )
777*
778 ELSE IF( itype.EQ.7 ) THEN
779*
780* Diagonal, random eigenvalues
781*
782 CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
783 $ 'T', 'N', work( n+1 ), 1, one,
784 $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
785 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
786*
787 ELSE IF( itype.EQ.8 ) THEN
788*
789* Symmetric, random eigenvalues
790*
791 CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
792 $ 'T', 'N', work( n+1 ), 1, one,
793 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
794 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
795*
796 ELSE IF( itype.EQ.9 ) THEN
797*
798* General, random eigenvalues
799*
800 CALL dlatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
801 $ 'T', 'N', work( n+1 ), 1, one,
802 $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
803 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
804 IF( n.GE.4 ) THEN
805 CALL dlaset( 'Full', 2, n, zero, zero, a, lda )
806 CALL dlaset( 'Full', n-3, 1, zero, zero, a( 3, 1 ),
807 $ lda )
808 CALL dlaset( 'Full', n-3, 2, zero, zero, a( 3, n-1 ),
809 $ lda )
810 CALL dlaset( 'Full', 1, n, zero, zero, a( n, 1 ),
811 $ lda )
812 END IF
813*
814 ELSE IF( itype.EQ.10 ) THEN
815*
816* Triangular, random eigenvalues
817*
818 CALL dlatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
819 $ 'T', 'N', work( n+1 ), 1, one,
820 $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
821 $ zero, anorm, 'NO', a, lda, iwork, iinfo )
822*
823 ELSE
824*
825 iinfo = 1
826 END IF
827*
828 IF( iinfo.NE.0 ) THEN
829 WRITE( nounit, fmt = 9992 )'Generator', iinfo, n, jtype,
830 $ ioldsd
831 info = abs( iinfo )
832 RETURN
833 END IF
834*
835 90 CONTINUE
836*
837* Test for minimal and generous workspace
838*
839 DO 130 iwk = 1, 3
840 IF( iwk.EQ.1 ) THEN
841 nnwork = 3*n
842 ELSE IF( iwk.EQ.2 ) THEN
843 nnwork = 6*n + n**2
844 ELSE
845 nnwork = 6*n + 2*n**2
846 END IF
847 nnwork = max( nnwork, 1 )
848*
849* Test for all balancing options
850*
851 DO 120 ibal = 1, 4
852 balanc = bal( ibal )
853*
854* Perform tests
855*
856 CALL dget23( .false., balanc, jtype, thresh, ioldsd,
857 $ nounit, n, a, lda, h, wr, wi, wr1, wi1,
858 $ vl, ldvl, vr, ldvr, lre, ldlre, rcondv,
859 $ rcndv1, rcdvin, rconde, rcnde1, rcdein,
860 $ scale, scale1, result, work, nnwork,
861 $ iwork, info )
862*
863* Check for RESULT(j) > THRESH
864*
865 ntest = 0
866 nfail = 0
867 DO 100 j = 1, 9
868 IF( result( j ).GE.zero )
869 $ ntest = ntest + 1
870 IF( result( j ).GE.thresh )
871 $ nfail = nfail + 1
872 100 CONTINUE
873*
874 IF( nfail.GT.0 )
875 $ ntestf = ntestf + 1
876 IF( ntestf.EQ.1 ) THEN
877 WRITE( nounit, fmt = 9999 )path
878 WRITE( nounit, fmt = 9998 )
879 WRITE( nounit, fmt = 9997 )
880 WRITE( nounit, fmt = 9996 )
881 WRITE( nounit, fmt = 9995 )thresh
882 ntestf = 2
883 END IF
884*
885 DO 110 j = 1, 9
886 IF( result( j ).GE.thresh ) THEN
887 WRITE( nounit, fmt = 9994 )balanc, n, iwk,
888 $ ioldsd, jtype, j, result( j )
889 END IF
890 110 CONTINUE
891*
892 nerrs = nerrs + nfail
893 ntestt = ntestt + ntest
894*
895 120 CONTINUE
896 130 CONTINUE
897 140 CONTINUE
898 150 CONTINUE
899*
900 160 CONTINUE
901*
902* Read in data from file to check accuracy of condition estimation.
903* Assume input eigenvalues are sorted lexicographically (increasing
904* by real part, then decreasing by imaginary part)
905*
906 jtype = 0
907 170 CONTINUE
908 READ( niunit, fmt = *, END = 220 )n
909*
910* Read input data until N=0
911*
912 IF( n.EQ.0 )
913 $ GO TO 220
914 jtype = jtype + 1
915 iseed( 1 ) = jtype
916 DO 180 i = 1, n
917 READ( niunit, fmt = * )( a( i, j ), j = 1, n )
918 180 CONTINUE
919 DO 190 i = 1, n
920 READ( niunit, fmt = * )wr1( i ), wi1( i ), rcdein( i ),
921 $ rcdvin( i )
922 190 CONTINUE
923 CALL dget23( .true., 'N', 22, thresh, iseed, nounit, n, a, lda, h,
924 $ wr, wi, wr1, wi1, vl, ldvl, vr, ldvr, lre, ldlre,
925 $ rcondv, rcndv1, rcdvin, rconde, rcnde1, rcdein,
926 $ scale, scale1, result, work, 6*n+2*n**2, iwork,
927 $ info )
928*
929* Check for RESULT(j) > THRESH
930*
931 ntest = 0
932 nfail = 0
933 DO 200 j = 1, 11
934 IF( result( j ).GE.zero )
935 $ ntest = ntest + 1
936 IF( result( j ).GE.thresh )
937 $ nfail = nfail + 1
938 200 CONTINUE
939*
940 IF( nfail.GT.0 )
941 $ ntestf = ntestf + 1
942 IF( ntestf.EQ.1 ) THEN
943 WRITE( nounit, fmt = 9999 )path
944 WRITE( nounit, fmt = 9998 )
945 WRITE( nounit, fmt = 9997 )
946 WRITE( nounit, fmt = 9996 )
947 WRITE( nounit, fmt = 9995 )thresh
948 ntestf = 2
949 END IF
950*
951 DO 210 j = 1, 11
952 IF( result( j ).GE.thresh ) THEN
953 WRITE( nounit, fmt = 9993 )n, jtype, j, result( j )
954 END IF
955 210 CONTINUE
956*
957 nerrs = nerrs + nfail
958 ntestt = ntestt + ntest
959 GO TO 170
960 220 CONTINUE
961*
962* Summary
963*
964 CALL dlasum( path, nounit, nerrs, ntestt )
965*
966 9999 FORMAT( / 1x, a3, ' -- Real Eigenvalue-Eigenvector Decomposition',
967 $ ' Expert Driver', /
968 $ ' Matrix types (see DDRVVX for details): ' )
969*
970 9998 FORMAT( / ' Special Matrices:', / ' 1=Zero matrix. ',
971 $ ' ', ' 5=Diagonal: geometr. spaced entries.',
972 $ / ' 2=Identity matrix. ', ' 6=Diagona',
973 $ 'l: clustered entries.', / ' 3=Transposed Jordan block. ',
974 $ ' ', ' 7=Diagonal: large, evenly spaced.', / ' ',
975 $ '4=Diagonal: evenly spaced entries. ', ' 8=Diagonal: s',
976 $ 'mall, evenly spaced.' )
977 9997 FORMAT( ' Dense, Non-Symmetric Matrices:', / ' 9=Well-cond., ev',
978 $ 'enly spaced eigenvals.', ' 14=Ill-cond., geomet. spaced e',
979 $ 'igenals.', / ' 10=Well-cond., geom. spaced eigenvals. ',
980 $ ' 15=Ill-conditioned, clustered e.vals.', / ' 11=Well-cond',
981 $ 'itioned, clustered e.vals. ', ' 16=Ill-cond., random comp',
982 $ 'lex ', / ' 12=Well-cond., random complex ', ' ',
983 $ ' 17=Ill-cond., large rand. complx ', / ' 13=Ill-condi',
984 $ 'tioned, evenly spaced. ', ' 18=Ill-cond., small rand.',
985 $ ' complx ' )
986 9996 FORMAT( ' 19=Matrix with random O(1) entries. ', ' 21=Matrix ',
987 $ 'with small random entries.', / ' 20=Matrix with large ran',
988 $ 'dom entries. ', ' 22=Matrix read from input file', / )
989 9995 FORMAT( ' Tests performed with test threshold =', f8.2,
990 $ / / ' 1 = | A VR - VR W | / ( n |A| ulp ) ',
991 $ / ' 2 = | transpose(A) VL - VL W | / ( n |A| ulp ) ',
992 $ / ' 3 = | |VR(i)| - 1 | / ulp ',
993 $ / ' 4 = | |VL(i)| - 1 | / ulp ',
994 $ / ' 5 = 0 if W same no matter if VR or VL computed,',
995 $ ' 1/ulp otherwise', /
996 $ ' 6 = 0 if VR same no matter what else computed,',
997 $ ' 1/ulp otherwise', /
998 $ ' 7 = 0 if VL same no matter what else computed,',
999 $ ' 1/ulp otherwise', /
1000 $ ' 8 = 0 if RCONDV same no matter what else computed,',
1001 $ ' 1/ulp otherwise', /
1002 $ ' 9 = 0 if SCALE, ILO, IHI, ABNRM same no matter what else',
1003 $ ' computed, 1/ulp otherwise',
1004 $ / ' 10 = | RCONDV - RCONDV(precomputed) | / cond(RCONDV),',
1005 $ / ' 11 = | RCONDE - RCONDE(precomputed) | / cond(RCONDE),' )
1006 9994 FORMAT( ' BALANC=''', a1, ''',N=', i4, ',IWK=', i1, ', seed=',
1007 $ 4( i4, ',' ), ' type ', i2, ', test(', i2, ')=', g10.3 )
1008 9993 FORMAT( ' N=', i5, ', input example =', i3, ', test(', i2, ')=',
1009 $ g10.3 )
1010 9992 FORMAT( ' DDRVVX: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
1011 $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1012*
1013 RETURN
1014*
1015* End of DDRVVX
1016*
1017 END
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
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:110
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:43
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:521
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:378
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:471
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:332
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:321