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