LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
alahd.f
Go to the documentation of this file.
1 *> \brief \b ALAHD
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 ALAHD( IOUNIT, PATH )
12 *
13 * .. Scalar Arguments ..
14 * CHARACTER*3 PATH
15 * INTEGER IOUNIT
16 * ..
17 *
18 *
19 *> \par Purpose:
20 * =============
21 *>
22 *> \verbatim
23 *>
24 *> ALAHD prints header information for the different test paths.
25 *> \endverbatim
26 *
27 * Arguments:
28 * ==========
29 *
30 *> \param[in] IOUNIT
31 *> \verbatim
32 *> IOUNIT is INTEGER
33 *> The unit number to which the header information should be
34 *> printed.
35 *> \endverbatim
36 *>
37 *> \param[in] PATH
38 *> \verbatim
39 *> PATH is CHARACTER*3
40 *> The name of the path for which the header information is to
41 *> be printed. Current paths are
42 *> _GE: General matrices
43 *> _GB: General band
44 *> _GT: General Tridiagonal
45 *> _PO: Symmetric or Hermitian positive definite
46 *> _PS: Symmetric or Hermitian positive semi-definite
47 *> _PP: Symmetric or Hermitian positive definite packed
48 *> _PB: Symmetric or Hermitian positive definite band
49 *> _PT: Symmetric or Hermitian positive definite tridiagonal
50 *> _SY: Symmetric indefinite,
51 *> with partial (Bunch-Kaufman) pivoting
52 *> _SR: Symmetric indefinite,
53 *> with "rook" (bounded Bunch-Kaufman) pivoting
54 *> _SP: Symmetric indefinite packed,
55 *> with partial (Bunch-Kaufman) pivoting
56 *> _HE: (complex) Hermitian indefinite,
57 *> with partial (Bunch-Kaufman) pivoting
58 *> _HR: Symmetric indefinite,
59 *> with "rook" (bounded Bunch-Kaufman) pivoting
60 *> _HP: (complex) Hermitian indefinite packed,
61 *> with partial (Bunch-Kaufman) pivoting
62 *> _TR: Triangular
63 *> _TP: Triangular packed
64 *> _TB: Triangular band
65 *> _QR: QR (general matrices)
66 *> _LQ: LQ (general matrices)
67 *> _QL: QL (general matrices)
68 *> _RQ: RQ (general matrices)
69 *> _QP: QR with column pivoting
70 *> _TZ: Trapezoidal
71 *> _LS: Least Squares driver routines
72 *> _LU: LU variants
73 *> _CH: Cholesky variants
74 *> _QS: QR variants
75 *> _QT: QRT (general matrices)
76 *> _QX: QRT (triangular-pentagonal matrices)
77 *> The first character must be one of S, D, C, or Z (C or Z only
78 *> if complex).
79 *> \endverbatim
80 *
81 * Authors:
82 * ========
83 *
84 *> \author Univ. of Tennessee
85 *> \author Univ. of California Berkeley
86 *> \author Univ. of Colorado Denver
87 *> \author NAG Ltd.
88 *
89 *> \date November 2015
90 *
91 *> \ingroup aux_lin
92 *
93 * =====================================================================
94  SUBROUTINE alahd( IOUNIT, PATH )
95 *
96 * -- LAPACK test routine (version 3.6.0) --
97 * -- LAPACK is a software package provided by Univ. of Tennessee, --
98 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
99 * November 2015
100 *
101 * .. Scalar Arguments ..
102  CHARACTER*3 PATH
103  INTEGER IOUNIT
104 * ..
105 *
106 * =====================================================================
107 *
108 * .. Local Scalars ..
109  LOGICAL CORZ, SORD
110  CHARACTER C1, C3
111  CHARACTER*2 P2
112  CHARACTER*4 EIGCNM
113  CHARACTER*32 SUBNAM
114  CHARACTER*9 SYM
115 * ..
116 * .. External Functions ..
117  LOGICAL LSAME, LSAMEN
118  EXTERNAL lsame, lsamen
119 * ..
120 * .. Intrinsic Functions ..
121  INTRINSIC len_trim
122 * ..
123 * .. Executable Statements ..
124 *
125  IF( iounit.LE.0 )
126  $ RETURN
127  c1 = path( 1: 1 )
128  c3 = path( 3: 3 )
129  p2 = path( 2: 3 )
130  sord = lsame( c1, 'S' ) .OR. lsame( c1, 'D' )
131  corz = lsame( c1, 'C' ) .OR. lsame( c1, 'Z' )
132  IF( .NOT.( sord .OR. corz ) )
133  $ RETURN
134 *
135  IF( lsamen( 2, p2, 'GE' ) ) THEN
136 *
137 * GE: General dense
138 *
139  WRITE( iounit, fmt = 9999 )path
140  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
141  WRITE( iounit, fmt = 9979 )
142  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
143  WRITE( iounit, fmt = 9962 )1
144  WRITE( iounit, fmt = 9961 )2
145  WRITE( iounit, fmt = 9960 )3
146  WRITE( iounit, fmt = 9959 )4
147  WRITE( iounit, fmt = 9958 )5
148  WRITE( iounit, fmt = 9957 )6
149  WRITE( iounit, fmt = 9956 )7
150  WRITE( iounit, fmt = 9955 )8
151  WRITE( iounit, fmt = '( '' Messages:'' )' )
152 *
153  ELSE IF( lsamen( 2, p2, 'GB' ) ) THEN
154 *
155 * GB: General band
156 *
157  WRITE( iounit, fmt = 9998 )path
158  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
159  WRITE( iounit, fmt = 9978 )
160  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
161  WRITE( iounit, fmt = 9962 )1
162  WRITE( iounit, fmt = 9960 )2
163  WRITE( iounit, fmt = 9959 )3
164  WRITE( iounit, fmt = 9958 )4
165  WRITE( iounit, fmt = 9957 )5
166  WRITE( iounit, fmt = 9956 )6
167  WRITE( iounit, fmt = 9955 )7
168  WRITE( iounit, fmt = '( '' Messages:'' )' )
169 *
170  ELSE IF( lsamen( 2, p2, 'GT' ) ) THEN
171 *
172 * GT: General tridiagonal
173 *
174  WRITE( iounit, fmt = 9997 )path
175  WRITE( iounit, fmt = 9977 )
176  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
177  WRITE( iounit, fmt = 9962 )1
178  WRITE( iounit, fmt = 9960 )2
179  WRITE( iounit, fmt = 9959 )3
180  WRITE( iounit, fmt = 9958 )4
181  WRITE( iounit, fmt = 9957 )5
182  WRITE( iounit, fmt = 9956 )6
183  WRITE( iounit, fmt = 9955 )7
184  WRITE( iounit, fmt = '( '' Messages:'' )' )
185 *
186  ELSE IF( lsamen( 2, p2, 'PO' ) .OR. lsamen( 2, p2, 'PP' ) ) THEN
187 *
188 * PO: Positive definite full
189 * PP: Positive definite packed
190 *
191  IF( sord ) THEN
192  sym = 'Symmetric'
193  ELSE
194  sym = 'Hermitian'
195  END IF
196  IF( lsame( c3, 'O' ) ) THEN
197  WRITE( iounit, fmt = 9996 )path, sym
198  ELSE
199  WRITE( iounit, fmt = 9995 )path, sym
200  END IF
201  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
202  WRITE( iounit, fmt = 9975 )path
203  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
204  WRITE( iounit, fmt = 9954 )1
205  WRITE( iounit, fmt = 9961 )2
206  WRITE( iounit, fmt = 9960 )3
207  WRITE( iounit, fmt = 9959 )4
208  WRITE( iounit, fmt = 9958 )5
209  WRITE( iounit, fmt = 9957 )6
210  WRITE( iounit, fmt = 9956 )7
211  WRITE( iounit, fmt = 9955 )8
212  WRITE( iounit, fmt = '( '' Messages:'' )' )
213 *
214  ELSE IF( lsamen( 2, p2, 'PS' ) ) THEN
215 *
216 * PS: Positive semi-definite full
217 *
218  IF( sord ) THEN
219  sym = 'Symmetric'
220  ELSE
221  sym = 'Hermitian'
222  END IF
223  IF( lsame( c1, 'S' ) .OR. lsame( c1, 'C' ) ) THEN
224  eigcnm = '1E04'
225  ELSE
226  eigcnm = '1D12'
227  END IF
228  WRITE( iounit, fmt = 9995 )path, sym
229  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
230  WRITE( iounit, fmt = 8973 )eigcnm, eigcnm, eigcnm
231  WRITE( iounit, fmt = '( '' Difference:'' )' )
232  WRITE( iounit, fmt = 8972 )c1
233  WRITE( iounit, fmt = '( '' Test ratio:'' )' )
234  WRITE( iounit, fmt = 8950 )
235  WRITE( iounit, fmt = '( '' Messages:'' )' )
236  ELSE IF( lsamen( 2, p2, 'PB' ) ) THEN
237 *
238 * PB: Positive definite band
239 *
240  IF( sord ) THEN
241  WRITE( iounit, fmt = 9994 )path, 'Symmetric'
242  ELSE
243  WRITE( iounit, fmt = 9994 )path, 'Hermitian'
244  END IF
245  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
246  WRITE( iounit, fmt = 9973 )path
247  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
248  WRITE( iounit, fmt = 9954 )1
249  WRITE( iounit, fmt = 9960 )2
250  WRITE( iounit, fmt = 9959 )3
251  WRITE( iounit, fmt = 9958 )4
252  WRITE( iounit, fmt = 9957 )5
253  WRITE( iounit, fmt = 9956 )6
254  WRITE( iounit, fmt = 9955 )7
255  WRITE( iounit, fmt = '( '' Messages:'' )' )
256 *
257  ELSE IF( lsamen( 2, p2, 'PT' ) ) THEN
258 *
259 * PT: Positive definite tridiagonal
260 *
261  IF( sord ) THEN
262  WRITE( iounit, fmt = 9993 )path, 'Symmetric'
263  ELSE
264  WRITE( iounit, fmt = 9993 )path, 'Hermitian'
265  END IF
266  WRITE( iounit, fmt = 9976 )
267  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
268  WRITE( iounit, fmt = 9952 )1
269  WRITE( iounit, fmt = 9960 )2
270  WRITE( iounit, fmt = 9959 )3
271  WRITE( iounit, fmt = 9958 )4
272  WRITE( iounit, fmt = 9957 )5
273  WRITE( iounit, fmt = 9956 )6
274  WRITE( iounit, fmt = 9955 )7
275  WRITE( iounit, fmt = '( '' Messages:'' )' )
276 *
277  ELSE IF( lsamen( 2, p2, 'SY' ) ) THEN
278 *
279 * SY: Symmetric indefinite full,
280 * with partial (Bunch-Kaufman) pivoting algorithm
281 *
282  IF( lsame( c3, 'Y' ) ) THEN
283  WRITE( iounit, fmt = 9992 )path, 'Symmetric'
284  ELSE
285  WRITE( iounit, fmt = 9991 )path, 'Symmetric'
286  END IF
287  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
288  IF( sord ) THEN
289  WRITE( iounit, fmt = 9972 )
290  ELSE
291  WRITE( iounit, fmt = 9971 )
292  END IF
293  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
294  WRITE( iounit, fmt = 9953 )1
295  WRITE( iounit, fmt = 9961 )2
296  WRITE( iounit, fmt = 9960 )3
297  WRITE( iounit, fmt = 9960 )4
298  WRITE( iounit, fmt = 9959 )5
299  WRITE( iounit, fmt = 9958 )6
300  WRITE( iounit, fmt = 9956 )7
301  WRITE( iounit, fmt = 9957 )8
302  WRITE( iounit, fmt = 9955 )9
303  WRITE( iounit, fmt = '( '' Messages:'' )' )
304 *
305  ELSE IF( lsamen( 2, p2, 'SR' ) ) THEN
306 *
307 * SR: Symmetric indefinite full,
308 * with "rook" (bounded Bunch-Kaufman) pivoting algorithm
309 *
310  WRITE( iounit, fmt = 9892 )path, 'Symmetric'
311 *
312  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
313  IF( sord ) THEN
314  WRITE( iounit, fmt = 9972 )
315  ELSE
316  WRITE( iounit, fmt = 9971 )
317  END IF
318 *
319  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
320  WRITE( iounit, fmt = 9953 )1
321  WRITE( iounit, fmt = 9961 )2
322  WRITE( iounit, fmt = 9927 )3
323  WRITE( iounit, fmt = 9928 )
324  WRITE( iounit, fmt = 9926 )4
325  WRITE( iounit, fmt = 9928 )
326  WRITE( iounit, fmt = 9960 )5
327  WRITE( iounit, fmt = 9959 )6
328  WRITE( iounit, fmt = 9955 )7
329  WRITE( iounit, fmt = '( '' Messages:'' )' )
330 *
331  ELSE IF( lsamen( 2, p2, 'SP' ) ) THEN
332 *
333 * SP: Symmetric indefinite packed,
334 * with partial (Bunch-Kaufman) pivoting algorithm
335 *
336  IF( lsame( c3, 'Y' ) ) THEN
337  WRITE( iounit, fmt = 9992 )path, 'Symmetric'
338  ELSE
339  WRITE( iounit, fmt = 9991 )path, 'Symmetric'
340  END IF
341  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
342  IF( sord ) THEN
343  WRITE( iounit, fmt = 9972 )
344  ELSE
345  WRITE( iounit, fmt = 9971 )
346  END IF
347  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
348  WRITE( iounit, fmt = 9953 )1
349  WRITE( iounit, fmt = 9961 )2
350  WRITE( iounit, fmt = 9960 )3
351  WRITE( iounit, fmt = 9959 )4
352  WRITE( iounit, fmt = 9958 )5
353  WRITE( iounit, fmt = 9956 )6
354  WRITE( iounit, fmt = 9957 )7
355  WRITE( iounit, fmt = 9955 )8
356  WRITE( iounit, fmt = '( '' Messages:'' )' )
357 *
358  ELSE IF( lsamen( 2, p2, 'HE' ) ) THEN
359 *
360 * HE: Hermitian indefinite full,
361 * with partial (Bunch-Kaufman) pivoting algorithm
362 *
363  WRITE( iounit, fmt = 9992 )path, 'Hermitian'
364 *
365  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
366  WRITE( iounit, fmt = 9972 )
367 *
368  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
369  WRITE( iounit, fmt = 9953 )1
370  WRITE( iounit, fmt = 9961 )2
371  WRITE( iounit, fmt = 9960 )3
372  WRITE( iounit, fmt = 9960 )4
373  WRITE( iounit, fmt = 9959 )5
374  WRITE( iounit, fmt = 9958 )6
375  WRITE( iounit, fmt = 9956 )7
376  WRITE( iounit, fmt = 9957 )8
377  WRITE( iounit, fmt = 9955 )9
378  WRITE( iounit, fmt = '( '' Messages:'' )' )
379 *
380  ELSE IF( lsamen( 2, p2, 'HR' ) ) THEN
381 *
382 * HR: Symmetric indefinite full,
383 * with "rook" (bounded Bunch-Kaufman) pivoting algorithm
384 *
385  WRITE( iounit, fmt = 9892 )path, 'Hermitian'
386 *
387  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
388  WRITE( iounit, fmt = 9972 )
389 *
390  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
391  WRITE( iounit, fmt = 9953 )1
392  WRITE( iounit, fmt = 9961 )2
393  WRITE( iounit, fmt = 9927 )3
394  WRITE( iounit, fmt = 9928 )
395  WRITE( iounit, fmt = 9926 )4
396  WRITE( iounit, fmt = 9928 )
397  WRITE( iounit, fmt = 9960 )5
398  WRITE( iounit, fmt = 9959 )6
399  WRITE( iounit, fmt = 9955 )7
400  WRITE( iounit, fmt = '( '' Messages:'' )' )
401 *
402  ELSE IF( lsamen( 2, p2, 'HP' ) ) THEN
403 *
404 * HP: Hermitian indefinite packed,
405 * with partial (Bunch-Kaufman) pivoting algorithm
406 *
407  IF( lsame( c3, 'E' ) ) THEN
408  WRITE( iounit, fmt = 9992 )path, 'Hermitian'
409  ELSE
410  WRITE( iounit, fmt = 9991 )path, 'Hermitian'
411  END IF
412  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
413  WRITE( iounit, fmt = 9972 )
414  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
415  WRITE( iounit, fmt = 9953 )1
416  WRITE( iounit, fmt = 9961 )2
417  WRITE( iounit, fmt = 9960 )3
418  WRITE( iounit, fmt = 9959 )4
419  WRITE( iounit, fmt = 9958 )5
420  WRITE( iounit, fmt = 9956 )6
421  WRITE( iounit, fmt = 9957 )7
422  WRITE( iounit, fmt = 9955 )8
423  WRITE( iounit, fmt = '( '' Messages:'' )' )
424 *
425  ELSE IF( lsamen( 2, p2, 'TR' ) .OR. lsamen( 2, p2, 'TP' ) ) THEN
426 *
427 * TR: Triangular full
428 * TP: Triangular packed
429 *
430  IF( lsame( c3, 'R' ) ) THEN
431  WRITE( iounit, fmt = 9990 )path
432  subnam = path( 1: 1 ) // 'LATRS'
433  ELSE
434  WRITE( iounit, fmt = 9989 )path
435  subnam = path( 1: 1 ) // 'LATPS'
436  END IF
437  WRITE( iounit, fmt = 9966 )path
438  WRITE( iounit, fmt = 9965 )subnam(1:len_trim( subnam ))
439  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
440  WRITE( iounit, fmt = 9961 )1
441  WRITE( iounit, fmt = 9960 )2
442  WRITE( iounit, fmt = 9959 )3
443  WRITE( iounit, fmt = 9958 )4
444  WRITE( iounit, fmt = 9957 )5
445  WRITE( iounit, fmt = 9956 )6
446  WRITE( iounit, fmt = 9955 )7
447  WRITE( iounit, fmt = 9951 )subnam(1:len_trim( subnam )), 8
448  WRITE( iounit, fmt = '( '' Messages:'' )' )
449 *
450  ELSE IF( lsamen( 2, p2, 'TB' ) ) THEN
451 *
452 * TB: Triangular band
453 *
454  WRITE( iounit, fmt = 9988 )path
455  subnam = path( 1: 1 ) // 'LATBS'
456  WRITE( iounit, fmt = 9964 )path
457  WRITE( iounit, fmt = 9963 )subnam(1:len_trim( subnam ))
458  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
459  WRITE( iounit, fmt = 9960 )1
460  WRITE( iounit, fmt = 9959 )2
461  WRITE( iounit, fmt = 9958 )3
462  WRITE( iounit, fmt = 9957 )4
463  WRITE( iounit, fmt = 9956 )5
464  WRITE( iounit, fmt = 9955 )6
465  WRITE( iounit, fmt = 9951 )subnam(1:len_trim( subnam )), 7
466  WRITE( iounit, fmt = '( '' Messages:'' )' )
467 *
468  ELSE IF( lsamen( 2, p2, 'QR' ) ) THEN
469 *
470 * QR decomposition of rectangular matrices
471 *
472  WRITE( iounit, fmt = 9987 )path, 'QR'
473  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
474  WRITE( iounit, fmt = 9970 )
475  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
476  WRITE( iounit, fmt = 9950 )1
477  WRITE( iounit, fmt = 6950 )8
478  WRITE( iounit, fmt = 9946 )2
479  WRITE( iounit, fmt = 9944 )3, 'M'
480  WRITE( iounit, fmt = 9943 )4, 'M'
481  WRITE( iounit, fmt = 9942 )5, 'M'
482  WRITE( iounit, fmt = 9941 )6, 'M'
483  WRITE( iounit, fmt = 9960 )7
484  WRITE( iounit, fmt = 6660 )9
485  WRITE( iounit, fmt = '( '' Messages:'' )' )
486 *
487  ELSE IF( lsamen( 2, p2, 'LQ' ) ) THEN
488 *
489 * LQ decomposition of rectangular matrices
490 *
491  WRITE( iounit, fmt = 9987 )path, 'LQ'
492  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
493  WRITE( iounit, fmt = 9970 )
494  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
495  WRITE( iounit, fmt = 9949 )1
496  WRITE( iounit, fmt = 9945 )2
497  WRITE( iounit, fmt = 9944 )3, 'N'
498  WRITE( iounit, fmt = 9943 )4, 'N'
499  WRITE( iounit, fmt = 9942 )5, 'N'
500  WRITE( iounit, fmt = 9941 )6, 'N'
501  WRITE( iounit, fmt = 9960 )7
502  WRITE( iounit, fmt = '( '' Messages:'' )' )
503 *
504  ELSE IF( lsamen( 2, p2, 'QL' ) ) THEN
505 *
506 * QL decomposition of rectangular matrices
507 *
508  WRITE( iounit, fmt = 9987 )path, 'QL'
509  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
510  WRITE( iounit, fmt = 9970 )
511  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
512  WRITE( iounit, fmt = 9948 )1
513  WRITE( iounit, fmt = 9946 )2
514  WRITE( iounit, fmt = 9944 )3, 'M'
515  WRITE( iounit, fmt = 9943 )4, 'M'
516  WRITE( iounit, fmt = 9942 )5, 'M'
517  WRITE( iounit, fmt = 9941 )6, 'M'
518  WRITE( iounit, fmt = 9960 )7
519  WRITE( iounit, fmt = '( '' Messages:'' )' )
520 *
521  ELSE IF( lsamen( 2, p2, 'RQ' ) ) THEN
522 *
523 * RQ decomposition of rectangular matrices
524 *
525  WRITE( iounit, fmt = 9987 )path, 'RQ'
526  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
527  WRITE( iounit, fmt = 9970 )
528  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
529  WRITE( iounit, fmt = 9947 )1
530  WRITE( iounit, fmt = 9945 )2
531  WRITE( iounit, fmt = 9944 )3, 'N'
532  WRITE( iounit, fmt = 9943 )4, 'N'
533  WRITE( iounit, fmt = 9942 )5, 'N'
534  WRITE( iounit, fmt = 9941 )6, 'N'
535  WRITE( iounit, fmt = 9960 )7
536  WRITE( iounit, fmt = '( '' Messages:'' )' )
537 *
538  ELSE IF( lsamen( 2, p2, 'QP' ) ) THEN
539 *
540 * QR decomposition with column pivoting
541 *
542  WRITE( iounit, fmt = 9986 )path
543  WRITE( iounit, fmt = 9969 )
544  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
545  WRITE( iounit, fmt = 9940 )1
546  WRITE( iounit, fmt = 9939 )2
547  WRITE( iounit, fmt = 9938 )3
548  WRITE( iounit, fmt = '( '' Messages:'' )' )
549 *
550  ELSE IF( lsamen( 2, p2, 'TZ' ) ) THEN
551 *
552 * TZ: Trapezoidal
553 *
554  WRITE( iounit, fmt = 9985 )path
555  WRITE( iounit, fmt = 9968 )
556  WRITE( iounit, fmt = 9929 )c1
557  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
558  WRITE( iounit, fmt = 9940 )1
559  WRITE( iounit, fmt = 9937 )2
560  WRITE( iounit, fmt = 9938 )3
561  WRITE( iounit, fmt = '( '' Messages:'' )' )
562 *
563  ELSE IF( lsamen( 2, p2, 'LS' ) ) THEN
564 *
565 * LS: Least Squares driver routines for
566 * LS, LSD, LSS, LSX and LSY.
567 *
568  WRITE( iounit, fmt = 9984 )path
569  WRITE( iounit, fmt = 9967 )
570  WRITE( iounit, fmt = 9921 )c1, c1, c1, c1
571  WRITE( iounit, fmt = 9935 )1
572  WRITE( iounit, fmt = 9931 )2
573  WRITE( iounit, fmt = 9933 )3
574  WRITE( iounit, fmt = 9935 )4
575  WRITE( iounit, fmt = 9934 )5
576  WRITE( iounit, fmt = 9932 )6
577  WRITE( iounit, fmt = 9920 )
578  WRITE( iounit, fmt = '( '' Messages:'' )' )
579 *
580  ELSE IF( lsamen( 2, p2, 'LU' ) ) THEN
581 *
582 * LU factorization variants
583 *
584  WRITE( iounit, fmt = 9983 )path
585  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
586  WRITE( iounit, fmt = 9979 )
587  WRITE( iounit, fmt = '( '' Test ratio:'' )' )
588  WRITE( iounit, fmt = 9962 )1
589  WRITE( iounit, fmt = '( '' Messages:'' )' )
590 *
591  ELSE IF( lsamen( 2, p2, 'CH' ) ) THEN
592 *
593 * Cholesky factorization variants
594 *
595  WRITE( iounit, fmt = 9982 )path
596  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
597  WRITE( iounit, fmt = 9974 )
598  WRITE( iounit, fmt = '( '' Test ratio:'' )' )
599  WRITE( iounit, fmt = 9954 )1
600  WRITE( iounit, fmt = '( '' Messages:'' )' )
601 *
602  ELSE IF( lsamen( 2, p2, 'QS' ) ) THEN
603 *
604 * QR factorization variants
605 *
606  WRITE( iounit, fmt = 9981 )path
607  WRITE( iounit, fmt = '( '' Matrix types:'' )' )
608  WRITE( iounit, fmt = 9970 )
609  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
610 *
611  ELSE IF( lsamen( 2, p2, 'QT' ) ) THEN
612 *
613 * QRT (general matrices)
614 *
615  WRITE( iounit, fmt = 8000 ) path
616  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
617  WRITE( iounit, fmt = 8011 ) 1
618  WRITE( iounit, fmt = 8012 ) 2
619  WRITE( iounit, fmt = 8013 ) 3
620  WRITE( iounit, fmt = 8014 ) 4
621  WRITE( iounit, fmt = 8015 ) 5
622  WRITE( iounit, fmt = 8016 ) 6
623 *
624  ELSE IF( lsamen( 2, p2, 'QX' ) ) THEN
625 *
626 * QRT (triangular-pentagonal)
627 *
628  WRITE( iounit, fmt = 8001 ) path
629  WRITE( iounit, fmt = '( '' Test ratios:'' )' )
630  WRITE( iounit, fmt = 8017 ) 1
631  WRITE( iounit, fmt = 8018 ) 2
632  WRITE( iounit, fmt = 8019 ) 3
633  WRITE( iounit, fmt = 8020 ) 4
634  WRITE( iounit, fmt = 8021 ) 5
635  WRITE( iounit, fmt = 8022 ) 6
636 *
637  ELSE
638 *
639 * Print error message if no header is available.
640 *
641  WRITE( iounit, fmt = 9980 )path
642  END IF
643 *
644 * First line of header
645 *
646  9999 FORMAT( / 1x, a3, ': General dense matrices' )
647  9998 FORMAT( / 1x, a3, ': General band matrices' )
648  9997 FORMAT( / 1x, a3, ': General tridiagonal' )
649  9996 FORMAT( / 1x, a3, ': ', a9, ' positive definite matrices' )
650  9995 FORMAT( / 1x, a3, ': ', a9, ' positive definite packed matrices'
651  $ )
652  9994 FORMAT( / 1x, a3, ': ', a9, ' positive definite band matrices' )
653  9993 FORMAT( / 1x, a3, ': ', a9, ' positive definite tridiagonal' )
654  9992 FORMAT( / 1x, a3, ': ', a9, ' indefinite matrices',
655  $ ', partial (Bunch-Kaufman) pivoting' )
656  9991 FORMAT( / 1x, a3, ': ', a9, ' indefinite packed matrices',
657  $ ', partial (Bunch-Kaufman) pivoting' )
658  9892 FORMAT( / 1x, a3, ': ', a9, ' indefinite matrices',
659  $ ', "rook" (bounded Bunch-Kaufman) pivoting' )
660  9891 FORMAT( / 1x, a3, ': ', a9, ' indefinite packed matrices',
661  $ ', "rook" (bounded Bunch-Kaufman) pivoting' )
662  9990 FORMAT( / 1x, a3, ': Triangular matrices' )
663  9989 FORMAT( / 1x, a3, ': Triangular packed matrices' )
664  9988 FORMAT( / 1x, a3, ': Triangular band matrices' )
665  9987 FORMAT( / 1x, a3, ': ', a2, ' factorization of general matrices'
666  $ )
667  9986 FORMAT( / 1x, a3, ': QR factorization with column pivoting' )
668  9985 FORMAT( / 1x, a3, ': RQ factorization of trapezoidal matrix' )
669  9984 FORMAT( / 1x, a3, ': Least squares driver routines' )
670  9983 FORMAT( / 1x, a3, ': LU factorization variants' )
671  9982 FORMAT( / 1x, a3, ': Cholesky factorization variants' )
672  9981 FORMAT( / 1x, a3, ': QR factorization variants' )
673  9980 FORMAT( / 1x, a3, ': No header available' )
674  8000 FORMAT( / 1x, a3, ': QRT factorization for general matrices' )
675  8001 FORMAT( / 1x, a3, ': QRT factorization for ',
676  $ 'triangular-pentagonal matrices' )
677 *
678 * GE matrix types
679 *
680  9979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
681  $ '2. Upper triangular', 16x,
682  $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
683  $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
684  $ / 4x, '4. Random, CNDNUM = 2', 13x,
685  $ '10. Scaled near underflow', / 4x, '5. First column zero',
686  $ 14x, '11. Scaled near overflow', / 4x,
687  $ '6. Last column zero' )
688 *
689 * GB matrix types
690 *
691  9978 FORMAT( 4x, '1. Random, CNDNUM = 2', 14x,
692  $ '5. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
693  $ '2. First column zero', 15x, '6. Random, CNDNUM = .01/EPS',
694  $ / 4x, '3. Last column zero', 16x,
695  $ '7. Scaled near underflow', / 4x,
696  $ '4. Last n/2 columns zero', 11x, '8. Scaled near overflow' )
697 *
698 * GT matrix types
699 *
700  9977 FORMAT( ' Matrix types (1-6 have specified condition numbers):',
701  $ / 4x, '1. Diagonal', 24x, '7. Random, unspecified CNDNUM',
702  $ / 4x, '2. Random, CNDNUM = 2', 14x, '8. First column zero',
703  $ / 4x, '3. Random, CNDNUM = sqrt(0.1/EPS)', 2x,
704  $ '9. Last column zero', / 4x, '4. Random, CNDNUM = 0.1/EPS',
705  $ 7x, '10. Last n/2 columns zero', / 4x,
706  $ '5. Scaled near underflow', 10x,
707  $ '11. Scaled near underflow', / 4x,
708  $ '6. Scaled near overflow', 11x, '12. Scaled near overflow' )
709 *
710 * PT matrix types
711 *
712  9976 FORMAT( ' Matrix types (1-6 have specified condition numbers):',
713  $ / 4x, '1. Diagonal', 24x, '7. Random, unspecified CNDNUM',
714  $ / 4x, '2. Random, CNDNUM = 2', 14x,
715  $ '8. First row and column zero', / 4x,
716  $ '3. Random, CNDNUM = sqrt(0.1/EPS)', 2x,
717  $ '9. Last row and column zero', / 4x,
718  $ '4. Random, CNDNUM = 0.1/EPS', 7x,
719  $ '10. Middle row and column zero', / 4x,
720  $ '5. Scaled near underflow', 10x,
721  $ '11. Scaled near underflow', / 4x,
722  $ '6. Scaled near overflow', 11x, '12. Scaled near overflow' )
723 *
724 * PO, PP matrix types
725 *
726  9975 FORMAT( 4x, '1. Diagonal', 24x,
727  $ '6. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
728  $ '2. Random, CNDNUM = 2', 14x, '7. Random, CNDNUM = 0.1/EPS',
729  $ / 3x, '*3. First row and column zero', 7x,
730  $ '8. Scaled near underflow', / 3x,
731  $ '*4. Last row and column zero', 8x,
732  $ '9. Scaled near overflow', / 3x,
733  $ '*5. Middle row and column zero', / 3x,
734  $ '(* - tests error exits from ', a3,
735  $ 'TRF, no test ratios are computed)' )
736 *
737 * CH matrix types
738 *
739  9974 FORMAT( 4x, '1. Diagonal', 24x,
740  $ '6. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
741  $ '2. Random, CNDNUM = 2', 14x, '7. Random, CNDNUM = 0.1/EPS',
742  $ / 3x, '*3. First row and column zero', 7x,
743  $ '8. Scaled near underflow', / 3x,
744  $ '*4. Last row and column zero', 8x,
745  $ '9. Scaled near overflow', / 3x,
746  $ '*5. Middle row and column zero', / 3x,
747  $ '(* - tests error exits, no test ratios are computed)' )
748 *
749 * PS matrix types
750 *
751  8973 FORMAT( 4x, '1. Diagonal', / 4x, '2. Random, CNDNUM = 2', 14x,
752  $ / 3x, '*3. Nonzero eigenvalues of: D(1:RANK-1)=1 and ',
753  $ 'D(RANK) = 1.0/', a4, / 3x,
754  $ '*4. Nonzero eigenvalues of: D(1)=1 and ',
755  $ ' D(2:RANK) = 1.0/', a4, / 3x,
756  $ '*5. Nonzero eigenvalues of: D(I) = ', a4,
757  $ '**(-(I-1)/(RANK-1)) ', ' I=1:RANK', / 4x,
758  $ '6. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
759  $ '7. Random, CNDNUM = 0.1/EPS', / 4x,
760  $ '8. Scaled near underflow', / 4x, '9. Scaled near overflow',
761  $ / 3x, '(* - Semi-definite tests )' )
762  8972 FORMAT( 3x, 'RANK minus computed rank, returned by ', a, 'PSTRF' )
763 *
764 * PB matrix types
765 *
766  9973 FORMAT( 4x, '1. Random, CNDNUM = 2', 14x,
767  $ '5. Random, CNDNUM = sqrt(0.1/EPS)', / 3x,
768  $ '*2. First row and column zero', 7x,
769  $ '6. Random, CNDNUM = 0.1/EPS', / 3x,
770  $ '*3. Last row and column zero', 8x,
771  $ '7. Scaled near underflow', / 3x,
772  $ '*4. Middle row and column zero', 6x,
773  $ '8. Scaled near overflow', / 3x,
774  $ '(* - tests error exits from ', a3,
775  $ 'TRF, no test ratios are computed)' )
776 *
777 * SSY, SSR, SSP, CHE, CHR, CHP matrix types
778 *
779  9972 FORMAT( 4x, '1. Diagonal', 24x,
780  $ '6. Last n/2 rows and columns zero', / 4x,
781  $ '2. Random, CNDNUM = 2', 14x,
782  $ '7. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
783  $ '3. First row and column zero', 7x,
784  $ '8. Random, CNDNUM = 0.1/EPS', / 4x,
785  $ '4. Last row and column zero', 8x,
786  $ '9. Scaled near underflow', / 4x,
787  $ '5. Middle row and column zero', 5x,
788  $ '10. Scaled near overflow' )
789 *
790 * CSY, CSR, CSP matrix types
791 *
792  9971 FORMAT( 4x, '1. Diagonal', 24x,
793  $ '7. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
794  $ '2. Random, CNDNUM = 2', 14x, '8. Random, CNDNUM = 0.1/EPS',
795  $ / 4x, '3. First row and column zero', 7x,
796  $ '9. Scaled near underflow', / 4x,
797  $ '4. Last row and column zero', 7x,
798  $ '10. Scaled near overflow', / 4x,
799  $ '5. Middle row and column zero', 5x,
800  $ '11. Block diagonal matrix', / 4x,
801  $ '6. Last n/2 rows and columns zero' )
802 *
803 * QR matrix types
804 *
805  9970 FORMAT( 4x, '1. Diagonal', 24x,
806  $ '5. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
807  $ '2. Upper triangular', 16x, '6. Random, CNDNUM = 0.1/EPS',
808  $ / 4x, '3. Lower triangular', 16x,
809  $ '7. Scaled near underflow', / 4x, '4. Random, CNDNUM = 2',
810  $ 14x, '8. Scaled near overflow' )
811 *
812 * QP matrix types
813 *
814  9969 FORMAT( ' Matrix types (2-6 have condition 1/EPS):', / 4x,
815  $ '1. Zero matrix', 21x, '4. First n/2 columns fixed', / 4x,
816  $ '2. One small eigenvalue', 12x, '5. Last n/2 columns fixed',
817  $ / 4x, '3. Geometric distribution', 10x,
818  $ '6. Every second column fixed' )
819 *
820 * TZ matrix types
821 *
822  9968 FORMAT( ' Matrix types (2-3 have condition 1/EPS):', / 4x,
823  $ '1. Zero matrix', / 4x, '2. One small eigenvalue', / 4x,
824  $ '3. Geometric distribution' )
825 *
826 * LS matrix types
827 *
828  9967 FORMAT( ' Matrix types (1-3: full rank, 4-6: rank deficient):',
829  $ / 4x, '1 and 4. Normal scaling', / 4x,
830  $ '2 and 5. Scaled near overflow', / 4x,
831  $ '3 and 6. Scaled near underflow' )
832 *
833 * TR, TP matrix types
834 *
835  9966 FORMAT( ' Matrix types for ', a3, ' routines:', / 4x,
836  $ '1. Diagonal', 24x, '6. Scaled near overflow', / 4x,
837  $ '2. Random, CNDNUM = 2', 14x, '7. Identity', / 4x,
838  $ '3. Random, CNDNUM = sqrt(0.1/EPS) ',
839  $ '8. Unit triangular, CNDNUM = 2', / 4x,
840  $ '4. Random, CNDNUM = 0.1/EPS', 8x,
841  $ '9. Unit, CNDNUM = sqrt(0.1/EPS)', / 4x,
842  $ '5. Scaled near underflow', 10x,
843  $ '10. Unit, CNDNUM = 0.1/EPS' )
844  9965 FORMAT( ' Special types for testing ', a, ':', / 3x,
845  $ '11. Matrix elements are O(1), large right hand side', / 3x,
846  $ '12. First diagonal causes overflow,',
847  $ ' offdiagonal column norms < 1', / 3x,
848  $ '13. First diagonal causes overflow,',
849  $ ' offdiagonal column norms > 1', / 3x,
850  $ '14. Growth factor underflows, solution does not overflow',
851  $ / 3x, '15. Small diagonal causes gradual overflow', / 3x,
852  $ '16. One zero diagonal element', / 3x,
853  $ '17. Large offdiagonals cause overflow when adding a column'
854  $ , / 3x, '18. Unit triangular with large right hand side' )
855 *
856 * TB matrix types
857 *
858  9964 FORMAT( ' Matrix types for ', a3, ' routines:', / 4x,
859  $ '1. Random, CNDNUM = 2', 14x, '6. Identity', / 4x,
860  $ '2. Random, CNDNUM = sqrt(0.1/EPS) ',
861  $ '7. Unit triangular, CNDNUM = 2', / 4x,
862  $ '3. Random, CNDNUM = 0.1/EPS', 8x,
863  $ '8. Unit, CNDNUM = sqrt(0.1/EPS)', / 4x,
864  $ '4. Scaled near underflow', 11x,
865  $ '9. Unit, CNDNUM = 0.1/EPS', / 4x,
866  $ '5. Scaled near overflow' )
867  9963 FORMAT( ' Special types for testing ', a, ':', / 3x,
868  $ '10. Matrix elements are O(1), large right hand side', / 3x,
869  $ '11. First diagonal causes overflow,',
870  $ ' offdiagonal column norms < 1', / 3x,
871  $ '12. First diagonal causes overflow,',
872  $ ' offdiagonal column norms > 1', / 3x,
873  $ '13. Growth factor underflows, solution does not overflow',
874  $ / 3x, '14. Small diagonal causes gradual overflow', / 3x,
875  $ '15. One zero diagonal element', / 3x,
876  $ '16. Large offdiagonals cause overflow when adding a column'
877  $ , / 3x, '17. Unit triangular with large right hand side' )
878 *
879 * Test ratios
880 *
881  9962 FORMAT( 3x, i2, ': norm( L * U - A ) / ( N * norm(A) * EPS )' )
882  9961 FORMAT( 3x, i2, ': norm( I - A*AINV ) / ',
883  $ '( N * norm(A) * norm(AINV) * EPS )' )
884  9960 FORMAT( 3x, i2, ': norm( B - A * X ) / ',
885  $ '( norm(A) * norm(X) * EPS )' )
886  6660 FORMAT( 3x, i2, ': diagonal is not non-negative')
887  9959 FORMAT( 3x, i2, ': norm( X - XACT ) / ',
888  $ '( norm(XACT) * CNDNUM * EPS )' )
889  9958 FORMAT( 3x, i2, ': norm( X - XACT ) / ',
890  $ '( norm(XACT) * CNDNUM * EPS ), refined' )
891  9957 FORMAT( 3x, i2, ': norm( X - XACT ) / ',
892  $ '( norm(XACT) * (error bound) )' )
893  9956 FORMAT( 3x, i2, ': (backward error) / EPS' )
894  9955 FORMAT( 3x, i2, ': RCOND * CNDNUM - 1.0' )
895  9954 FORMAT( 3x, i2, ': norm( U'' * U - A ) / ( N * norm(A) * EPS )',
896  $ ', or', / 7x, 'norm( L * L'' - A ) / ( N * norm(A) * EPS )'
897  $ )
898  8950 FORMAT( 3x,
899  $ 'norm( P * U'' * U * P'' - A ) / ( N * norm(A) * EPS )',
900  $ ', or', / 3x,
901  $ 'norm( P * L * L'' * P'' - A ) / ( N * norm(A) * EPS )' )
902  9953 FORMAT( 3x, i2, ': norm( U*D*U'' - A ) / ( N * norm(A) * EPS )',
903  $ ', or', / 7x, 'norm( L*D*L'' - A ) / ( N * norm(A) * EPS )'
904  $ )
905  9952 FORMAT( 3x, i2, ': norm( U''*D*U - A ) / ( N * norm(A) * EPS )',
906  $ ', or', / 7x, 'norm( L*D*L'' - A ) / ( N * norm(A) * EPS )'
907  $ )
908  9951 FORMAT( ' Test ratio for ', a, ':', / 3x, i2,
909  $ ': norm( s*b - A*x ) / ( norm(A) * norm(x) * EPS )' )
910  9950 FORMAT( 3x, i2, ': norm( R - Q'' * A ) / ( M * norm(A) * EPS )' )
911  6950 FORMAT( 3x, i2, ': norm( R - Q'' * A ) / ( M * norm(A) * EPS )
912  $ [RFPG]' )
913  9949 FORMAT( 3x, i2, ': norm( L - A * Q'' ) / ( N * norm(A) * EPS )' )
914  9948 FORMAT( 3x, i2, ': norm( L - Q'' * A ) / ( M * norm(A) * EPS )' )
915  9947 FORMAT( 3x, i2, ': norm( R - A * Q'' ) / ( N * norm(A) * EPS )' )
916  9946 FORMAT( 3x, i2, ': norm( I - Q''*Q ) / ( M * EPS )' )
917  9945 FORMAT( 3x, i2, ': norm( I - Q*Q'' ) / ( N * EPS )' )
918  9944 FORMAT( 3x, i2, ': norm( Q*C - Q*C ) / ', '( ', a1,
919  $ ' * norm(C) * EPS )' )
920  9943 FORMAT( 3x, i2, ': norm( C*Q - C*Q ) / ', '( ', a1,
921  $ ' * norm(C) * EPS )' )
922  9942 FORMAT( 3x, i2, ': norm( Q''*C - Q''*C )/ ', '( ', a1,
923  $ ' * norm(C) * EPS )' )
924  9941 FORMAT( 3x, i2, ': norm( C*Q'' - C*Q'' )/ ', '( ', a1,
925  $ ' * norm(C) * EPS )' )
926  9940 FORMAT( 3x, i2, ': norm(svd(A) - svd(R)) / ',
927  $ '( M * norm(svd(R)) * EPS )' )
928  9939 FORMAT( 3x, i2, ': norm( A*P - Q*R ) / ( M * norm(A) * EPS )'
929  $ )
930  9938 FORMAT( 3x, i2, ': norm( I - Q''*Q ) / ( M * EPS )' )
931  9937 FORMAT( 3x, i2, ': norm( A - R*Q ) / ( M * norm(A) * EPS )'
932  $ )
933  9935 FORMAT( 3x, i2, ': norm( B - A * X ) / ',
934  $ '( max(M,N) * norm(A) * norm(X) * EPS )' )
935  9934 FORMAT( 3x, i2, ': norm( (A*X-B)'' *A ) / ',
936  $ '( max(M,N,NRHS) * norm(A) * norm(B) * EPS )' )
937  9933 FORMAT( 3x, i2, ': norm(svd(A)-svd(R)) / ',
938  $ '( min(M,N) * norm(svd(R)) * EPS )' )
939  9932 FORMAT( 3x, i2, ': Check if X is in the row space of A or A''' )
940  9931 FORMAT( 3x, i2, ': norm( (A*X-B)'' *A ) / ',
941  $ '( max(M,N,NRHS) * norm(A) * norm(B) * EPS )', / 7x,
942  $ 'if TRANS=''N'.GE.' and MN or TRANS=''T'.LT.' and MN, ',
943  $ 'otherwise', / 7x,
944  $ 'check if X is in the row space of A or A'' ',
945  $ '(overdetermined case)' )
946  9929 FORMAT( ' Test ratios (1-3: ', a1, 'TZRZF):' )
947  9920 FORMAT( 3x, ' 7-10: same as 3-6', 3x, ' 11-14: same as 3-6' )
948  9921 FORMAT( ' Test ratios:', / ' (1-2: ', a1, 'GELS, 3-6: ', a1,
949  $ 'GELSY, 7-10: ', a1, 'GELSS, 11-14: ', a1, 'GELSD)' )
950  9928 FORMAT( 7x, 'where ALPHA = ( 1 + SQRT( 17 ) ) / 8' )
951  9927 FORMAT( 3x, i2, ': ABS( Largest element in L )', / 12x,
952  $ ' - ( 1 / ( 1 - ALPHA ) ) + THRESH' )
953  9926 FORMAT( 3x, i2, ': Largest 2-Norm of 2-by-2 pivots', / 12x,
954  $ ' - ( ( 1 + ALPHA ) / ( 1 - ALPHA ) ) + THRESH' )
955  8011 FORMAT(3x,i2,': norm( R - Q''*A ) / ( M * norm(A) * EPS )' )
956  8012 FORMAT(3x,i2,': norm( I - Q''*Q ) / ( M * EPS )' )
957  8013 FORMAT(3x,i2,': norm( Q*C - Q*C ) / ( M * norm(C) * EPS )' )
958  8014 FORMAT(3x,i2,': norm( Q''*C - Q''*C ) / ( M * norm(C) * EPS )')
959  8015 FORMAT(3x,i2,': norm( C*Q - C*Q ) / ( M * norm(C) * EPS )' )
960  8016 FORMAT(3x,i2,': norm( C*Q'' - C*Q'' ) / ( M * norm(C) * EPS )')
961  8017 FORMAT(3x,i2,': norm( R - Q''*A ) / ( (M+N) * norm(A) * EPS )' )
962  8018 FORMAT(3x,i2,': norm( I - Q''*Q ) / ( (M+N) * EPS )' )
963  8019 FORMAT(3x,i2,': norm( Q*C - Q*C ) / ( (M+N) * norm(C) * EPS )' )
964  8020 FORMAT(3x,i2,
965  $ ': norm( Q''*C - Q''*C ) / ( (M+N) * norm(C) * EPS )')
966  8021 FORMAT(3x,i2,': norm( C*Q - C*Q ) / ( (M+N) * norm(C) * EPS )' )
967  8022 FORMAT(3x,i2,
968  $ ': norm( C*Q'' - C*Q'' ) / ( (M+N) * norm(C) * EPS )')
969 *
970  RETURN
971 *
972 * End of ALAHD
973 *
974  END
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95