LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dchkhs.f
Go to the documentation of this file.
1 *> \brief \b DCHKHS
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 DCHKHS( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
12 * NOUNIT, A, LDA, H, T1, T2, U, LDU, Z, UZ, WR1,
13 * WI1, WR2, WI2, WR3, WI3, EVECTL, EVECTR, EVECTY,
14 * EVECTX, UU, TAU, WORK, NWORK, IWORK, SELECT,
15 * RESULT, INFO )
16 *
17 * .. Scalar Arguments ..
18 * INTEGER INFO, LDA, LDU, NOUNIT, NSIZES, NTYPES, NWORK
19 * DOUBLE PRECISION THRESH
20 * ..
21 * .. Array Arguments ..
22 * LOGICAL DOTYPE( * ), SELECT( * )
23 * INTEGER ISEED( 4 ), IWORK( * ), NN( * )
24 * DOUBLE PRECISION A( LDA, * ), EVECTL( LDU, * ),
25 * $ EVECTR( LDU, * ), EVECTX( LDU, * ),
26 * $ EVECTY( LDU, * ), H( LDA, * ), RESULT( 14 ),
27 * $ T1( LDA, * ), T2( LDA, * ), TAU( * ),
28 * $ U( LDU, * ), UU( LDU, * ), UZ( LDU, * ),
29 * $ WI1( * ), WI2( * ), WI3( * ), WORK( * ),
30 * $ WR1( * ), WR2( * ), WR3( * ), Z( LDU, * )
31 * ..
32 *
33 *
34 *> \par Purpose:
35 * =============
36 *>
37 *> \verbatim
38 *>
39 *> DCHKHS checks the nonsymmetric eigenvalue problem routines.
40 *>
41 *> DGEHRD factors A as U H U' , where ' means transpose,
42 *> H is hessenberg, and U is an orthogonal matrix.
43 *>
44 *> DORGHR generates the orthogonal matrix U.
45 *>
46 *> DORMHR multiplies a matrix by the orthogonal matrix U.
47 *>
48 *> DHSEQR factors H as Z T Z' , where Z is orthogonal and
49 *> T is "quasi-triangular", and the eigenvalue vector W.
50 *>
51 *> DTREVC computes the left and right eigenvector matrices
52 *> L and R for T.
53 *>
54 *> DHSEIN computes the left and right eigenvector matrices
55 *> Y and X for H, using inverse iteration.
56 *>
57 *> When DCHKHS is called, a number of matrix "sizes" ("n's") and a
58 *> number of matrix "types" are specified. For each size ("n")
59 *> and each type of matrix, one matrix will be generated and used
60 *> to test the nonsymmetric eigenroutines. For each matrix, 14
61 *> tests will be performed:
62 *>
63 *> (1) | A - U H U**T | / ( |A| n ulp )
64 *>
65 *> (2) | I - UU**T | / ( n ulp )
66 *>
67 *> (3) | H - Z T Z**T | / ( |H| n ulp )
68 *>
69 *> (4) | I - ZZ**T | / ( n ulp )
70 *>
71 *> (5) | A - UZ H (UZ)**T | / ( |A| n ulp )
72 *>
73 *> (6) | I - UZ (UZ)**T | / ( n ulp )
74 *>
75 *> (7) | T(Z computed) - T(Z not computed) | / ( |T| ulp )
76 *>
77 *> (8) | W(Z computed) - W(Z not computed) | / ( |W| ulp )
78 *>
79 *> (9) | TR - RW | / ( |T| |R| ulp )
80 *>
81 *> (10) | L**H T - W**H L | / ( |T| |L| ulp )
82 *>
83 *> (11) | HX - XW | / ( |H| |X| ulp )
84 *>
85 *> (12) | Y**H H - W**H Y | / ( |H| |Y| ulp )
86 *>
87 *> (13) | AX - XW | / ( |A| |X| ulp )
88 *>
89 *> (14) | Y**H A - W**H Y | / ( |A| |Y| ulp )
90 *>
91 *> The "sizes" are specified by an array NN(1:NSIZES); the value of
92 *> each element NN(j) specifies one size.
93 *> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
94 *> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
95 *> Currently, the list of possible types is:
96 *>
97 *> (1) The zero matrix.
98 *> (2) The identity matrix.
99 *> (3) A (transposed) Jordan block, with 1's on the diagonal.
100 *>
101 *> (4) A diagonal matrix with evenly spaced entries
102 *> 1, ..., ULP and random signs.
103 *> (ULP = (first number larger than 1) - 1 )
104 *> (5) A diagonal matrix with geometrically spaced entries
105 *> 1, ..., ULP and random signs.
106 *> (6) A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
107 *> and random signs.
108 *>
109 *> (7) Same as (4), but multiplied by SQRT( overflow threshold )
110 *> (8) Same as (4), but multiplied by SQRT( underflow threshold )
111 *>
112 *> (9) A matrix of the form U' T U, where U is orthogonal and
113 *> T has evenly spaced entries 1, ..., ULP with random signs
114 *> on the diagonal and random O(1) entries in the upper
115 *> triangle.
116 *>
117 *> (10) A matrix of the form U' T U, where U is orthogonal and
118 *> T has geometrically spaced entries 1, ..., ULP with random
119 *> signs on the diagonal and random O(1) entries in the upper
120 *> triangle.
121 *>
122 *> (11) A matrix of the form U' T U, where U is orthogonal and
123 *> T has "clustered" entries 1, ULP,..., ULP with random
124 *> signs on the diagonal and random O(1) entries in the upper
125 *> triangle.
126 *>
127 *> (12) A matrix of the form U' T U, where U is orthogonal and
128 *> T has real or complex conjugate paired eigenvalues randomly
129 *> chosen from ( ULP, 1 ) and random O(1) entries in the upper
130 *> triangle.
131 *>
132 *> (13) A matrix of the form X' T X, where X has condition
133 *> SQRT( ULP ) and T has evenly spaced entries 1, ..., ULP
134 *> with random signs on the diagonal and random O(1) entries
135 *> in the upper triangle.
136 *>
137 *> (14) A matrix of the form X' T X, where X has condition
138 *> SQRT( ULP ) and T has geometrically spaced entries
139 *> 1, ..., ULP with random signs on the diagonal and random
140 *> O(1) entries in the upper triangle.
141 *>
142 *> (15) A matrix of the form X' T X, where X has condition
143 *> SQRT( ULP ) and T has "clustered" entries 1, ULP,..., ULP
144 *> with random signs on the diagonal and random O(1) entries
145 *> in the upper triangle.
146 *>
147 *> (16) A matrix of the form X' T X, where X has condition
148 *> SQRT( ULP ) and T has real or complex conjugate paired
149 *> eigenvalues randomly chosen from ( ULP, 1 ) and random
150 *> O(1) entries in the upper triangle.
151 *>
152 *> (17) Same as (16), but multiplied by SQRT( overflow threshold )
153 *> (18) Same as (16), but multiplied by SQRT( underflow threshold )
154 *>
155 *> (19) Nonsymmetric matrix with random entries chosen from (-1,1).
156 *> (20) Same as (19), but multiplied by SQRT( overflow threshold )
157 *> (21) Same as (19), but multiplied by SQRT( underflow threshold )
158 *> \endverbatim
159 *
160 * Arguments:
161 * ==========
162 *
163 *> \verbatim
164 *> NSIZES - INTEGER
165 *> The number of sizes of matrices to use. If it is zero,
166 *> DCHKHS does nothing. It must be at least zero.
167 *> Not modified.
168 *>
169 *> NN - INTEGER array, dimension (NSIZES)
170 *> An array containing the sizes to be used for the matrices.
171 *> Zero values will be skipped. The values must be at least
172 *> zero.
173 *> Not modified.
174 *>
175 *> NTYPES - INTEGER
176 *> The number of elements in DOTYPE. If it is zero, DCHKHS
177 *> does nothing. It must be at least zero. If it is MAXTYP+1
178 *> and NSIZES is 1, then an additional type, MAXTYP+1 is
179 *> defined, which is to use whatever matrix is in A. This
180 *> is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
181 *> DOTYPE(MAXTYP+1) is .TRUE. .
182 *> Not modified.
183 *>
184 *> DOTYPE - LOGICAL array, dimension (NTYPES)
185 *> If DOTYPE(j) is .TRUE., then for each size in NN a
186 *> matrix of that size and of type j will be generated.
187 *> If NTYPES is smaller than the maximum number of types
188 *> defined (PARAMETER MAXTYP), then types NTYPES+1 through
189 *> MAXTYP will not be generated. If NTYPES is larger
190 *> than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
191 *> will be ignored.
192 *> Not modified.
193 *>
194 *> ISEED - INTEGER array, dimension (4)
195 *> On entry ISEED specifies the seed of the random number
196 *> generator. The array elements should be between 0 and 4095;
197 *> if not they will be reduced mod 4096. Also, ISEED(4) must
198 *> be odd. The random number generator uses a linear
199 *> congruential sequence limited to small integers, and so
200 *> should produce machine independent random numbers. The
201 *> values of ISEED are changed on exit, and can be used in the
202 *> next call to DCHKHS to continue the same random number
203 *> sequence.
204 *> Modified.
205 *>
206 *> THRESH - DOUBLE PRECISION
207 *> A test will count as "failed" if the "error", computed as
208 *> described above, exceeds THRESH. Note that the error
209 *> is scaled to be O(1), so THRESH should be a reasonably
210 *> small multiple of 1, e.g., 10 or 100. In particular,
211 *> it should not depend on the precision (single vs. double)
212 *> or the size of the matrix. It must be at least zero.
213 *> Not modified.
214 *>
215 *> NOUNIT - INTEGER
216 *> The FORTRAN unit number for printing out error messages
217 *> (e.g., if a routine returns IINFO not equal to 0.)
218 *> Not modified.
219 *>
220 *> A - DOUBLE PRECISION array, dimension (LDA,max(NN))
221 *> Used to hold the matrix whose eigenvalues are to be
222 *> computed. On exit, A contains the last matrix actually
223 *> used.
224 *> Modified.
225 *>
226 *> LDA - INTEGER
227 *> The leading dimension of A, H, T1 and T2. It must be at
228 *> least 1 and at least max( NN ).
229 *> Not modified.
230 *>
231 *> H - DOUBLE PRECISION array, dimension (LDA,max(NN))
232 *> The upper hessenberg matrix computed by DGEHRD. On exit,
233 *> H contains the Hessenberg form of the matrix in A.
234 *> Modified.
235 *>
236 *> T1 - DOUBLE PRECISION array, dimension (LDA,max(NN))
237 *> The Schur (="quasi-triangular") matrix computed by DHSEQR
238 *> if Z is computed. On exit, T1 contains the Schur form of
239 *> the matrix in A.
240 *> Modified.
241 *>
242 *> T2 - DOUBLE PRECISION array, dimension (LDA,max(NN))
243 *> The Schur matrix computed by DHSEQR when Z is not computed.
244 *> This should be identical to T1.
245 *> Modified.
246 *>
247 *> LDU - INTEGER
248 *> The leading dimension of U, Z, UZ and UU. It must be at
249 *> least 1 and at least max( NN ).
250 *> Not modified.
251 *>
252 *> U - DOUBLE PRECISION array, dimension (LDU,max(NN))
253 *> The orthogonal matrix computed by DGEHRD.
254 *> Modified.
255 *>
256 *> Z - DOUBLE PRECISION array, dimension (LDU,max(NN))
257 *> The orthogonal matrix computed by DHSEQR.
258 *> Modified.
259 *>
260 *> UZ - DOUBLE PRECISION array, dimension (LDU,max(NN))
261 *> The product of U times Z.
262 *> Modified.
263 *>
264 *> WR1 - DOUBLE PRECISION array, dimension (max(NN))
265 *> WI1 - DOUBLE PRECISION array, dimension (max(NN))
266 *> The real and imaginary parts of the eigenvalues of A,
267 *> as computed when Z is computed.
268 *> On exit, WR1 + WI1*i are the eigenvalues of the matrix in A.
269 *> Modified.
270 *>
271 *> WR2 - DOUBLE PRECISION array, dimension (max(NN))
272 *> WI2 - DOUBLE PRECISION array, dimension (max(NN))
273 *> The real and imaginary parts of the eigenvalues of A,
274 *> as computed when T is computed but not Z.
275 *> On exit, WR2 + WI2*i are the eigenvalues of the matrix in A.
276 *> Modified.
277 *>
278 *> WR3 - DOUBLE PRECISION array, dimension (max(NN))
279 *> WI3 - DOUBLE PRECISION array, dimension (max(NN))
280 *> Like WR1, WI1, these arrays contain the eigenvalues of A,
281 *> but those computed when DHSEQR only computes the
282 *> eigenvalues, i.e., not the Schur vectors and no more of the
283 *> Schur form than is necessary for computing the
284 *> eigenvalues.
285 *> Modified.
286 *>
287 *> EVECTL - DOUBLE PRECISION array, dimension (LDU,max(NN))
288 *> The (upper triangular) left eigenvector matrix for the
289 *> matrix in T1. For complex conjugate pairs, the real part
290 *> is stored in one row and the imaginary part in the next.
291 *> Modified.
292 *>
293 *> EVEZTR - DOUBLE PRECISION array, dimension (LDU,max(NN))
294 *> The (upper triangular) right eigenvector matrix for the
295 *> matrix in T1. For complex conjugate pairs, the real part
296 *> is stored in one column and the imaginary part in the next.
297 *> Modified.
298 *>
299 *> EVECTY - DOUBLE PRECISION array, dimension (LDU,max(NN))
300 *> The left eigenvector matrix for the
301 *> matrix in H. For complex conjugate pairs, the real part
302 *> is stored in one row and the imaginary part in the next.
303 *> Modified.
304 *>
305 *> EVECTX - DOUBLE PRECISION array, dimension (LDU,max(NN))
306 *> The right eigenvector matrix for the
307 *> matrix in H. For complex conjugate pairs, the real part
308 *> is stored in one column and the imaginary part in the next.
309 *> Modified.
310 *>
311 *> UU - DOUBLE PRECISION array, dimension (LDU,max(NN))
312 *> Details of the orthogonal matrix computed by DGEHRD.
313 *> Modified.
314 *>
315 *> TAU - DOUBLE PRECISION array, dimension(max(NN))
316 *> Further details of the orthogonal matrix computed by DGEHRD.
317 *> Modified.
318 *>
319 *> WORK - DOUBLE PRECISION array, dimension (NWORK)
320 *> Workspace.
321 *> Modified.
322 *>
323 *> NWORK - INTEGER
324 *> The number of entries in WORK. NWORK >= 4*NN(j)*NN(j) + 2.
325 *>
326 *> IWORK - INTEGER array, dimension (max(NN))
327 *> Workspace.
328 *> Modified.
329 *>
330 *> SELECT - LOGICAL array, dimension (max(NN))
331 *> Workspace.
332 *> Modified.
333 *>
334 *> RESULT - DOUBLE PRECISION array, dimension (14)
335 *> The values computed by the fourteen tests described above.
336 *> The values are currently limited to 1/ulp, to avoid
337 *> overflow.
338 *> Modified.
339 *>
340 *> INFO - INTEGER
341 *> If 0, then everything ran OK.
342 *> -1: NSIZES < 0
343 *> -2: Some NN(j) < 0
344 *> -3: NTYPES < 0
345 *> -6: THRESH < 0
346 *> -9: LDA < 1 or LDA < NMAX, where NMAX is max( NN(j) ).
347 *> -14: LDU < 1 or LDU < NMAX.
348 *> -28: NWORK too small.
349 *> If DLATMR, SLATMS, or SLATME returns an error code, the
350 *> absolute value of it is returned.
351 *> If 1, then DHSEQR could not find all the shifts.
352 *> If 2, then the EISPACK code (for small blocks) failed.
353 *> If >2, then 30*N iterations were not enough to find an
354 *> eigenvalue or to decompose the problem.
355 *> Modified.
356 *>
357 *>-----------------------------------------------------------------------
358 *>
359 *> Some Local Variables and Parameters:
360 *> ---- ----- --------- --- ----------
361 *>
362 *> ZERO, ONE Real 0 and 1.
363 *> MAXTYP The number of types defined.
364 *> MTEST The number of tests defined: care must be taken
365 *> that (1) the size of RESULT, (2) the number of
366 *> tests actually performed, and (3) MTEST agree.
367 *> NTEST The number of tests performed on this matrix
368 *> so far. This should be less than MTEST, and
369 *> equal to it by the last test. It will be less
370 *> if any of the routines being tested indicates
371 *> that it could not compute the matrices that
372 *> would be tested.
373 *> NMAX Largest value in NN.
374 *> NMATS The number of matrices generated so far.
375 *> NERRS The number of tests which have exceeded THRESH
376 *> so far (computed by DLAFTS).
377 *> COND, CONDS,
378 *> IMODE Values to be passed to the matrix generators.
379 *> ANORM Norm of A; passed to matrix generators.
380 *>
381 *> OVFL, UNFL Overflow and underflow thresholds.
382 *> ULP, ULPINV Finest relative precision and its inverse.
383 *> RTOVFL, RTUNFL,
384 *> RTULP, RTULPI Square roots of the previous 4 values.
385 *>
386 *> The following four arrays decode JTYPE:
387 *> KTYPE(j) The general type (1-10) for type "j".
388 *> KMODE(j) The MODE value to be passed to the matrix
389 *> generator for type "j".
390 *> KMAGN(j) The order of magnitude ( O(1),
391 *> O(overflow^(1/2) ), O(underflow^(1/2) )
392 *> KCONDS(j) Selects whether CONDS is to be 1 or
393 *> 1/sqrt(ulp). (0 means irrelevant.)
394 *> \endverbatim
395 *
396 * Authors:
397 * ========
398 *
399 *> \author Univ. of Tennessee
400 *> \author Univ. of California Berkeley
401 *> \author Univ. of Colorado Denver
402 *> \author NAG Ltd.
403 *
404 *> \date November 2015
405 *
406 *> \ingroup double_eig
407 *
408 * =====================================================================
409  SUBROUTINE dchkhs( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
410  $ nounit, a, lda, h, t1, t2, u, ldu, z, uz, wr1,
411  $ wi1, wr2, wi2, wr3, wi3, evectl, evectr,
412  $ evecty, evectx, uu, tau, work, nwork, iwork,
413  $ SELECT, result, info )
414 *
415 * -- LAPACK test routine (version 3.6.0) --
416 * -- LAPACK is a software package provided by Univ. of Tennessee, --
417 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
418 * November 2015
419 *
420 * .. Scalar Arguments ..
421  INTEGER INFO, LDA, LDU, NOUNIT, NSIZES, NTYPES, NWORK
422  DOUBLE PRECISION THRESH
423 * ..
424 * .. Array Arguments ..
425  LOGICAL DOTYPE( * ), SELECT( * )
426  INTEGER ISEED( 4 ), IWORK( * ), NN( * )
427  DOUBLE PRECISION A( lda, * ), EVECTL( ldu, * ),
428  $ evectr( ldu, * ), evectx( ldu, * ),
429  $ evecty( ldu, * ), h( lda, * ), result( 14 ),
430  $ t1( lda, * ), t2( lda, * ), tau( * ),
431  $ u( ldu, * ), uu( ldu, * ), uz( ldu, * ),
432  $ wi1( * ), wi2( * ), wi3( * ), work( * ),
433  $ wr1( * ), wr2( * ), wr3( * ), z( ldu, * )
434 * ..
435 *
436 * =====================================================================
437 *
438 * .. Parameters ..
439  DOUBLE PRECISION ZERO, ONE
440  parameter ( zero = 0.0d0, one = 1.0d0 )
441  INTEGER MAXTYP
442  parameter ( maxtyp = 21 )
443 * ..
444 * .. Local Scalars ..
445  LOGICAL BADNN, MATCH
446  INTEGER I, IHI, IINFO, ILO, IMODE, IN, ITYPE, J, JCOL,
447  $ jj, jsize, jtype, k, mtypes, n, n1, nerrs,
448  $ nmats, nmax, nselc, nselr, ntest, ntestt
449  DOUBLE PRECISION ANINV, ANORM, COND, CONDS, OVFL, RTOVFL, RTULP,
450  $ rtulpi, rtunfl, temp1, temp2, ulp, ulpinv, unfl
451 * ..
452 * .. Local Arrays ..
453  CHARACTER ADUMMA( 1 )
454  INTEGER IDUMMA( 1 ), IOLDSD( 4 ), KCONDS( maxtyp ),
455  $ kmagn( maxtyp ), kmode( maxtyp ),
456  $ ktype( maxtyp )
457  DOUBLE PRECISION DUMMA( 6 )
458 * ..
459 * .. External Functions ..
460  DOUBLE PRECISION DLAMCH
461  EXTERNAL dlamch
462 * ..
463 * .. External Subroutines ..
464  EXTERNAL dcopy, dgehrd, dgemm, dget10, dget22, dhsein,
467  $ dtrevc, xerbla
468 * ..
469 * .. Intrinsic Functions ..
470  INTRINSIC abs, dble, max, min, sqrt
471 * ..
472 * .. Data statements ..
473  DATA ktype / 1, 2, 3, 5*4, 4*6, 6*6, 3*9 /
474  DATA kmagn / 3*1, 1, 1, 1, 2, 3, 4*1, 1, 1, 1, 1, 2,
475  $ 3, 1, 2, 3 /
476  DATA kmode / 3*0, 4, 3, 1, 4, 4, 4, 3, 1, 5, 4, 3,
477  $ 1, 5, 5, 5, 4, 3, 1 /
478  DATA kconds / 3*0, 5*0, 4*1, 6*2, 3*0 /
479 * ..
480 * .. Executable Statements ..
481 *
482 * Check for errors
483 *
484  ntestt = 0
485  info = 0
486 *
487  badnn = .false.
488  nmax = 0
489  DO 10 j = 1, nsizes
490  nmax = max( nmax, nn( j ) )
491  IF( nn( j ).LT.0 )
492  $ badnn = .true.
493  10 CONTINUE
494 *
495 * Check for errors
496 *
497  IF( nsizes.LT.0 ) THEN
498  info = -1
499  ELSE IF( badnn ) THEN
500  info = -2
501  ELSE IF( ntypes.LT.0 ) THEN
502  info = -3
503  ELSE IF( thresh.LT.zero ) THEN
504  info = -6
505  ELSE IF( lda.LE.1 .OR. lda.LT.nmax ) THEN
506  info = -9
507  ELSE IF( ldu.LE.1 .OR. ldu.LT.nmax ) THEN
508  info = -14
509  ELSE IF( 4*nmax*nmax+2.GT.nwork ) THEN
510  info = -28
511  END IF
512 *
513  IF( info.NE.0 ) THEN
514  CALL xerbla( 'DCHKHS', -info )
515  RETURN
516  END IF
517 *
518 * Quick return if possible
519 *
520  IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
521  $ RETURN
522 *
523 * More important constants
524 *
525  unfl = dlamch( 'Safe minimum' )
526  ovfl = dlamch( 'Overflow' )
527  CALL dlabad( unfl, ovfl )
528  ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
529  ulpinv = one / ulp
530  rtunfl = sqrt( unfl )
531  rtovfl = sqrt( ovfl )
532  rtulp = sqrt( ulp )
533  rtulpi = one / rtulp
534 *
535 * Loop over sizes, types
536 *
537  nerrs = 0
538  nmats = 0
539 *
540  DO 270 jsize = 1, nsizes
541  n = nn( jsize )
542  IF( n.EQ.0 )
543  $ GO TO 270
544  n1 = max( 1, n )
545  aninv = one / dble( n1 )
546 *
547  IF( nsizes.NE.1 ) THEN
548  mtypes = min( maxtyp, ntypes )
549  ELSE
550  mtypes = min( maxtyp+1, ntypes )
551  END IF
552 *
553  DO 260 jtype = 1, mtypes
554  IF( .NOT.dotype( jtype ) )
555  $ GO TO 260
556  nmats = nmats + 1
557  ntest = 0
558 *
559 * Save ISEED in case of an error.
560 *
561  DO 20 j = 1, 4
562  ioldsd( j ) = iseed( j )
563  20 CONTINUE
564 *
565 * Initialize RESULT
566 *
567  DO 30 j = 1, 14
568  result( j ) = zero
569  30 CONTINUE
570 *
571 * Compute "A"
572 *
573 * Control parameters:
574 *
575 * KMAGN KCONDS KMODE KTYPE
576 * =1 O(1) 1 clustered 1 zero
577 * =2 large large clustered 2 identity
578 * =3 small exponential Jordan
579 * =4 arithmetic diagonal, (w/ eigenvalues)
580 * =5 random log symmetric, w/ eigenvalues
581 * =6 random general, w/ eigenvalues
582 * =7 random diagonal
583 * =8 random symmetric
584 * =9 random general
585 * =10 random triangular
586 *
587  IF( mtypes.GT.maxtyp )
588  $ GO TO 100
589 *
590  itype = ktype( jtype )
591  imode = kmode( jtype )
592 *
593 * Compute norm
594 *
595  GO TO ( 40, 50, 60 )kmagn( jtype )
596 *
597  40 CONTINUE
598  anorm = one
599  GO TO 70
600 *
601  50 CONTINUE
602  anorm = ( rtovfl*ulp )*aninv
603  GO TO 70
604 *
605  60 CONTINUE
606  anorm = rtunfl*n*ulpinv
607  GO TO 70
608 *
609  70 CONTINUE
610 *
611  CALL dlaset( 'Full', lda, n, zero, zero, a, lda )
612  iinfo = 0
613  cond = ulpinv
614 *
615 * Special Matrices
616 *
617  IF( itype.EQ.1 ) THEN
618 *
619 * Zero
620 *
621  iinfo = 0
622 *
623  ELSE IF( itype.EQ.2 ) THEN
624 *
625 * Identity
626 *
627  DO 80 jcol = 1, n
628  a( jcol, jcol ) = anorm
629  80 CONTINUE
630 *
631  ELSE IF( itype.EQ.3 ) THEN
632 *
633 * Jordan Block
634 *
635  DO 90 jcol = 1, n
636  a( jcol, jcol ) = anorm
637  IF( jcol.GT.1 )
638  $ a( jcol, jcol-1 ) = one
639  90 CONTINUE
640 *
641  ELSE IF( itype.EQ.4 ) THEN
642 *
643 * Diagonal Matrix, [Eigen]values Specified
644 *
645  CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
646  $ anorm, 0, 0, 'N', a, lda, work( n+1 ),
647  $ iinfo )
648 *
649  ELSE IF( itype.EQ.5 ) THEN
650 *
651 * Symmetric, eigenvalues specified
652 *
653  CALL dlatms( n, n, 'S', iseed, 'S', work, imode, cond,
654  $ anorm, n, n, 'N', a, lda, work( n+1 ),
655  $ iinfo )
656 *
657  ELSE IF( itype.EQ.6 ) THEN
658 *
659 * General, eigenvalues specified
660 *
661  IF( kconds( jtype ).EQ.1 ) THEN
662  conds = one
663  ELSE IF( kconds( jtype ).EQ.2 ) THEN
664  conds = rtulpi
665  ELSE
666  conds = zero
667  END IF
668 *
669  adumma( 1 ) = ' '
670  CALL dlatme( n, 'S', iseed, work, imode, cond, one,
671  $ adumma, 'T', 'T', 'T', work( n+1 ), 4,
672  $ conds, n, n, anorm, a, lda, work( 2*n+1 ),
673  $ iinfo )
674 *
675  ELSE IF( itype.EQ.7 ) THEN
676 *
677 * Diagonal, random eigenvalues
678 *
679  CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
680  $ 'T', 'N', work( n+1 ), 1, one,
681  $ work( 2*n+1 ), 1, one, 'N', idumma, 0, 0,
682  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
683 *
684  ELSE IF( itype.EQ.8 ) THEN
685 *
686 * Symmetric, random eigenvalues
687 *
688  CALL dlatmr( n, n, 'S', iseed, 'S', work, 6, one, one,
689  $ 'T', 'N', work( n+1 ), 1, one,
690  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
691  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
692 *
693  ELSE IF( itype.EQ.9 ) THEN
694 *
695 * General, random eigenvalues
696 *
697  CALL dlatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
698  $ 'T', 'N', work( n+1 ), 1, one,
699  $ work( 2*n+1 ), 1, one, 'N', idumma, n, n,
700  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
701 *
702  ELSE IF( itype.EQ.10 ) THEN
703 *
704 * Triangular, random eigenvalues
705 *
706  CALL dlatmr( n, n, 'S', iseed, 'N', work, 6, one, one,
707  $ 'T', 'N', work( n+1 ), 1, one,
708  $ work( 2*n+1 ), 1, one, 'N', idumma, n, 0,
709  $ zero, anorm, 'NO', a, lda, iwork, iinfo )
710 *
711  ELSE
712 *
713  iinfo = 1
714  END IF
715 *
716  IF( iinfo.NE.0 ) THEN
717  WRITE( nounit, fmt = 9999 )'Generator', iinfo, n, jtype,
718  $ ioldsd
719  info = abs( iinfo )
720  RETURN
721  END IF
722 *
723  100 CONTINUE
724 *
725 * Call DGEHRD to compute H and U, do tests.
726 *
727  CALL dlacpy( ' ', n, n, a, lda, h, lda )
728 *
729  ntest = 1
730 *
731  ilo = 1
732  ihi = n
733 *
734  CALL dgehrd( n, ilo, ihi, h, lda, work, work( n+1 ),
735  $ nwork-n, iinfo )
736 *
737  IF( iinfo.NE.0 ) THEN
738  result( 1 ) = ulpinv
739  WRITE( nounit, fmt = 9999 )'DGEHRD', iinfo, n, jtype,
740  $ ioldsd
741  info = abs( iinfo )
742  GO TO 250
743  END IF
744 *
745  DO 120 j = 1, n - 1
746  uu( j+1, j ) = zero
747  DO 110 i = j + 2, n
748  u( i, j ) = h( i, j )
749  uu( i, j ) = h( i, j )
750  h( i, j ) = zero
751  110 CONTINUE
752  120 CONTINUE
753  CALL dcopy( n-1, work, 1, tau, 1 )
754  CALL dorghr( n, ilo, ihi, u, ldu, work, work( n+1 ),
755  $ nwork-n, iinfo )
756  ntest = 2
757 *
758  CALL dhst01( n, ilo, ihi, a, lda, h, lda, u, ldu, work,
759  $ nwork, result( 1 ) )
760 *
761 * Call DHSEQR to compute T1, T2 and Z, do tests.
762 *
763 * Eigenvalues only (WR3,WI3)
764 *
765  CALL dlacpy( ' ', n, n, h, lda, t2, lda )
766  ntest = 3
767  result( 3 ) = ulpinv
768 *
769  CALL dhseqr( 'E', 'N', n, ilo, ihi, t2, lda, wr3, wi3, uz,
770  $ ldu, work, nwork, iinfo )
771  IF( iinfo.NE.0 ) THEN
772  WRITE( nounit, fmt = 9999 )'DHSEQR(E)', iinfo, n, jtype,
773  $ ioldsd
774  IF( iinfo.LE.n+2 ) THEN
775  info = abs( iinfo )
776  GO TO 250
777  END IF
778  END IF
779 *
780 * Eigenvalues (WR2,WI2) and Full Schur Form (T2)
781 *
782  CALL dlacpy( ' ', n, n, h, lda, t2, lda )
783 *
784  CALL dhseqr( 'S', 'N', n, ilo, ihi, t2, lda, wr2, wi2, uz,
785  $ ldu, work, nwork, iinfo )
786  IF( iinfo.NE.0 .AND. iinfo.LE.n+2 ) THEN
787  WRITE( nounit, fmt = 9999 )'DHSEQR(S)', iinfo, n, jtype,
788  $ ioldsd
789  info = abs( iinfo )
790  GO TO 250
791  END IF
792 *
793 * Eigenvalues (WR1,WI1), Schur Form (T1), and Schur vectors
794 * (UZ)
795 *
796  CALL dlacpy( ' ', n, n, h, lda, t1, lda )
797  CALL dlacpy( ' ', n, n, u, ldu, uz, ldu )
798 *
799  CALL dhseqr( 'S', 'V', n, ilo, ihi, t1, lda, wr1, wi1, uz,
800  $ ldu, work, nwork, iinfo )
801  IF( iinfo.NE.0 .AND. iinfo.LE.n+2 ) THEN
802  WRITE( nounit, fmt = 9999 )'DHSEQR(V)', iinfo, n, jtype,
803  $ ioldsd
804  info = abs( iinfo )
805  GO TO 250
806  END IF
807 *
808 * Compute Z = U' UZ
809 *
810  CALL dgemm( 'T', 'N', n, n, n, one, u, ldu, uz, ldu, zero,
811  $ z, ldu )
812  ntest = 8
813 *
814 * Do Tests 3: | H - Z T Z' | / ( |H| n ulp )
815 * and 4: | I - Z Z' | / ( n ulp )
816 *
817  CALL dhst01( n, ilo, ihi, h, lda, t1, lda, z, ldu, work,
818  $ nwork, result( 3 ) )
819 *
820 * Do Tests 5: | A - UZ T (UZ)' | / ( |A| n ulp )
821 * and 6: | I - UZ (UZ)' | / ( n ulp )
822 *
823  CALL dhst01( n, ilo, ihi, a, lda, t1, lda, uz, ldu, work,
824  $ nwork, result( 5 ) )
825 *
826 * Do Test 7: | T2 - T1 | / ( |T| n ulp )
827 *
828  CALL dget10( n, n, t2, lda, t1, lda, work, result( 7 ) )
829 *
830 * Do Test 8: | W2 - W1 | / ( max(|W1|,|W2|) ulp )
831 *
832  temp1 = zero
833  temp2 = zero
834  DO 130 j = 1, n
835  temp1 = max( temp1, abs( wr1( j ) )+abs( wi1( j ) ),
836  $ abs( wr2( j ) )+abs( wi2( j ) ) )
837  temp2 = max( temp2, abs( wr1( j )-wr2( j ) )+
838  & abs( wi1( j )-wi2( j ) ) )
839  130 CONTINUE
840 *
841  result( 8 ) = temp2 / max( unfl, ulp*max( temp1, temp2 ) )
842 *
843 * Compute the Left and Right Eigenvectors of T
844 *
845 * Compute the Right eigenvector Matrix:
846 *
847  ntest = 9
848  result( 9 ) = ulpinv
849 *
850 * Select last max(N/4,1) real, max(N/4,1) complex eigenvectors
851 *
852  nselc = 0
853  nselr = 0
854  j = n
855  140 CONTINUE
856  IF( wi1( j ).EQ.zero ) THEN
857  IF( nselr.LT.max( n / 4, 1 ) ) THEN
858  nselr = nselr + 1
859  SELECT( j ) = .true.
860  ELSE
861  SELECT( j ) = .false.
862  END IF
863  j = j - 1
864  ELSE
865  IF( nselc.LT.max( n / 4, 1 ) ) THEN
866  nselc = nselc + 1
867  SELECT( j ) = .true.
868  SELECT( j-1 ) = .false.
869  ELSE
870  SELECT( j ) = .false.
871  SELECT( j-1 ) = .false.
872  END IF
873  j = j - 2
874  END IF
875  IF( j.GT.0 )
876  $ GO TO 140
877 *
878  CALL dtrevc( 'Right', 'All', SELECT, n, t1, lda, dumma, ldu,
879  $ evectr, ldu, n, in, work, iinfo )
880  IF( iinfo.NE.0 ) THEN
881  WRITE( nounit, fmt = 9999 )'DTREVC(R,A)', iinfo, n,
882  $ jtype, ioldsd
883  info = abs( iinfo )
884  GO TO 250
885  END IF
886 *
887 * Test 9: | TR - RW | / ( |T| |R| ulp )
888 *
889  CALL dget22( 'N', 'N', 'N', n, t1, lda, evectr, ldu, wr1,
890  $ wi1, work, dumma( 1 ) )
891  result( 9 ) = dumma( 1 )
892  IF( dumma( 2 ).GT.thresh ) THEN
893  WRITE( nounit, fmt = 9998 )'Right', 'DTREVC',
894  $ dumma( 2 ), n, jtype, ioldsd
895  END IF
896 *
897 * Compute selected right eigenvectors and confirm that
898 * they agree with previous right eigenvectors
899 *
900  CALL dtrevc( 'Right', 'Some', SELECT, n, t1, lda, dumma,
901  $ ldu, evectl, ldu, n, in, work, iinfo )
902  IF( iinfo.NE.0 ) THEN
903  WRITE( nounit, fmt = 9999 )'DTREVC(R,S)', iinfo, n,
904  $ jtype, ioldsd
905  info = abs( iinfo )
906  GO TO 250
907  END IF
908 *
909  k = 1
910  match = .true.
911  DO 170 j = 1, n
912  IF( SELECT( j ) .AND. wi1( j ).EQ.zero ) THEN
913  DO 150 jj = 1, n
914  IF( evectr( jj, j ).NE.evectl( jj, k ) ) THEN
915  match = .false.
916  GO TO 180
917  END IF
918  150 CONTINUE
919  k = k + 1
920  ELSE IF( SELECT( j ) .AND. wi1( j ).NE.zero ) THEN
921  DO 160 jj = 1, n
922  IF( evectr( jj, j ).NE.evectl( jj, k ) .OR.
923  $ evectr( jj, j+1 ).NE.evectl( jj, k+1 ) ) THEN
924  match = .false.
925  GO TO 180
926  END IF
927  160 CONTINUE
928  k = k + 2
929  END IF
930  170 CONTINUE
931  180 CONTINUE
932  IF( .NOT.match )
933  $ WRITE( nounit, fmt = 9997 )'Right', 'DTREVC', n, jtype,
934  $ ioldsd
935 *
936 * Compute the Left eigenvector Matrix:
937 *
938  ntest = 10
939  result( 10 ) = ulpinv
940  CALL dtrevc( 'Left', 'All', SELECT, n, t1, lda, evectl, ldu,
941  $ dumma, ldu, n, in, work, iinfo )
942  IF( iinfo.NE.0 ) THEN
943  WRITE( nounit, fmt = 9999 )'DTREVC(L,A)', iinfo, n,
944  $ jtype, ioldsd
945  info = abs( iinfo )
946  GO TO 250
947  END IF
948 *
949 * Test 10: | LT - WL | / ( |T| |L| ulp )
950 *
951  CALL dget22( 'Trans', 'N', 'Conj', n, t1, lda, evectl, ldu,
952  $ wr1, wi1, work, dumma( 3 ) )
953  result( 10 ) = dumma( 3 )
954  IF( dumma( 4 ).GT.thresh ) THEN
955  WRITE( nounit, fmt = 9998 )'Left', 'DTREVC', dumma( 4 ),
956  $ n, jtype, ioldsd
957  END IF
958 *
959 * Compute selected left eigenvectors and confirm that
960 * they agree with previous left eigenvectors
961 *
962  CALL dtrevc( 'Left', 'Some', SELECT, n, t1, lda, evectr,
963  $ ldu, dumma, ldu, n, in, work, iinfo )
964  IF( iinfo.NE.0 ) THEN
965  WRITE( nounit, fmt = 9999 )'DTREVC(L,S)', iinfo, n,
966  $ jtype, ioldsd
967  info = abs( iinfo )
968  GO TO 250
969  END IF
970 *
971  k = 1
972  match = .true.
973  DO 210 j = 1, n
974  IF( SELECT( j ) .AND. wi1( j ).EQ.zero ) THEN
975  DO 190 jj = 1, n
976  IF( evectl( jj, j ).NE.evectr( jj, k ) ) THEN
977  match = .false.
978  GO TO 220
979  END IF
980  190 CONTINUE
981  k = k + 1
982  ELSE IF( SELECT( j ) .AND. wi1( j ).NE.zero ) THEN
983  DO 200 jj = 1, n
984  IF( evectl( jj, j ).NE.evectr( jj, k ) .OR.
985  $ evectl( jj, j+1 ).NE.evectr( jj, k+1 ) ) THEN
986  match = .false.
987  GO TO 220
988  END IF
989  200 CONTINUE
990  k = k + 2
991  END IF
992  210 CONTINUE
993  220 CONTINUE
994  IF( .NOT.match )
995  $ WRITE( nounit, fmt = 9997 )'Left', 'DTREVC', n, jtype,
996  $ ioldsd
997 *
998 * Call DHSEIN for Right eigenvectors of H, do test 11
999 *
1000  ntest = 11
1001  result( 11 ) = ulpinv
1002  DO 230 j = 1, n
1003  SELECT( j ) = .true.
1004  230 CONTINUE
1005 *
1006  CALL dhsein( 'Right', 'Qr', 'Ninitv', SELECT, n, h, lda,
1007  $ wr3, wi3, dumma, ldu, evectx, ldu, n1, in,
1008  $ work, iwork, iwork, iinfo )
1009  IF( iinfo.NE.0 ) THEN
1010  WRITE( nounit, fmt = 9999 )'DHSEIN(R)', iinfo, n, jtype,
1011  $ ioldsd
1012  info = abs( iinfo )
1013  IF( iinfo.LT.0 )
1014  $ GO TO 250
1015  ELSE
1016 *
1017 * Test 11: | HX - XW | / ( |H| |X| ulp )
1018 *
1019 * (from inverse iteration)
1020 *
1021  CALL dget22( 'N', 'N', 'N', n, h, lda, evectx, ldu, wr3,
1022  $ wi3, work, dumma( 1 ) )
1023  IF( dumma( 1 ).LT.ulpinv )
1024  $ result( 11 ) = dumma( 1 )*aninv
1025  IF( dumma( 2 ).GT.thresh ) THEN
1026  WRITE( nounit, fmt = 9998 )'Right', 'DHSEIN',
1027  $ dumma( 2 ), n, jtype, ioldsd
1028  END IF
1029  END IF
1030 *
1031 * Call DHSEIN for Left eigenvectors of H, do test 12
1032 *
1033  ntest = 12
1034  result( 12 ) = ulpinv
1035  DO 240 j = 1, n
1036  SELECT( j ) = .true.
1037  240 CONTINUE
1038 *
1039  CALL dhsein( 'Left', 'Qr', 'Ninitv', SELECT, n, h, lda, wr3,
1040  $ wi3, evecty, ldu, dumma, ldu, n1, in, work,
1041  $ iwork, iwork, iinfo )
1042  IF( iinfo.NE.0 ) THEN
1043  WRITE( nounit, fmt = 9999 )'DHSEIN(L)', iinfo, n, jtype,
1044  $ ioldsd
1045  info = abs( iinfo )
1046  IF( iinfo.LT.0 )
1047  $ GO TO 250
1048  ELSE
1049 *
1050 * Test 12: | YH - WY | / ( |H| |Y| ulp )
1051 *
1052 * (from inverse iteration)
1053 *
1054  CALL dget22( 'C', 'N', 'C', n, h, lda, evecty, ldu, wr3,
1055  $ wi3, work, dumma( 3 ) )
1056  IF( dumma( 3 ).LT.ulpinv )
1057  $ result( 12 ) = dumma( 3 )*aninv
1058  IF( dumma( 4 ).GT.thresh ) THEN
1059  WRITE( nounit, fmt = 9998 )'Left', 'DHSEIN',
1060  $ dumma( 4 ), n, jtype, ioldsd
1061  END IF
1062  END IF
1063 *
1064 * Call DORMHR for Right eigenvectors of A, do test 13
1065 *
1066  ntest = 13
1067  result( 13 ) = ulpinv
1068 *
1069  CALL dormhr( 'Left', 'No transpose', n, n, ilo, ihi, uu,
1070  $ ldu, tau, evectx, ldu, work, nwork, iinfo )
1071  IF( iinfo.NE.0 ) THEN
1072  WRITE( nounit, fmt = 9999 )'DORMHR(R)', iinfo, n, jtype,
1073  $ ioldsd
1074  info = abs( iinfo )
1075  IF( iinfo.LT.0 )
1076  $ GO TO 250
1077  ELSE
1078 *
1079 * Test 13: | AX - XW | / ( |A| |X| ulp )
1080 *
1081 * (from inverse iteration)
1082 *
1083  CALL dget22( 'N', 'N', 'N', n, a, lda, evectx, ldu, wr3,
1084  $ wi3, work, dumma( 1 ) )
1085  IF( dumma( 1 ).LT.ulpinv )
1086  $ result( 13 ) = dumma( 1 )*aninv
1087  END IF
1088 *
1089 * Call DORMHR for Left eigenvectors of A, do test 14
1090 *
1091  ntest = 14
1092  result( 14 ) = ulpinv
1093 *
1094  CALL dormhr( 'Left', 'No transpose', n, n, ilo, ihi, uu,
1095  $ ldu, tau, evecty, ldu, work, nwork, iinfo )
1096  IF( iinfo.NE.0 ) THEN
1097  WRITE( nounit, fmt = 9999 )'DORMHR(L)', iinfo, n, jtype,
1098  $ ioldsd
1099  info = abs( iinfo )
1100  IF( iinfo.LT.0 )
1101  $ GO TO 250
1102  ELSE
1103 *
1104 * Test 14: | YA - WY | / ( |A| |Y| ulp )
1105 *
1106 * (from inverse iteration)
1107 *
1108  CALL dget22( 'C', 'N', 'C', n, a, lda, evecty, ldu, wr3,
1109  $ wi3, work, dumma( 3 ) )
1110  IF( dumma( 3 ).LT.ulpinv )
1111  $ result( 14 ) = dumma( 3 )*aninv
1112  END IF
1113 *
1114 * End of Loop -- Check for RESULT(j) > THRESH
1115 *
1116  250 CONTINUE
1117 *
1118  ntestt = ntestt + ntest
1119  CALL dlafts( 'DHS', n, n, jtype, ntest, result, ioldsd,
1120  $ thresh, nounit, nerrs )
1121 *
1122  260 CONTINUE
1123  270 CONTINUE
1124 *
1125 * Summary
1126 *
1127  CALL dlasum( 'DHS', nounit, nerrs, ntestt )
1128 *
1129  RETURN
1130 *
1131  9999 FORMAT( ' DCHKHS: ', a, ' returned INFO=', i6, '.', / 9x, 'N=',
1132  $ i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1133  9998 FORMAT( ' DCHKHS: ', a, ' Eigenvectors from ', a, ' incorrectly ',
1134  $ 'normalized.', / ' Bits of error=', 0p, g10.3, ',', 9x,
1135  $ 'N=', i6, ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5,
1136  $ ')' )
1137  9997 FORMAT( ' DCHKHS: Selected ', a, ' Eigenvectors from ', a,
1138  $ ' do not match other eigenvectors ', 9x, 'N=', i6,
1139  $ ', JTYPE=', i6, ', ISEED=(', 3( i5, ',' ), i5, ')' )
1140 *
1141 * End of DCHKHS
1142 *
1143  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 dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dget22(TRANSA, TRANSE, TRANSW, N, A, LDA, E, LDE, WR, WI, WORK, RESULT)
DGET22
Definition: dget22.f:169
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 dhst01(N, ILO, IHI, A, LDA, H, LDH, Q, LDQ, WORK, LWORK, RESULT)
DHST01
Definition: dhst01.f:136
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DGEHRD
Definition: dgehrd.f:169
subroutine dget10(M, N, A, LDA, B, LDB, WORK, RESULT)
DGET10
Definition: dget10.f:95
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine dchkhs(NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH, NOUNIT, A, LDA, H, T1, T2, U, LDU, Z, UZ, WR1, WI1, WR2, WI2, WR3, WI3, EVECTL, EVECTR, EVECTY, EVECTX, UU, TAU, WORK, NWORK, IWORK, SELECT, RESULT, INFO)
DCHKHS
Definition: dchkhs.f:414
subroutine dtrevc(SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR, LDVR, MM, M, WORK, INFO)
DTREVC
Definition: dtrevc.f:224
subroutine dorghr(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DORGHR
Definition: dorghr.f:128
subroutine dlafts(TYPE, M, N, IMAT, NTESTS, RESULT, ISEED, THRESH, IOUNIT, IE)
DLAFTS
Definition: dlafts.f:101
subroutine dlasum(TYPE, IOUNIT, IE, NRUN)
DLASUM
Definition: dlasum.f:45
subroutine dormhr(SIDE, TRANS, M, N, ILO, IHI, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMHR
Definition: dormhr.f:180
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dhseqr(JOB, COMPZ, N, ILO, IHI, H, LDH, WR, WI, Z, LDZ, WORK, LWORK, INFO)
DHSEQR
Definition: dhseqr.f:318
subroutine dhsein(SIDE, EIGSRC, INITV, SELECT, N, H, LDH, WR, WI, VL, LDVL, VR, LDVR, MM, M, WORK, IFAILL, IFAILR, INFO)
DHSEIN
Definition: dhsein.f:265