OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
derrgex.f
Go to the documentation of this file.
1*> \brief \b DERRGEX
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 DERRGE( PATH, NUNIT )
12*
13* .. Scalar Arguments ..
14* CHARACTER*3 PATH
15* INTEGER NUNIT
16* ..
17*
18*
19*> \par Purpose:
20* =============
21*>
22*> \verbatim
23*>
24*> DERRGE tests the error exits for the DOUBLE PRECISION routines
25*> for general matrices.
26*>
27*> Note that this file is used only when the XBLAS are available,
28*> otherwise derrge.f defines this subroutine.
29*> \endverbatim
30*
31* Arguments:
32* ==========
33*
34*> \param[in] PATH
35*> \verbatim
36*> PATH is CHARACTER*3
37*> The LAPACK path name for the routines to be tested.
38*> \endverbatim
39*>
40*> \param[in] NUNIT
41*> \verbatim
42*> NUNIT is INTEGER
43*> The unit number for output.
44*> \endverbatim
45*
46* Authors:
47* ========
48*
49*> \author Univ. of Tennessee
50*> \author Univ. of California Berkeley
51*> \author Univ. of Colorado Denver
52*> \author NAG Ltd.
53*
54*> \ingroup double_lin
55*
56* =====================================================================
57 SUBROUTINE derrge( PATH, NUNIT )
58*
59* -- LAPACK test routine --
60* -- LAPACK is a software package provided by Univ. of Tennessee, --
61* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
62*
63* .. Scalar Arguments ..
64 CHARACTER*3 PATH
65 INTEGER NUNIT
66* ..
67*
68* =====================================================================
69*
70* .. Parameters ..
71 INTEGER NMAX, LW
72 parameter( nmax = 4, lw = 3*nmax )
73* ..
74* .. Local Scalars ..
75 CHARACTER EQ
76 CHARACTER*2 C2
77 INTEGER I, INFO, J, N_ERR_BNDS, NPARAMS
78 DOUBLE PRECISION ANRM, CCOND, RCOND, BERR
79* ..
80* .. Local Arrays ..
81 INTEGER IP( NMAX ), IW( NMAX )
82 DOUBLE PRECISION A( NMAX, NMAX ), AF( NMAX, NMAX ), B( NMAX ),
83 $ C( NMAX ), R( NMAX ), R1( NMAX ), R2( NMAX ),
84 $ W( LW ), X( NMAX ), ERR_BNDS_N( NMAX, 3 ),
85 $ ERR_BNDS_C( NMAX, 3 ), PARAMS( 1 )
86* ..
87* .. External Functions ..
88 LOGICAL LSAMEN
89 EXTERNAL lsamen
90* ..
91* .. External Subroutines ..
92 EXTERNAL alaesm, chkxer, dgbcon, dgbequ, dgbrfs, dgbtf2,
96* ..
97* .. Scalars in Common ..
98 LOGICAL LERR, OK
99 CHARACTER*32 SRNAMT
100 INTEGER INFOT, NOUT
101* ..
102* .. Common blocks ..
103 COMMON / infoc / infot, nout, ok, lerr
104 COMMON / srnamc / srnamt
105* ..
106* .. Intrinsic Functions ..
107 INTRINSIC dble
108* ..
109* .. Executable Statements ..
110*
111 nout = nunit
112 WRITE( nout, fmt = * )
113 c2 = path( 2: 3 )
114*
115* Set the variables to innocuous values.
116*
117 DO 20 j = 1, nmax
118 DO 10 i = 1, nmax
119 a( i, j ) = 1.d0 / dble( i+j )
120 af( i, j ) = 1.d0 / dble( i+j )
121 10 CONTINUE
122 b( j ) = 0.d0
123 r1( j ) = 0.d0
124 r2( j ) = 0.d0
125 w( j ) = 0.d0
126 x( j ) = 0.d0
127 c( j ) = 0.d0
128 r( j ) = 0.d0
129 ip( j ) = j
130 iw( j ) = j
131 20 CONTINUE
132 ok = .true.
133*
134 IF( lsamen( 2, c2, 'GE' ) ) THEN
135*
136* Test error exits of the routines that use the LU decomposition
137* of a general matrix.
138*
139* DGETRF
140*
141 srnamt = 'DGETRF'
142 infot = 1
143 CALL dgetrf( -1, 0, a, 1, ip, info )
144 CALL chkxer( 'DGETRF', infot, nout, lerr, ok )
145 infot = 2
146 CALL dgetrf( 0, -1, a, 1, ip, info )
147 CALL chkxer( 'DGETRF', infot, nout, lerr, ok )
148 infot = 4
149 CALL dgetrf( 2, 1, a, 1, ip, info )
150 CALL chkxer( 'DGETRF', infot, nout, lerr, ok )
151*
152* DGETF2
153*
154 srnamt = 'DGETF2'
155 infot = 1
156 CALL dgetf2( -1, 0, a, 1, ip, info )
157 CALL chkxer( 'DGETF2', infot, nout, lerr, ok )
158 infot = 2
159 CALL dgetf2( 0, -1, a, 1, ip, info )
160 CALL chkxer( 'DGETF2', infot, nout, lerr, ok )
161 infot = 4
162 CALL dgetf2( 2, 1, a, 1, ip, info )
163 CALL chkxer( 'DGETF2', infot, nout, lerr, ok )
164*
165* DGETRI
166*
167 srnamt = 'DGETRI'
168 infot = 1
169 CALL dgetri( -1, a, 1, ip, w, lw, info )
170 CALL chkxer( 'DGETRI', infot, nout, lerr, ok )
171 infot = 3
172 CALL dgetri( 2, a, 1, ip, w, lw, info )
173 CALL chkxer( 'DGETRI', infot, nout, lerr, ok )
174*
175* DGETRS
176*
177 srnamt = 'DGETRS'
178 infot = 1
179 CALL dgetrs( '/', 0, 0, a, 1, ip, b, 1, info )
180 CALL chkxer( 'DGETRS', infot, nout, lerr, ok )
181 infot = 2
182 CALL dgetrs( 'N', -1, 0, a, 1, ip, b, 1, info )
183 CALL chkxer( 'DGETRS', infot, nout, lerr, ok )
184 infot = 3
185 CALL dgetrs( 'N', 0, -1, a, 1, ip, b, 1, info )
186 CALL chkxer( 'dgetrs', INFOT, NOUT, LERR, OK )
187 INFOT = 5
188 CALL DGETRS( 'n', 2, 1, A, 1, IP, B, 2, INFO )
189 CALL CHKXER( 'dgetrs', INFOT, NOUT, LERR, OK )
190 INFOT = 8
191 CALL DGETRS( 'n', 2, 1, A, 2, IP, B, 1, INFO )
192 CALL CHKXER( 'dgetrs', INFOT, NOUT, LERR, OK )
193*
194* DGERFS
195*
196 SRNAMT = 'dgerfs'
197 INFOT = 1
198 CALL DGERFS( '/', 0, 0, A, 1, AF, 1, IP, B, 1, X, 1, R1, R2, W,
199 $ IW, INFO )
200 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
201 INFOT = 2
202 CALL DGERFS( 'n', -1, 0, A, 1, AF, 1, IP, B, 1, X, 1, R1, R2,
203 $ W, IW, INFO )
204 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
205 INFOT = 3
206 CALL DGERFS( 'n', 0, -1, A, 1, AF, 1, IP, B, 1, X, 1, R1, R2,
207 $ W, IW, INFO )
208 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
209 INFOT = 5
210 CALL DGERFS( 'n', 2, 1, A, 1, AF, 2, IP, B, 2, X, 2, R1, R2, W,
211 $ IW, INFO )
212 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
213 INFOT = 7
214 CALL DGERFS( 'n', 2, 1, A, 2, AF, 1, IP, B, 2, X, 2, R1, R2, W,
215 $ IW, INFO )
216 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
217 INFOT = 10
218 CALL DGERFS( 'n', 2, 1, A, 2, AF, 2, IP, B, 1, X, 2, R1, R2, W,
219 $ IW, INFO )
220 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
221 INFOT = 12
222 CALL DGERFS( 'n', 2, 1, A, 2, AF, 2, IP, B, 2, X, 1, R1, R2, W,
223 $ IW, INFO )
224 CALL CHKXER( 'dgerfs', INFOT, NOUT, LERR, OK )
225*
226* DGERFSX
227*
228 N_ERR_BNDS = 3
229 NPARAMS = 0
230 SRNAMT = 'dgerfsx'
231 INFOT = 1
232 CALL DGERFSX( '/', EQ, 0, 0, A, 1, AF, 1, IP, R, C, B, 1, X,
233 $ 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
234 $ NPARAMS, PARAMS, W, IW, INFO )
235 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
236 INFOT = 2
237 EQ = '/'
238 CALL DGERFSX( 'n', EQ, 2, 1, A, 1, AF, 2, IP, R, C, B, 2, X,
239 $ 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
240 $ NPARAMS, PARAMS, W, IW, INFO )
241 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
242 INFOT = 3
243 EQ = 'r'
244 CALL DGERFSX( 'n', EQ, -1, 0, A, 1, AF, 1, IP, R, C, B, 1, X,
245 $ 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
246 $ NPARAMS, PARAMS, W, IW, INFO )
247 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
248 INFOT = 4
249 CALL DGERFSX( 'n', EQ, 0, -1, A, 1, AF, 1, IP, R, C, B, 1, X,
250 $ 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
251 $ NPARAMS, PARAMS, W, IW, INFO )
252 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
253 INFOT = 6
254 CALL DGERFSX( 'n', EQ, 2, 1, A, 1, AF, 2, IP, R, C, B, 2, X,
255 $ 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
256 $ NPARAMS, PARAMS, W, IW, INFO )
257 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
258 INFOT = 8
259 CALL DGERFSX( 'n', EQ, 2, 1, A, 2, AF, 1, IP, R, C, B, 2, X,
260 $ 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
261 $ NPARAMS, PARAMS, W, IW, INFO )
262 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
263 INFOT = 13
264 EQ = 'c'
265 CALL DGERFSX( 'n', EQ, 2, 1, A, 2, AF, 2, IP, R, C, B, 1, X,
266 $ 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
267 $ NPARAMS, PARAMS, W, IW, INFO )
268 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
269 INFOT = 15
270 CALL DGERFSX( 'n', EQ, 2, 1, A, 2, AF, 2, IP, R, C, B, 2, X,
271 $ 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
272 $ NPARAMS, PARAMS, W, IW, INFO )
273 CALL CHKXER( 'dgerfsx', INFOT, NOUT, LERR, OK )
274*
275* DGECON
276*
277 SRNAMT = 'dgecon'
278 INFOT = 1
279 CALL DGECON( '/', 0, A, 1, ANRM, RCOND, W, IW, INFO )
280 CALL CHKXER( 'dgecon', INFOT, NOUT, LERR, OK )
281 INFOT = 2
282 CALL DGECON( '1', -1, A, 1, ANRM, RCOND, W, IW, INFO )
283 CALL CHKXER( 'dgecon', INFOT, NOUT, LERR, OK )
284 INFOT = 4
285 CALL DGECON( '1', 2, A, 1, ANRM, RCOND, W, IW, INFO )
286 CALL CHKXER( 'dgecon', INFOT, NOUT, LERR, OK )
287*
288* DGEEQU
289*
290 SRNAMT = 'dgeequ'
291 INFOT = 1
292 CALL DGEEQU( -1, 0, A, 1, R1, R2, RCOND, CCOND, ANRM, INFO )
293 CALL CHKXER( 'dgeequ', INFOT, NOUT, LERR, OK )
294 INFOT = 2
295 CALL DGEEQU( 0, -1, A, 1, R1, R2, RCOND, CCOND, ANRM, INFO )
296 CALL CHKXER( 'dgeequ', INFOT, NOUT, LERR, OK )
297 INFOT = 4
298 CALL DGEEQU( 2, 2, A, 1, R1, R2, RCOND, CCOND, ANRM, INFO )
299 CALL CHKXER( 'dgeequ', INFOT, NOUT, LERR, OK )
300*
301* DGEEQUB
302*
303 SRNAMT = 'dgeequb'
304 INFOT = 1
305 CALL DGEEQUB( -1, 0, A, 1, R1, R2, RCOND, CCOND, ANRM, INFO )
306 CALL CHKXER( 'dgeequb', INFOT, NOUT, LERR, OK )
307 INFOT = 2
308 CALL DGEEQUB( 0, -1, A, 1, R1, R2, RCOND, CCOND, ANRM, INFO )
309 CALL CHKXER( 'dgeequb', INFOT, NOUT, LERR, OK )
310 INFOT = 4
311 CALL DGEEQUB( 2, 2, A, 1, R1, R2, RCOND, CCOND, ANRM, INFO )
312 CALL CHKXER( 'dgeequb', INFOT, NOUT, LERR, OK )
313*
314 ELSE IF( LSAMEN( 2, C2, 'gb' ) ) THEN
315*
316* Test error exits of the routines that use the LU decomposition
317* of a general band matrix.
318*
319* DGBTRF
320*
321 SRNAMT = 'dgbtrf'
322 INFOT = 1
323 CALL DGBTRF( -1, 0, 0, 0, A, 1, IP, INFO )
324 CALL CHKXER( 'dgbtrf', INFOT, NOUT, LERR, OK )
325 INFOT = 2
326 CALL DGBTRF( 0, -1, 0, 0, A, 1, IP, INFO )
327 CALL CHKXER( 'dgbtrf', INFOT, NOUT, LERR, OK )
328 INFOT = 3
329 CALL DGBTRF( 1, 1, -1, 0, A, 1, IP, INFO )
330 CALL CHKXER( 'dgbtrf', INFOT, NOUT, LERR, OK )
331 INFOT = 4
332 CALL DGBTRF( 1, 1, 0, -1, A, 1, IP, INFO )
333 CALL CHKXER( 'dgbtrf', INFOT, NOUT, LERR, OK )
334 INFOT = 6
335 CALL DGBTRF( 2, 2, 1, 1, A, 3, IP, INFO )
336 CALL CHKXER( 'dgbtrf', INFOT, NOUT, LERR, OK )
337*
338* DGBTF2
339*
340 SRNAMT = 'dgbtf2'
341 INFOT = 1
342 CALL DGBTF2( -1, 0, 0, 0, A, 1, IP, INFO )
343 CALL CHKXER( 'dgbtf2', INFOT, NOUT, LERR, OK )
344 INFOT = 2
345 CALL DGBTF2( 0, -1, 0, 0, A, 1, IP, INFO )
346 CALL CHKXER( 'dgbtf2', INFOT, NOUT, LERR, OK )
347 INFOT = 3
348 CALL DGBTF2( 1, 1, -1, 0, A, 1, IP, INFO )
349 CALL CHKXER( 'dgbtf2', INFOT, NOUT, LERR, OK )
350 INFOT = 4
351 CALL DGBTF2( 1, 1, 0, -1, A, 1, IP, INFO )
352 CALL CHKXER( 'dgbtf2', INFOT, NOUT, LERR, OK )
353 INFOT = 6
354 CALL DGBTF2( 2, 2, 1, 1, A, 3, IP, INFO )
355 CALL CHKXER( 'dgbtf2', INFOT, NOUT, LERR, OK )
356*
357* DGBTRS
358*
359 SRNAMT = 'dgbtrs'
360 INFOT = 1
361 CALL DGBTRS( '/', 0, 0, 0, 1, A, 1, IP, B, 1, INFO )
362 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
363 INFOT = 2
364 CALL DGBTRS( 'n', -1, 0, 0, 1, A, 1, IP, B, 1, INFO )
365 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
366 INFOT = 3
367 CALL DGBTRS( 'n', 1, -1, 0, 1, A, 1, IP, B, 1, INFO )
368 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
369 INFOT = 4
370 CALL DGBTRS( 'n', 1, 0, -1, 1, A, 1, IP, B, 1, INFO )
371 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
372 INFOT = 5
373 CALL DGBTRS( 'n', 1, 0, 0, -1, A, 1, IP, B, 1, INFO )
374 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
375 INFOT = 7
376 CALL DGBTRS( 'n', 2, 1, 1, 1, A, 3, IP, B, 2, INFO )
377 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
378 INFOT = 10
379 CALL DGBTRS( 'n', 2, 0, 0, 1, A, 1, IP, B, 1, INFO )
380 CALL CHKXER( 'dgbtrs', INFOT, NOUT, LERR, OK )
381*
382* DGBRFS
383*
384 SRNAMT = 'dgbrfs'
385 INFOT = 1
386 CALL DGBRFS( '/', 0, 0, 0, 0, A, 1, AF, 1, IP, B, 1, X, 1, R1,
387 $ R2, W, IW, INFO )
388 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
389 INFOT = 2
390 CALL DGBRFS( 'n', -1, 0, 0, 0, A, 1, AF, 1, IP, B, 1, X, 1, R1,
391 $ R2, W, IW, INFO )
392 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
393 INFOT = 3
394 CALL DGBRFS( 'n', 1, -1, 0, 0, A, 1, AF, 1, IP, B, 1, X, 1, R1,
395 $ R2, W, IW, INFO )
396 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
397 INFOT = 4
398 CALL DGBRFS( 'n', 1, 0, -1, 0, A, 1, AF, 1, IP, B, 1, X, 1, R1,
399 $ R2, W, IW, INFO )
400 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
401 INFOT = 5
402 CALL DGBRFS( 'n', 1, 0, 0, -1, A, 1, AF, 1, IP, B, 1, X, 1, R1,
403 $ R2, W, IW, INFO )
404 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
405 INFOT = 7
406 CALL DGBRFS( 'n', 2, 1, 1, 1, A, 2, AF, 4, IP, B, 2, X, 2, R1,
407 $ R2, W, IW, INFO )
408 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
409 INFOT = 9
410 CALL DGBRFS( 'n', 2, 1, 1, 1, A, 3, AF, 3, IP, B, 2, X, 2, R1,
411 $ R2, W, IW, INFO )
412 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
413 INFOT = 12
414 CALL DGBRFS( 'n', 2, 0, 0, 1, A, 1, AF, 1, IP, B, 1, X, 2, R1,
415 $ R2, W, IW, INFO )
416 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
417 INFOT = 14
418 CALL DGBRFS( 'n', 2, 0, 0, 1, A, 1, AF, 1, IP, B, 2, X, 1, R1,
419 $ R2, W, IW, INFO )
420 CALL CHKXER( 'dgbrfs', INFOT, NOUT, LERR, OK )
421*
422* DGBRFSX
423*
424 N_ERR_BNDS = 3
425 NPARAMS = 0
426 SRNAMT = 'dgbrfsx'
427 INFOT = 1
428 CALL DGBRFSX( '/', EQ, 0, 0, 0, 0, A, 1, AF, 1, IP, R, C, B, 1,
429 $ X, 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
430 $ NPARAMS, PARAMS, W, IW, INFO )
431 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
432 INFOT = 2
433 EQ = '/'
434 CALL DGBRFSX( 'n', EQ, 2, 1, 1, 1, A, 1, AF, 2, IP, R, C, B, 2,
435 $ X, 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
436 $ NPARAMS, PARAMS, W, IW, INFO )
437 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
438 INFOT = 3
439 EQ = 'r'
440 CALL DGBRFSX( 'n', EQ, -1, 1, 1, 0, A, 1, AF, 1, IP, R, C, B,
441 $ 1, X, 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
442 $ NPARAMS, PARAMS, W, IW, INFO )
443 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
444 INFOT = 4
445 EQ = 'r'
446 CALL DGBRFSX( 'n', EQ, 2, -1, 1, 1, A, 3, AF, 4, IP, R, C, B,
447 $ 1, X, 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
448 $ NPARAMS, PARAMS, W, IW, INFO )
449 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
450 INFOT = 5
451 EQ = 'r'
452 CALL DGBRFSX( 'n', EQ, 2, 1, -1, 1, A, 3, AF, 4, IP, R, C, B,
453 $ 1, X, 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
454 $ NPARAMS, PARAMS, W, IW, INFO )
455 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
456 INFOT = 6
457 CALL DGBRFSX( 'n', EQ, 0, 0, 0, -1, A, 1, AF, 1, IP, R, C, B,
458 $ 1, X, 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
459 $ NPARAMS, PARAMS, W, IW, INFO )
460 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
461 INFOT = 8
462 CALL DGBRFSX( 'n', EQ, 2, 1, 1, 1, A, 1, AF, 2, IP, R, C, B,
463 $ 2, X, 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
464 $ NPARAMS, PARAMS, W, IW, INFO )
465 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
466 INFOT = 10
467 CALL DGBRFSX( 'n', EQ, 2, 1, 1, 1, A, 3, AF, 3, IP, R, C, B, 2,
468 $ X, 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
469 $ NPARAMS, PARAMS, W, IW, INFO )
470 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
471 INFOT = 13
472 EQ = 'c'
473 CALL DGBRFSX( 'n', EQ, 2, 1, 1, 1, A, 3, AF, 5, IP, R, C, B,
474 $ 1, X, 2, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
475 $ NPARAMS, PARAMS, W, IW, INFO )
476 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
477 INFOT = 15
478 CALL DGBRFSX( 'n', EQ, 2, 1, 1, 1, A, 3, AF, 5, IP, R, C, B, 2,
479 $ X, 1, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_N, ERR_BNDS_C,
480 $ NPARAMS, PARAMS, W, IW, INFO )
481 CALL CHKXER( 'dgbrfsx', INFOT, NOUT, LERR, OK )
482*
483* DGBCON
484*
485 SRNAMT = 'dgbcon'
486 INFOT = 1
487 CALL DGBCON( '/', 0, 0, 0, A, 1, IP, ANRM, RCOND, W, IW, INFO )
488 CALL CHKXER( 'dgbcon', INFOT, NOUT, LERR, OK )
489 INFOT = 2
490 CALL DGBCON( '1', -1, 0, 0, A, 1, IP, ANRM, RCOND, W, IW,
491 $ INFO )
492 CALL CHKXER( 'dgbcon', INFOT, NOUT, LERR, OK )
493 INFOT = 3
494 CALL DGBCON( '1', 1, -1, 0, A, 1, IP, ANRM, RCOND, W, IW,
495 $ INFO )
496 CALL CHKXER( 'dgbcon', INFOT, NOUT, LERR, OK )
497 INFOT = 4
498 CALL DGBCON( '1', 1, 0, -1, A, 1, IP, ANRM, RCOND, W, IW,
499 $ INFO )
500 CALL CHKXER( 'dgbcon', INFOT, NOUT, LERR, OK )
501 INFOT = 6
502 CALL DGBCON( '1', 2, 1, 1, A, 3, IP, ANRM, RCOND, W, IW, INFO )
503 CALL CHKXER( 'dgbcon', INFOT, NOUT, LERR, OK )
504*
505* DGBEQU
506*
507 SRNAMT = 'dgbequ'
508 INFOT = 1
509 CALL DGBEQU( -1, 0, 0, 0, A, 1, R1, R2, RCOND, CCOND, ANRM,
510 $ INFO )
511 CALL CHKXER( 'dgbequ', INFOT, NOUT, LERR, OK )
512 INFOT = 2
513 CALL DGBEQU( 0, -1, 0, 0, A, 1, R1, R2, RCOND, CCOND, ANRM,
514 $ INFO )
515 CALL CHKXER( 'dgbequ', INFOT, NOUT, LERR, OK )
516 INFOT = 3
517 CALL DGBEQU( 1, 1, -1, 0, A, 1, R1, R2, RCOND, CCOND, ANRM,
518 $ INFO )
519 CALL CHKXER( 'dgbequ', INFOT, NOUT, LERR, OK )
520 INFOT = 4
521 CALL DGBEQU( 1, 1, 0, -1, A, 1, R1, R2, RCOND, CCOND, ANRM,
522 $ INFO )
523 CALL CHKXER( 'dgbequ', INFOT, NOUT, LERR, OK )
524 INFOT = 6
525 CALL DGBEQU( 2, 2, 1, 1, A, 2, R1, R2, RCOND, CCOND, ANRM,
526 $ INFO )
527 CALL CHKXER( 'dgbequ', INFOT, NOUT, LERR, OK )
528*
529* DGBEQUB
530*
531 SRNAMT = 'dgbequb'
532 INFOT = 1
533 CALL DGBEQUB( -1, 0, 0, 0, A, 1, R1, R2, RCOND, CCOND, ANRM,
534 $ INFO )
535 CALL CHKXER( 'dgbequb', INFOT, NOUT, LERR, OK )
536 INFOT = 2
537 CALL DGBEQUB( 0, -1, 0, 0, A, 1, R1, R2, RCOND, CCOND, ANRM,
538 $ INFO )
539 CALL CHKXER( 'dgbequb', INFOT, NOUT, LERR, OK )
540 INFOT = 3
541 CALL DGBEQUB( 1, 1, -1, 0, A, 1, R1, R2, RCOND, CCOND, ANRM,
542 $ INFO )
543 CALL CHKXER( 'dgbequb', INFOT, NOUT, LERR, OK )
544 INFOT = 4
545 CALL DGBEQUB( 1, 1, 0, -1, A, 1, R1, R2, RCOND, CCOND, ANRM,
546 $ INFO )
547 CALL CHKXER( 'dgbequb', INFOT, NOUT, LERR, OK )
548 INFOT = 6
549 CALL DGBEQUB( 2, 2, 1, 1, A, 2, R1, R2, RCOND, CCOND, ANRM,
550 $ INFO )
551 CALL CHKXER( 'dgbequb', INFOT, NOUT, LERR, OK )
552 END IF
553*
554* Print a summary line.
555*
556 CALL ALAESM( PATH, OK, NOUT )
557*
558 RETURN
559*
560* End of DERRGEX
561*
562 END
subroutine chkxer(srnamt, infot, nout, lerr, ok)
Definition cblat2.f:3196
subroutine alaesm(path, ok, nout)
ALAESM
Definition alaesm.f:63
subroutine dgbtf2(m, n, kl, ku, ab, ldab, ipiv, info)
DGBTF2 computes the LU factorization of a general band matrix using the unblocked version of the algo...
Definition dgbtf2.f:145
subroutine dgbrfsx(trans, equed, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, r, c, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
DGBRFSX
Definition dgbrfsx.f:440
subroutine dgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
DGBTRS
Definition dgbtrs.f:138
subroutine dgbequb(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, info)
DGBEQUB
Definition dgbequb.f:160
subroutine dgbequ(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, info)
DGBEQU
Definition dgbequ.f:153
subroutine dgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
DGBTRF
Definition dgbtrf.f:144
subroutine dgbrfs(trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DGBRFS
Definition dgbrfs.f:205
subroutine dgbcon(norm, n, kl, ku, ab, ldab, ipiv, anorm, rcond, work, iwork, info)
DGBCON
Definition dgbcon.f:146
subroutine dgecon(norm, n, a, lda, anorm, rcond, work, iwork, info)
DGECON
Definition dgecon.f:124
subroutine dgetf2(m, n, a, lda, ipiv, info)
DGETF2 computes the LU factorization of a general m-by-n matrix using partial pivoting with row inter...
Definition dgetf2.f:108
subroutine dgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
DGEEQU
Definition dgeequ.f:139
subroutine dgetri(n, a, lda, ipiv, work, lwork, info)
DGETRI
Definition dgetri.f:114
subroutine dgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
DGETRS
Definition dgetrs.f:121
subroutine dgerfsx(trans, equed, n, nrhs, a, lda, af, ldaf, ipiv, r, c, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
DGERFSX
Definition dgerfsx.f:414
subroutine dgetrf(m, n, a, lda, ipiv, info)
DGETRF
Definition dgetrf.f:108
subroutine dgeequb(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
DGEEQUB
Definition dgeequb.f:146
subroutine dgerfs(trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DGERFS
Definition dgerfs.f:185
subroutine derrge(path, nunit)
DERRGE
Definition derrge.f:55