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