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