LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvgbx.f
Go to the documentation of this file.
1*> \brief \b SDRVGBX
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 SDRVGB( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, LA,
12* AFB, LAFB, ASAV, B, BSAV, X, XACT, S, WORK,
13* RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER LA, LAFB, NN, NOUT, NRHS
18* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NVAL( * )
23* REAL A( * ), AFB( * ), ASAV( * ), B( * ), BSAV( * ),
24* $ RWORK( * ), S( * ), WORK( * ), X( * ),
25* $ XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> SDRVGB tests the driver routines SGBSV, -SVX, and -SVXX.
35*>
36*> Note that this file is used only when the XBLAS are available,
37*> otherwise sdrvgb.f defines this subroutine.
38*> \endverbatim
39*
40* Arguments:
41* ==========
42*
43*> \param[in] DOTYPE
44*> \verbatim
45*> DOTYPE is LOGICAL array, dimension (NTYPES)
46*> The matrix types to be used for testing. Matrices of type j
47*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
48*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
49*> \endverbatim
50*>
51*> \param[in] NN
52*> \verbatim
53*> NN is INTEGER
54*> The number of values of N contained in the vector NVAL.
55*> \endverbatim
56*>
57*> \param[in] NVAL
58*> \verbatim
59*> NVAL is INTEGER array, dimension (NN)
60*> The values of the matrix column dimension N.
61*> \endverbatim
62*>
63*> \param[in] NRHS
64*> \verbatim
65*> NRHS is INTEGER
66*> The number of right hand side vectors to be generated for
67*> each linear system.
68*> \endverbatim
69*>
70*> \param[in] THRESH
71*> \verbatim
72*> THRESH is REAL
73*> The threshold value for the test ratios. A result is
74*> included in the output file if RESULT >= THRESH. To have
75*> every test ratio printed, use THRESH = 0.
76*> \endverbatim
77*>
78*> \param[in] TSTERR
79*> \verbatim
80*> TSTERR is LOGICAL
81*> Flag that indicates whether error exits are to be tested.
82*> \endverbatim
83*>
84*> \param[out] A
85*> \verbatim
86*> A is REAL array, dimension (LA)
87*> \endverbatim
88*>
89*> \param[in] LA
90*> \verbatim
91*> LA is INTEGER
92*> The length of the array A. LA >= (2*NMAX-1)*NMAX
93*> where NMAX is the largest entry in NVAL.
94*> \endverbatim
95*>
96*> \param[out] AFB
97*> \verbatim
98*> AFB is REAL array, dimension (LAFB)
99*> \endverbatim
100*>
101*> \param[in] LAFB
102*> \verbatim
103*> LAFB is INTEGER
104*> The length of the array AFB. LAFB >= (3*NMAX-2)*NMAX
105*> where NMAX is the largest entry in NVAL.
106*> \endverbatim
107*>
108*> \param[out] ASAV
109*> \verbatim
110*> ASAV is REAL array, dimension (LA)
111*> \endverbatim
112*>
113*> \param[out] B
114*> \verbatim
115*> B is REAL array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] BSAV
119*> \verbatim
120*> BSAV is REAL array, dimension (NMAX*NRHS)
121*> \endverbatim
122*>
123*> \param[out] X
124*> \verbatim
125*> X is REAL array, dimension (NMAX*NRHS)
126*> \endverbatim
127*>
128*> \param[out] XACT
129*> \verbatim
130*> XACT is REAL array, dimension (NMAX*NRHS)
131*> \endverbatim
132*>
133*> \param[out] S
134*> \verbatim
135*> S is REAL array, dimension (2*NMAX)
136*> \endverbatim
137*>
138*> \param[out] WORK
139*> \verbatim
140*> WORK is REAL array, dimension
141*> (NMAX*max(3,NRHS,NMAX))
142*> \endverbatim
143*>
144*> \param[out] RWORK
145*> \verbatim
146*> RWORK is REAL array, dimension
147*> (max(2*NMAX,NMAX+2*NRHS))
148*> \endverbatim
149*>
150*> \param[out] IWORK
151*> \verbatim
152*> IWORK is INTEGER array, dimension (2*NMAX)
153*> \endverbatim
154*>
155*> \param[in] NOUT
156*> \verbatim
157*> NOUT is INTEGER
158*> The unit number for output.
159*> \endverbatim
160*
161* Authors:
162* ========
163*
164*> \author Univ. of Tennessee
165*> \author Univ. of California Berkeley
166*> \author Univ. of Colorado Denver
167*> \author NAG Ltd.
168*
169*> \ingroup single_lin
170*
171* =====================================================================
172 SUBROUTINE sdrvgb( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, LA,
173 $ AFB, LAFB, ASAV, B, BSAV, X, XACT, S, WORK,
174 $ RWORK, IWORK, NOUT )
175*
176* -- LAPACK test routine --
177* -- LAPACK is a software package provided by Univ. of Tennessee, --
178* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
179*
180* .. Scalar Arguments ..
181 LOGICAL TSTERR
182 INTEGER LA, LAFB, NN, NOUT, NRHS
183 REAL THRESH
184* ..
185* .. Array Arguments ..
186 LOGICAL DOTYPE( * )
187 INTEGER IWORK( * ), NVAL( * )
188 REAL A( * ), AFB( * ), ASAV( * ), B( * ), BSAV( * ),
189 $ rwork( * ), s( * ), work( * ), x( * ),
190 $ xact( * )
191* ..
192*
193* =====================================================================
194*
195* .. Parameters ..
196 REAL ONE, ZERO
197 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
198 INTEGER NTYPES
199 parameter( ntypes = 8 )
200 INTEGER NTESTS
201 parameter( ntests = 7 )
202 INTEGER NTRAN
203 parameter( ntran = 3 )
204* ..
205* .. Local Scalars ..
206 LOGICAL EQUIL, NOFACT, PREFAC, TRFCON, ZEROT
207 CHARACTER DIST, EQUED, FACT, TRANS, TYPE, XTYPE
208 CHARACTER*3 PATH
209 INTEGER I, I1, I2, IEQUED, IFACT, IKL, IKU, IMAT, IN,
210 $ info, ioff, itran, izero, j, k, k1, kl, ku,
211 $ lda, ldafb, ldb, mode, n, nb, nbmin, nerrs,
212 $ nfact, nfail, nimat, nkl, nku, nrun, nt,
213 $ n_err_bnds
214 REAL AINVNM, AMAX, ANORM, ANORMI, ANORMO, ANRMPV,
215 $ CNDNUM, COLCND, RCOND, RCONDC, RCONDI, RCONDO,
216 $ roldc, roldi, roldo, rowcnd, rpvgrw,
217 $ rpvgrw_svxx
218* ..
219* .. Local Arrays ..
220 CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN )
221 INTEGER ISEED( 4 ), ISEEDY( 4 )
222 REAL RESULT( NTESTS ), BERR( NRHS ),
223 $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
224* ..
225* .. External Functions ..
226 LOGICAL LSAME
227 REAL SGET06, SLAMCH, SLANGB, SLANGE, SLANTB,
229 EXTERNAL lsame, sget06, slamch, slangb, slange, slantb,
231* ..
232* .. External Subroutines ..
233 EXTERNAL aladhd, alaerh, alasvm, serrvx, sgbequ, sgbsv,
237* ..
238* .. Intrinsic Functions ..
239 INTRINSIC abs, max, min
240* ..
241* .. Scalars in Common ..
242 LOGICAL LERR, OK
243 CHARACTER*32 SRNAMT
244 INTEGER INFOT, NUNIT
245* ..
246* .. Common blocks ..
247 COMMON / infoc / infot, nunit, ok, lerr
248 COMMON / srnamc / srnamt
249* ..
250* .. Data statements ..
251 DATA iseedy / 1988, 1989, 1990, 1991 /
252 DATA transs / 'N', 'T', 'C' /
253 DATA facts / 'F', 'N', 'E' /
254 DATA equeds / 'N', 'R', 'C', 'B' /
255* ..
256* .. Executable Statements ..
257*
258* Initialize constants and the random number seed.
259*
260 path( 1: 1 ) = 'Single precision'
261 path( 2: 3 ) = 'GB'
262 nrun = 0
263 nfail = 0
264 nerrs = 0
265 DO 10 i = 1, 4
266 iseed( i ) = iseedy( i )
267 10 CONTINUE
268*
269* Test the error exits
270*
271 IF( tsterr )
272 $ CALL serrvx( path, nout )
273 infot = 0
274*
275* Set the block size and minimum block size for testing.
276*
277 nb = 1
278 nbmin = 2
279 CALL xlaenv( 1, nb )
280 CALL xlaenv( 2, nbmin )
281*
282* Do for each value of N in NVAL
283*
284 DO 150 in = 1, nn
285 n = nval( in )
286 ldb = max( n, 1 )
287 xtype = 'N'
288*
289* Set limits on the number of loop iterations.
290*
291 nkl = max( 1, min( n, 4 ) )
292 IF( n.EQ.0 )
293 $ nkl = 1
294 nku = nkl
295 nimat = ntypes
296 IF( n.LE.0 )
297 $ nimat = 1
298*
299 DO 140 ikl = 1, nkl
300*
301* Do for KL = 0, N-1, (3N-1)/4, and (N+1)/4. This order makes
302* it easier to skip redundant values for small values of N.
303*
304 IF( ikl.EQ.1 ) THEN
305 kl = 0
306 ELSE IF( ikl.EQ.2 ) THEN
307 kl = max( n-1, 0 )
308 ELSE IF( ikl.EQ.3 ) THEN
309 kl = ( 3*n-1 ) / 4
310 ELSE IF( ikl.EQ.4 ) THEN
311 kl = ( n+1 ) / 4
312 END IF
313 DO 130 iku = 1, nku
314*
315* Do for KU = 0, N-1, (3N-1)/4, and (N+1)/4. This order
316* makes it easier to skip redundant values for small
317* values of N.
318*
319 IF( iku.EQ.1 ) THEN
320 ku = 0
321 ELSE IF( iku.EQ.2 ) THEN
322 ku = max( n-1, 0 )
323 ELSE IF( iku.EQ.3 ) THEN
324 ku = ( 3*n-1 ) / 4
325 ELSE IF( iku.EQ.4 ) THEN
326 ku = ( n+1 ) / 4
327 END IF
328*
329* Check that A and AFB are big enough to generate this
330* matrix.
331*
332 lda = kl + ku + 1
333 ldafb = 2*kl + ku + 1
334 IF( lda*n.GT.la .OR. ldafb*n.GT.lafb ) THEN
335 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
336 $ CALL aladhd( nout, path )
337 IF( lda*n.GT.la ) THEN
338 WRITE( nout, fmt = 9999 )la, n, kl, ku,
339 $ n*( kl+ku+1 )
340 nerrs = nerrs + 1
341 END IF
342 IF( ldafb*n.GT.lafb ) THEN
343 WRITE( nout, fmt = 9998 )lafb, n, kl, ku,
344 $ n*( 2*kl+ku+1 )
345 nerrs = nerrs + 1
346 END IF
347 GO TO 130
348 END IF
349*
350 DO 120 imat = 1, nimat
351*
352* Do the tests only if DOTYPE( IMAT ) is true.
353*
354 IF( .NOT.dotype( imat ) )
355 $ GO TO 120
356*
357* Skip types 2, 3, or 4 if the matrix is too small.
358*
359 zerot = imat.GE.2 .AND. imat.LE.4
360 IF( zerot .AND. n.LT.imat-1 )
361 $ GO TO 120
362*
363* Set up parameters with SLATB4 and generate a
364* test matrix with SLATMS.
365*
366 CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm,
367 $ mode, cndnum, dist )
368 rcondc = one / cndnum
369*
370 srnamt = 'SLATMS'
371 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
372 $ cndnum, anorm, kl, ku, 'Z', a, lda, work,
373 $ info )
374*
375* Check the error code from SLATMS.
376*
377 IF( info.NE.0 ) THEN
378 CALL alaerh( path, 'SLATMS', info, 0, ' ', n, n,
379 $ kl, ku, -1, imat, nfail, nerrs, nout )
380 GO TO 120
381 END IF
382*
383* For types 2, 3, and 4, zero one or more columns of
384* the matrix to test that INFO is returned correctly.
385*
386 izero = 0
387 IF( zerot ) THEN
388 IF( imat.EQ.2 ) THEN
389 izero = 1
390 ELSE IF( imat.EQ.3 ) THEN
391 izero = n
392 ELSE
393 izero = n / 2 + 1
394 END IF
395 ioff = ( izero-1 )*lda
396 IF( imat.LT.4 ) THEN
397 i1 = max( 1, ku+2-izero )
398 i2 = min( kl+ku+1, ku+1+( n-izero ) )
399 DO 20 i = i1, i2
400 a( ioff+i ) = zero
401 20 CONTINUE
402 ELSE
403 DO 40 j = izero, n
404 DO 30 i = max( 1, ku+2-j ),
405 $ min( kl+ku+1, ku+1+( n-j ) )
406 a( ioff+i ) = zero
407 30 CONTINUE
408 ioff = ioff + lda
409 40 CONTINUE
410 END IF
411 END IF
412*
413* Save a copy of the matrix A in ASAV.
414*
415 CALL slacpy( 'Full', kl+ku+1, n, a, lda, asav, lda )
416*
417 DO 110 iequed = 1, 4
418 equed = equeds( iequed )
419 IF( iequed.EQ.1 ) THEN
420 nfact = 3
421 ELSE
422 nfact = 1
423 END IF
424*
425 DO 100 ifact = 1, nfact
426 fact = facts( ifact )
427 prefac = lsame( fact, 'F' )
428 nofact = lsame( fact, 'N' )
429 equil = lsame( fact, 'E' )
430*
431 IF( zerot ) THEN
432 IF( prefac )
433 $ GO TO 100
434 rcondo = zero
435 rcondi = zero
436*
437 ELSE IF( .NOT.nofact ) THEN
438*
439* Compute the condition number for comparison
440* with the value returned by SGESVX (FACT =
441* 'N' reuses the condition number from the
442* previous iteration with FACT = 'F').
443*
444 CALL slacpy( 'Full', kl+ku+1, n, asav, lda,
445 $ afb( kl+1 ), ldafb )
446 IF( equil .OR. iequed.GT.1 ) THEN
447*
448* Compute row and column scale factors to
449* equilibrate the matrix A.
450*
451 CALL sgbequ( n, n, kl, ku, afb( kl+1 ),
452 $ ldafb, s, s( n+1 ), rowcnd,
453 $ colcnd, amax, info )
454 IF( info.EQ.0 .AND. n.GT.0 ) THEN
455 IF( lsame( equed, 'R' ) ) THEN
456 rowcnd = zero
457 colcnd = one
458 ELSE IF( lsame( equed, 'C' ) ) THEN
459 rowcnd = one
460 colcnd = zero
461 ELSE IF( lsame( equed, 'B' ) ) THEN
462 rowcnd = zero
463 colcnd = zero
464 END IF
465*
466* Equilibrate the matrix.
467*
468 CALL slaqgb( n, n, kl, ku, afb( kl+1 ),
469 $ ldafb, s, s( n+1 ),
470 $ rowcnd, colcnd, amax,
471 $ equed )
472 END IF
473 END IF
474*
475* Save the condition number of the
476* non-equilibrated system for use in SGET04.
477*
478 IF( equil ) THEN
479 roldo = rcondo
480 roldi = rcondi
481 END IF
482*
483* Compute the 1-norm and infinity-norm of A.
484*
485 anormo = slangb( '1', n, kl, ku, afb( kl+1 ),
486 $ ldafb, rwork )
487 anormi = slangb( 'I', n, kl, ku, afb( kl+1 ),
488 $ ldafb, rwork )
489*
490* Factor the matrix A.
491*
492 CALL sgbtrf( n, n, kl, ku, afb, ldafb, iwork,
493 $ info )
494*
495* Form the inverse of A.
496*
497 CALL slaset( 'Full', n, n, zero, one, work,
498 $ ldb )
499 srnamt = 'SGBTRS'
500 CALL sgbtrs( 'No transpose', n, kl, ku, n,
501 $ afb, ldafb, iwork, work, ldb,
502 $ info )
503*
504* Compute the 1-norm condition number of A.
505*
506 ainvnm = slange( '1', n, n, work, ldb,
507 $ rwork )
508 IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
509 rcondo = one
510 ELSE
511 rcondo = ( one / anormo ) / ainvnm
512 END IF
513*
514* Compute the infinity-norm condition number
515* of A.
516*
517 ainvnm = slange( 'I', n, n, work, ldb,
518 $ rwork )
519 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
520 rcondi = one
521 ELSE
522 rcondi = ( one / anormi ) / ainvnm
523 END IF
524 END IF
525*
526 DO 90 itran = 1, ntran
527*
528* Do for each value of TRANS.
529*
530 trans = transs( itran )
531 IF( itran.EQ.1 ) THEN
532 rcondc = rcondo
533 ELSE
534 rcondc = rcondi
535 END IF
536*
537* Restore the matrix A.
538*
539 CALL slacpy( 'Full', kl+ku+1, n, asav, lda,
540 $ a, lda )
541*
542* Form an exact solution and set the right hand
543* side.
544*
545 srnamt = 'SLARHS'
546 CALL slarhs( path, xtype, 'Full', trans, n,
547 $ n, kl, ku, nrhs, a, lda, xact,
548 $ ldb, b, ldb, iseed, info )
549 xtype = 'C'
550 CALL slacpy( 'Full', n, nrhs, b, ldb, bsav,
551 $ ldb )
552*
553 IF( nofact .AND. itran.EQ.1 ) THEN
554*
555* --- Test SGBSV ---
556*
557* Compute the LU factorization of the matrix
558* and solve the system.
559*
560 CALL slacpy( 'Full', kl+ku+1, n, a, lda,
561 $ afb( kl+1 ), ldafb )
562 CALL slacpy( 'Full', n, nrhs, b, ldb, x,
563 $ ldb )
564*
565 srnamt = 'SGBSV '
566 CALL sgbsv( n, kl, ku, nrhs, afb, ldafb,
567 $ iwork, x, ldb, info )
568*
569* Check error code from SGBSV .
570*
571 IF( info.NE.izero )
572 $ CALL alaerh( path, 'SGBSV ', info,
573 $ izero, ' ', n, n, kl, ku,
574 $ nrhs, imat, nfail, nerrs,
575 $ nout )
576*
577* Reconstruct matrix from factors and
578* compute residual.
579*
580 CALL sgbt01( n, n, kl, ku, a, lda, afb,
581 $ ldafb, iwork, work,
582 $ result( 1 ) )
583 nt = 1
584 IF( izero.EQ.0 ) THEN
585*
586* Compute residual of the computed
587* solution.
588*
589 CALL slacpy( 'Full', n, nrhs, b, ldb,
590 $ work, ldb )
591 CALL sgbt02( 'No transpose', n, n, kl,
592 $ ku, nrhs, a, lda, x, ldb,
593 $ work, ldb, rwork,
594 $ result( 2 ) )
595*
596* Check solution from generated exact
597* solution.
598*
599 CALL sget04( n, nrhs, x, ldb, xact,
600 $ ldb, rcondc, result( 3 ) )
601 nt = 3
602 END IF
603*
604* Print information about the tests that did
605* not pass the threshold.
606*
607 DO 50 k = 1, nt
608 IF( result( k ).GE.thresh ) THEN
609 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
610 $ CALL aladhd( nout, path )
611 WRITE( nout, fmt = 9997 )'SGBSV ',
612 $ n, kl, ku, imat, k, result( k )
613 nfail = nfail + 1
614 END IF
615 50 CONTINUE
616 nrun = nrun + nt
617 END IF
618*
619* --- Test SGBSVX ---
620*
621 IF( .NOT.prefac )
622 $ CALL slaset( 'Full', 2*kl+ku+1, n, zero,
623 $ zero, afb, ldafb )
624 CALL slaset( 'Full', n, nrhs, zero, zero, x,
625 $ ldb )
626 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
627*
628* Equilibrate the matrix if FACT = 'F' and
629* EQUED = 'R', 'C', or 'B'.
630*
631 CALL slaqgb( n, n, kl, ku, a, lda, s,
632 $ s( n+1 ), rowcnd, colcnd,
633 $ amax, equed )
634 END IF
635*
636* Solve the system and compute the condition
637* number and error bounds using SGBSVX.
638*
639 srnamt = 'SGBSVX'
640 CALL sgbsvx( fact, trans, n, kl, ku, nrhs, a,
641 $ lda, afb, ldafb, iwork, equed,
642 $ s, s( n+1 ), b, ldb, x, ldb,
643 $ rcond, rwork, rwork( nrhs+1 ),
644 $ work, iwork( n+1 ), info )
645*
646* Check the error code from SGBSVX.
647*
648 IF( info.NE.izero )
649 $ CALL alaerh( path, 'SGBSVX', info, izero,
650 $ fact // trans, n, n, kl, ku,
651 $ nrhs, imat, nfail, nerrs,
652 $ nout )
653*
654* Compare WORK(1) from SGBSVX with the computed
655* reciprocal pivot growth factor RPVGRW
656*
657 IF( info.NE.0 ) THEN
658 anrmpv = zero
659 DO 70 j = 1, info
660 DO 60 i = max( ku+2-j, 1 ),
661 $ min( n+ku+1-j, kl+ku+1 )
662 anrmpv = max( anrmpv,
663 $ abs( a( i+( j-1 )*lda ) ) )
664 60 CONTINUE
665 70 CONTINUE
666 rpvgrw = slantb( 'M', 'U', 'N', info,
667 $ min( info-1, kl+ku ),
668 $ afb( max( 1, kl+ku+2-info ) ),
669 $ ldafb, work )
670 IF( rpvgrw.EQ.zero ) THEN
671 rpvgrw = one
672 ELSE
673 rpvgrw = anrmpv / rpvgrw
674 END IF
675 ELSE
676 rpvgrw = slantb( 'M', 'U', 'N', n, kl+ku,
677 $ afb, ldafb, work )
678 IF( rpvgrw.EQ.zero ) THEN
679 rpvgrw = one
680 ELSE
681 rpvgrw = slangb( 'M', n, kl, ku, a,
682 $ lda, work ) / rpvgrw
683 END IF
684 END IF
685 result( 7 ) = abs( rpvgrw-work( 1 ) ) /
686 $ max( work( 1 ), rpvgrw ) /
687 $ slamch( 'E' )
688*
689 IF( .NOT.prefac ) THEN
690*
691* Reconstruct matrix from factors and
692* compute residual.
693*
694 CALL sgbt01( n, n, kl, ku, a, lda, afb,
695 $ ldafb, iwork, work,
696 $ result( 1 ) )
697 k1 = 1
698 ELSE
699 k1 = 2
700 END IF
701*
702 IF( info.EQ.0 ) THEN
703 trfcon = .false.
704*
705* Compute residual of the computed solution.
706*
707 CALL slacpy( 'Full', n, nrhs, bsav, ldb,
708 $ work, ldb )
709 CALL sgbt02( trans, n, n, kl, ku, nrhs,
710 $ asav, lda, x, ldb, work, ldb,
711 $ rwork( 2*nrhs+1 ),
712 $ result( 2 ) )
713*
714* Check solution from generated exact
715* solution.
716*
717 IF( nofact .OR. ( prefac .AND.
718 $ lsame( equed, 'N' ) ) ) THEN
719 CALL sget04( n, nrhs, x, ldb, xact,
720 $ ldb, rcondc, result( 3 ) )
721 ELSE
722 IF( itran.EQ.1 ) THEN
723 roldc = roldo
724 ELSE
725 roldc = roldi
726 END IF
727 CALL sget04( n, nrhs, x, ldb, xact,
728 $ ldb, roldc, result( 3 ) )
729 END IF
730*
731* Check the error bounds from iterative
732* refinement.
733*
734 CALL sgbt05( trans, n, kl, ku, nrhs, asav,
735 $ lda, b, ldb, x, ldb, xact,
736 $ ldb, rwork, rwork( nrhs+1 ),
737 $ result( 4 ) )
738 ELSE
739 trfcon = .true.
740 END IF
741*
742* Compare RCOND from SGBSVX with the computed
743* value in RCONDC.
744*
745 result( 6 ) = sget06( rcond, rcondc )
746*
747* Print information about the tests that did
748* not pass the threshold.
749*
750 IF( .NOT.trfcon ) THEN
751 DO 80 k = k1, ntests
752 IF( result( k ).GE.thresh ) THEN
753 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
754 $ CALL aladhd( nout, path )
755 IF( prefac ) THEN
756 WRITE( nout, fmt = 9995 )
757 $ 'SGBSVX', fact, trans, n, kl,
758 $ ku, equed, imat, k,
759 $ result( k )
760 ELSE
761 WRITE( nout, fmt = 9996 )
762 $ 'SGBSVX', fact, trans, n, kl,
763 $ ku, imat, k, result( k )
764 END IF
765 nfail = nfail + 1
766 END IF
767 80 CONTINUE
768 nrun = nrun + 7 - k1
769 ELSE
770 IF( result( 1 ).GE.thresh .AND. .NOT.
771 $ prefac ) THEN
772 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
773 $ CALL aladhd( nout, path )
774 IF( prefac ) THEN
775 WRITE( nout, fmt = 9995 )'SGBSVX',
776 $ fact, trans, n, kl, ku, equed,
777 $ imat, 1, result( 1 )
778 ELSE
779 WRITE( nout, fmt = 9996 )'SGBSVX',
780 $ fact, trans, n, kl, ku, imat, 1,
781 $ result( 1 )
782 END IF
783 nfail = nfail + 1
784 nrun = nrun + 1
785 END IF
786 IF( result( 6 ).GE.thresh ) THEN
787 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
788 $ CALL aladhd( nout, path )
789 IF( prefac ) THEN
790 WRITE( nout, fmt = 9995 )'SGBSVX',
791 $ fact, trans, n, kl, ku, equed,
792 $ imat, 6, result( 6 )
793 ELSE
794 WRITE( nout, fmt = 9996 )'SGBSVX',
795 $ fact, trans, n, kl, ku, imat, 6,
796 $ result( 6 )
797 END IF
798 nfail = nfail + 1
799 nrun = nrun + 1
800 END IF
801 IF( result( 7 ).GE.thresh ) THEN
802 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
803 $ CALL aladhd( nout, path )
804 IF( prefac ) THEN
805 WRITE( nout, fmt = 9995 )'SGBSVX',
806 $ fact, trans, n, kl, ku, equed,
807 $ imat, 7, result( 7 )
808 ELSE
809 WRITE( nout, fmt = 9996 )'SGBSVX',
810 $ fact, trans, n, kl, ku, imat, 7,
811 $ result( 7 )
812 END IF
813 nfail = nfail + 1
814 nrun = nrun + 1
815 END IF
816*
817 END IF
818*
819* --- Test SGBSVXX ---
820*
821* Restore the matrices A and B.
822*
823 CALL slacpy( 'Full', kl+ku+1, n, asav, lda, a,
824 $ lda )
825 CALL slacpy( 'Full', n, nrhs, bsav, ldb, b, ldb )
826
827 IF( .NOT.prefac )
828 $ CALL slaset( 'Full', 2*kl+ku+1, n, zero, zero,
829 $ afb, ldafb )
830 CALL slaset( 'Full', n, nrhs, zero, zero, x, ldb )
831 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
832*
833* Equilibrate the matrix if FACT = 'F' and
834* EQUED = 'R', 'C', or 'B'.
835*
836 CALL slaqgb( n, n, kl, ku, a, lda, s,
837 $ s( n+1 ), rowcnd, colcnd, amax, equed )
838 END IF
839*
840* Solve the system and compute the condition number
841* and error bounds using SGBSVXX.
842*
843 srnamt = 'SGBSVXX'
844 n_err_bnds = 3
845 CALL sgbsvxx( fact, trans, n, kl, ku, nrhs, a, lda,
846 $ afb, ldafb, iwork, equed, s, s( n+1 ), b, ldb,
847 $ x, ldb, rcond, rpvgrw_svxx, berr, n_err_bnds,
848 $ errbnds_n, errbnds_c, 0, zero, work,
849 $ iwork( n+1 ), info )
850
851* Check the error code from SGBSVXX.
852*
853 IF( info.EQ.n+1 ) GOTO 90
854 IF( info.NE.izero ) THEN
855 CALL alaerh( path, 'SGBSVXX', info, izero,
856 $ fact // trans, n, n, -1, -1, nrhs,
857 $ imat, nfail, nerrs, nout )
858 GOTO 90
859 END IF
860*
861* Compare rpvgrw_svxx from SGBSVXX with the computed
862* reciprocal pivot growth factor RPVGRW
863*
864
865 IF ( info .GT. 0 .AND. info .LT. n+1 ) THEN
866 rpvgrw = sla_gbrpvgrw(n, kl, ku, info, a, lda,
867 $ afb, ldafb )
868 ELSE
869 rpvgrw = sla_gbrpvgrw(n, kl, ku, n, a, lda,
870 $ afb, ldafb )
871 ENDIF
872
873 result( 7 ) = abs( rpvgrw-rpvgrw_svxx ) /
874 $ max( rpvgrw_svxx, rpvgrw ) /
875 $ slamch( 'E' )
876*
877 IF( .NOT.prefac ) THEN
878*
879* Reconstruct matrix from factors and compute
880* residual.
881*
882 CALL sgbt01( n, n, kl, ku, a, lda, afb, ldafb,
883 $ iwork, work,
884 $ result( 1 ) )
885 k1 = 1
886 ELSE
887 k1 = 2
888 END IF
889*
890 IF( info.EQ.0 ) THEN
891 trfcon = .false.
892*
893* Compute residual of the computed solution.
894*
895 CALL slacpy( 'Full', n, nrhs, bsav, ldb, work,
896 $ ldb )
897 CALL sgbt02( trans, n, n, kl, ku, nrhs, asav,
898 $ lda, x, ldb, work, ldb, rwork,
899 $ result( 2 ) )
900*
901* Check solution from generated exact solution.
902*
903 IF( nofact .OR. ( prefac .AND. lsame( equed,
904 $ 'N' ) ) ) THEN
905 CALL sget04( n, nrhs, x, ldb, xact, ldb,
906 $ rcondc, result( 3 ) )
907 ELSE
908 IF( itran.EQ.1 ) THEN
909 roldc = roldo
910 ELSE
911 roldc = roldi
912 END IF
913 CALL sget04( n, nrhs, x, ldb, xact, ldb,
914 $ roldc, result( 3 ) )
915 END IF
916 ELSE
917 trfcon = .true.
918 END IF
919*
920* Compare RCOND from SGBSVXX with the computed value
921* in RCONDC.
922*
923 result( 6 ) = sget06( rcond, rcondc )
924*
925* Print information about the tests that did not pass
926* the threshold.
927*
928 IF( .NOT.trfcon ) THEN
929 DO 45 k = k1, ntests
930 IF( result( k ).GE.thresh ) THEN
931 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
932 $ CALL aladhd( nout, path )
933 IF( prefac ) THEN
934 WRITE( nout, fmt = 9995 )'SGBSVXX',
935 $ fact, trans, n, kl, ku, equed,
936 $ imat, k, result( k )
937 ELSE
938 WRITE( nout, fmt = 9996 )'SGBSVXX',
939 $ fact, trans, n, kl, ku, imat, k,
940 $ result( k )
941 END IF
942 nfail = nfail + 1
943 END IF
944 45 CONTINUE
945 nrun = nrun + 7 - k1
946 ELSE
947 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
948 $ THEN
949 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
950 $ CALL aladhd( nout, path )
951 IF( prefac ) THEN
952 WRITE( nout, fmt = 9995 )'SGBSVXX', fact,
953 $ trans, n, kl, ku, equed, imat, 1,
954 $ result( 1 )
955 ELSE
956 WRITE( nout, fmt = 9996 )'SGBSVXX', fact,
957 $ trans, n, kl, ku, imat, 1,
958 $ result( 1 )
959 END IF
960 nfail = nfail + 1
961 nrun = nrun + 1
962 END IF
963 IF( result( 6 ).GE.thresh ) THEN
964 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
965 $ CALL aladhd( nout, path )
966 IF( prefac ) THEN
967 WRITE( nout, fmt = 9995 )'SGBSVXX', fact,
968 $ trans, n, kl, ku, equed, imat, 6,
969 $ result( 6 )
970 ELSE
971 WRITE( nout, fmt = 9996 )'SGBSVXX', fact,
972 $ trans, n, kl, ku, imat, 6,
973 $ result( 6 )
974 END IF
975 nfail = nfail + 1
976 nrun = nrun + 1
977 END IF
978 IF( result( 7 ).GE.thresh ) THEN
979 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
980 $ CALL aladhd( nout, path )
981 IF( prefac ) THEN
982 WRITE( nout, fmt = 9995 )'SGBSVXX', fact,
983 $ trans, n, kl, ku, equed, imat, 7,
984 $ result( 7 )
985 ELSE
986 WRITE( nout, fmt = 9996 )'SGBSVXX', fact,
987 $ trans, n, kl, ku, imat, 7,
988 $ result( 7 )
989 END IF
990 nfail = nfail + 1
991 nrun = nrun + 1
992 END IF
993
994 END IF
995*
996 90 CONTINUE
997 100 CONTINUE
998 110 CONTINUE
999 120 CONTINUE
1000 130 CONTINUE
1001 140 CONTINUE
1002 150 CONTINUE
1003*
1004* Print a summary of the results.
1005*
1006 CALL alasvm( path, nout, nfail, nrun, nerrs )
1007*
1008
1009* Test Error Bounds from SGBSVXX
1010
1011 CALL sebchvxx(thresh, path)
1012
1013 9999 FORMAT( ' *** In SDRVGB, LA=', i5, ' is too small for N=', i5,
1014 $ ', KU=', i5, ', KL=', i5, / ' ==> Increase LA to at least ',
1015 $ i5 )
1016 9998 FORMAT( ' *** In SDRVGB, LAFB=', i5, ' is too small for N=', i5,
1017 $ ', KU=', i5, ', KL=', i5, /
1018 $ ' ==> Increase LAFB to at least ', i5 )
1019 9997 FORMAT( 1x, a, ', N=', i5, ', KL=', i5, ', KU=', i5, ', type ',
1020 $ i1, ', test(', i1, ')=', g12.5 )
1021 9996 FORMAT( 1x, a, '( ''', a1, ''',''', a1, ''',', i5, ',', i5, ',',
1022 $ i5, ',...), type ', i1, ', test(', i1, ')=', g12.5 )
1023 9995 FORMAT( 1x, a, '( ''', a1, ''',''', a1, ''',', i5, ',', i5, ',',
1024 $ i5, ',...), EQUED=''', a1, ''', type ', i1, ', test(', i1,
1025 $ ')=', g12.5 )
1026*
1027 RETURN
1028*
1029* End of SDRVGBX
1030*
1031 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine slarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
SLARHS
Definition slarhs.f:205
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine sgbequ(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, info)
SGBEQU
Definition sgbequ.f:153
subroutine sgbsv(n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
SGBSV computes the solution to system of linear equations A * X = B for GB matrices (simple driver)
Definition sgbsv.f:162
subroutine sgbsvx(fact, trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, equed, r, c, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
SGBSVX computes the solution to system of linear equations A * X = B for GB matrices
Definition sgbsvx.f:368
subroutine sgbsvxx(fact, trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, equed, r, c, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
SGBSVXX computes the solution to system of linear equations A * X = B for GB matrices
Definition sgbsvxx.f:563
subroutine sgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
SGBTRF
Definition sgbtrf.f:144
subroutine sgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
SGBTRS
Definition sgbtrs.f:138
real function sla_gbrpvgrw(n, kl, ku, ncols, ab, ldab, afb, ldafb)
SLA_GBRPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a general banded matrix.
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slaqgb(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, equed)
SLAQGB scales a general band matrix, using row and column scaling factors computed by sgbequ.
Definition slaqgb.f:159
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sdrvgb(dotype, nn, nval, nrhs, thresh, tsterr, a, la, afb, lafb, asav, b, bsav, x, xact, s, work, rwork, iwork, nout)
SDRVGB
Definition sdrvgb.f:172
subroutine sebchvxx(thresh, path)
SEBCHVXX
Definition sebchvxx.f:96
subroutine serrvx(path, nunit)
SERRVX
Definition serrvx.f:55
subroutine sgbt01(m, n, kl, ku, a, lda, afac, ldafac, ipiv, work, resid)
SGBT01
Definition sgbt01.f:126
subroutine sgbt02(trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
SGBT02
Definition sgbt02.f:149
subroutine sgbt05(trans, n, kl, ku, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
SGBT05
Definition sgbt05.f:176
subroutine sget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
SGET04
Definition sget04.f:102
subroutine slatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
SLATB4
Definition slatb4.f:120
subroutine slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321