LAPACK 3.3.1
Linear Algebra PACKage
|
00001 SUBROUTINE DDRVGB( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, LA, 00002 $ AFB, LAFB, ASAV, B, BSAV, X, XACT, S, WORK, 00003 $ RWORK, IWORK, NOUT ) 00004 * 00005 * -- LAPACK test routine (version 3.2.2) -- 00006 * Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. 00007 * April 2009 00008 * 00009 * .. Scalar Arguments .. 00010 LOGICAL TSTERR 00011 INTEGER LA, LAFB, NN, NOUT, NRHS 00012 DOUBLE PRECISION THRESH 00013 * .. 00014 * .. Array Arguments .. 00015 LOGICAL DOTYPE( * ) 00016 INTEGER IWORK( * ), NVAL( * ) 00017 DOUBLE PRECISION A( * ), AFB( * ), ASAV( * ), B( * ), BSAV( * ), 00018 $ RWORK( * ), S( * ), WORK( * ), X( * ), 00019 $ XACT( * ) 00020 * .. 00021 * 00022 * Purpose 00023 * ======= 00024 * 00025 * DDRVGB tests the driver routines DGBSV, -SVX, and -SVXX. 00026 * 00027 * Note that this file is used only when the XBLAS are available, 00028 * otherwise ddrvgb.f defines this subroutine. 00029 * 00030 * Arguments 00031 * ========= 00032 * 00033 * DOTYPE (input) LOGICAL array, dimension (NTYPES) 00034 * The matrix types to be used for testing. Matrices of type j 00035 * (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) = 00036 * .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used. 00037 * 00038 * NN (input) INTEGER 00039 * The number of values of N contained in the vector NVAL. 00040 * 00041 * NVAL (input) INTEGER array, dimension (NN) 00042 * The values of the matrix column dimension N. 00043 * 00044 * NRHS (input) INTEGER 00045 * The number of right hand side vectors to be generated for 00046 * each linear system. 00047 * 00048 * THRESH (input) DOUBLE PRECISION 00049 * The threshold value for the test ratios. A result is 00050 * included in the output file if RESULT >= THRESH. To have 00051 * every test ratio printed, use THRESH = 0. 00052 * 00053 * TSTERR (input) LOGICAL 00054 * Flag that indicates whether error exits are to be tested. 00055 * 00056 * A (workspace) DOUBLE PRECISION array, dimension (LA) 00057 * 00058 * LA (input) INTEGER 00059 * The length of the array A. LA >= (2*NMAX-1)*NMAX 00060 * where NMAX is the largest entry in NVAL. 00061 * 00062 * AFB (workspace) DOUBLE PRECISION array, dimension (LAFB) 00063 * 00064 * LAFB (input) INTEGER 00065 * The length of the array AFB. LAFB >= (3*NMAX-2)*NMAX 00066 * where NMAX is the largest entry in NVAL. 00067 * 00068 * ASAV (workspace) DOUBLE PRECISION array, dimension (LA) 00069 * 00070 * B (workspace) DOUBLE PRECISION array, dimension (NMAX*NRHS) 00071 * 00072 * BSAV (workspace) DOUBLE PRECISION array, dimension (NMAX*NRHS) 00073 * 00074 * X (workspace) DOUBLE PRECISION array, dimension (NMAX*NRHS) 00075 * 00076 * XACT (workspace) DOUBLE PRECISION array, dimension (NMAX*NRHS) 00077 * 00078 * S (workspace) DOUBLE PRECISION array, dimension (2*NMAX) 00079 * 00080 * WORK (workspace) DOUBLE PRECISION array, dimension 00081 * (NMAX*max(3,NRHS,NMAX)) 00082 * 00083 * RWORK (workspace) DOUBLE PRECISION array, dimension 00084 * (max(NMAX,2*NRHS)) 00085 * 00086 * IWORK (workspace) INTEGER array, dimension (2*NMAX) 00087 * 00088 * NOUT (input) INTEGER 00089 * The unit number for output. 00090 * 00091 * ===================================================================== 00092 * 00093 * .. Parameters .. 00094 DOUBLE PRECISION ONE, ZERO 00095 PARAMETER ( ONE = 1.0D+0, ZERO = 0.0D+0 ) 00096 INTEGER NTYPES 00097 PARAMETER ( NTYPES = 8 ) 00098 INTEGER NTESTS 00099 PARAMETER ( NTESTS = 7 ) 00100 INTEGER NTRAN 00101 PARAMETER ( NTRAN = 3 ) 00102 * .. 00103 * .. Local Scalars .. 00104 LOGICAL EQUIL, NOFACT, PREFAC, TRFCON, ZEROT 00105 CHARACTER DIST, EQUED, FACT, TRANS, TYPE, XTYPE 00106 CHARACTER*3 PATH 00107 INTEGER I, I1, I2, IEQUED, IFACT, IKL, IKU, IMAT, IN, 00108 $ INFO, IOFF, ITRAN, IZERO, J, K, K1, KL, KU, 00109 $ LDA, LDAFB, LDB, MODE, N, NB, NBMIN, NERRS, 00110 $ NFACT, NFAIL, NIMAT, NKL, NKU, NRUN, NT, 00111 $ N_ERR_BNDS 00112 DOUBLE PRECISION AINVNM, AMAX, ANORM, ANORMI, ANORMO, ANRMPV, 00113 $ CNDNUM, COLCND, RCOND, RCONDC, RCONDI, RCONDO, 00114 $ ROLDC, ROLDI, ROLDO, ROWCND, RPVGRW, 00115 $ RPVGRW_SVXX 00116 * .. 00117 * .. Local Arrays .. 00118 CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN ) 00119 INTEGER ISEED( 4 ), ISEEDY( 4 ) 00120 DOUBLE PRECISION RESULT( NTESTS ), BERR( NRHS ), 00121 $ ERRBNDS_N( NRHS, 3 ), ERRBNDS_C( NRHS, 3 ) 00122 * .. 00123 * .. External Functions .. 00124 LOGICAL LSAME 00125 DOUBLE PRECISION DGET06, DLAMCH, DLANGB, DLANGE, DLANTB, 00126 $ DLA_GBRPVGRW 00127 EXTERNAL LSAME, DGET06, DLAMCH, DLANGB, DLANGE, DLANTB, 00128 $ DLA_GBRPVGRW 00129 * .. 00130 * .. External Subroutines .. 00131 EXTERNAL ALADHD, ALAERH, ALASVM, DERRVX, DGBEQU, DGBSV, 00132 $ DGBSVX, DGBT01, DGBT02, DGBT05, DGBTRF, DGBTRS, 00133 $ DGET04, DLACPY, DLAQGB, DLARHS, DLASET, DLATB4, 00134 $ DLATMS, XLAENV, DGBSVXX, DGBEQUB 00135 * .. 00136 * .. Intrinsic Functions .. 00137 INTRINSIC ABS, MAX, MIN 00138 * .. 00139 * .. Scalars in Common .. 00140 LOGICAL LERR, OK 00141 CHARACTER*32 SRNAMT 00142 INTEGER INFOT, NUNIT 00143 * .. 00144 * .. Common blocks .. 00145 COMMON / INFOC / INFOT, NUNIT, OK, LERR 00146 COMMON / SRNAMC / SRNAMT 00147 * .. 00148 * .. Data statements .. 00149 DATA ISEEDY / 1988, 1989, 1990, 1991 / 00150 DATA TRANSS / 'N', 'T', 'C' / 00151 DATA FACTS / 'F', 'N', 'E' / 00152 DATA EQUEDS / 'N', 'R', 'C', 'B' / 00153 * .. 00154 * .. Executable Statements .. 00155 * 00156 * Initialize constants and the random number seed. 00157 * 00158 PATH( 1: 1 ) = 'Double precision' 00159 PATH( 2: 3 ) = 'GB' 00160 NRUN = 0 00161 NFAIL = 0 00162 NERRS = 0 00163 DO 10 I = 1, 4 00164 ISEED( I ) = ISEEDY( I ) 00165 10 CONTINUE 00166 * 00167 * Test the error exits 00168 * 00169 IF( TSTERR ) 00170 $ CALL DERRVX( PATH, NOUT ) 00171 INFOT = 0 00172 * 00173 * Set the block size and minimum block size for testing. 00174 * 00175 NB = 1 00176 NBMIN = 2 00177 CALL XLAENV( 1, NB ) 00178 CALL XLAENV( 2, NBMIN ) 00179 * 00180 * Do for each value of N in NVAL 00181 * 00182 DO 150 IN = 1, NN 00183 N = NVAL( IN ) 00184 LDB = MAX( N, 1 ) 00185 XTYPE = 'N' 00186 * 00187 * Set limits on the number of loop iterations. 00188 * 00189 NKL = MAX( 1, MIN( N, 4 ) ) 00190 IF( N.EQ.0 ) 00191 $ NKL = 1 00192 NKU = NKL 00193 NIMAT = NTYPES 00194 IF( N.LE.0 ) 00195 $ NIMAT = 1 00196 * 00197 DO 140 IKL = 1, NKL 00198 * 00199 * Do for KL = 0, N-1, (3N-1)/4, and (N+1)/4. This order makes 00200 * it easier to skip redundant values for small values of N. 00201 * 00202 IF( IKL.EQ.1 ) THEN 00203 KL = 0 00204 ELSE IF( IKL.EQ.2 ) THEN 00205 KL = MAX( N-1, 0 ) 00206 ELSE IF( IKL.EQ.3 ) THEN 00207 KL = ( 3*N-1 ) / 4 00208 ELSE IF( IKL.EQ.4 ) THEN 00209 KL = ( N+1 ) / 4 00210 END IF 00211 DO 130 IKU = 1, NKU 00212 * 00213 * Do for KU = 0, N-1, (3N-1)/4, and (N+1)/4. This order 00214 * makes it easier to skip redundant values for small 00215 * values of N. 00216 * 00217 IF( IKU.EQ.1 ) THEN 00218 KU = 0 00219 ELSE IF( IKU.EQ.2 ) THEN 00220 KU = MAX( N-1, 0 ) 00221 ELSE IF( IKU.EQ.3 ) THEN 00222 KU = ( 3*N-1 ) / 4 00223 ELSE IF( IKU.EQ.4 ) THEN 00224 KU = ( N+1 ) / 4 00225 END IF 00226 * 00227 * Check that A and AFB are big enough to generate this 00228 * matrix. 00229 * 00230 LDA = KL + KU + 1 00231 LDAFB = 2*KL + KU + 1 00232 IF( LDA*N.GT.LA .OR. LDAFB*N.GT.LAFB ) THEN 00233 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00234 $ CALL ALADHD( NOUT, PATH ) 00235 IF( LDA*N.GT.LA ) THEN 00236 WRITE( NOUT, FMT = 9999 )LA, N, KL, KU, 00237 $ N*( KL+KU+1 ) 00238 NERRS = NERRS + 1 00239 END IF 00240 IF( LDAFB*N.GT.LAFB ) THEN 00241 WRITE( NOUT, FMT = 9998 )LAFB, N, KL, KU, 00242 $ N*( 2*KL+KU+1 ) 00243 NERRS = NERRS + 1 00244 END IF 00245 GO TO 130 00246 END IF 00247 * 00248 DO 120 IMAT = 1, NIMAT 00249 * 00250 * Do the tests only if DOTYPE( IMAT ) is true. 00251 * 00252 IF( .NOT.DOTYPE( IMAT ) ) 00253 $ GO TO 120 00254 * 00255 * Skip types 2, 3, or 4 if the matrix is too small. 00256 * 00257 ZEROT = IMAT.GE.2 .AND. IMAT.LE.4 00258 IF( ZEROT .AND. N.LT.IMAT-1 ) 00259 $ GO TO 120 00260 * 00261 * Set up parameters with DLATB4 and generate a 00262 * test matrix with DLATMS. 00263 * 00264 CALL DLATB4( PATH, IMAT, N, N, TYPE, KL, KU, ANORM, 00265 $ MODE, CNDNUM, DIST ) 00266 RCONDC = ONE / CNDNUM 00267 * 00268 SRNAMT = 'DLATMS' 00269 CALL DLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE, 00270 $ CNDNUM, ANORM, KL, KU, 'Z', A, LDA, WORK, 00271 $ INFO ) 00272 * 00273 * Check the error code from DLATMS. 00274 * 00275 IF( INFO.NE.0 ) THEN 00276 CALL ALAERH( PATH, 'DLATMS', INFO, 0, ' ', N, N, 00277 $ KL, KU, -1, IMAT, NFAIL, NERRS, NOUT ) 00278 GO TO 120 00279 END IF 00280 * 00281 * For types 2, 3, and 4, zero one or more columns of 00282 * the matrix to test that INFO is returned correctly. 00283 * 00284 IZERO = 0 00285 IF( ZEROT ) THEN 00286 IF( IMAT.EQ.2 ) THEN 00287 IZERO = 1 00288 ELSE IF( IMAT.EQ.3 ) THEN 00289 IZERO = N 00290 ELSE 00291 IZERO = N / 2 + 1 00292 END IF 00293 IOFF = ( IZERO-1 )*LDA 00294 IF( IMAT.LT.4 ) THEN 00295 I1 = MAX( 1, KU+2-IZERO ) 00296 I2 = MIN( KL+KU+1, KU+1+( N-IZERO ) ) 00297 DO 20 I = I1, I2 00298 A( IOFF+I ) = ZERO 00299 20 CONTINUE 00300 ELSE 00301 DO 40 J = IZERO, N 00302 DO 30 I = MAX( 1, KU+2-J ), 00303 $ MIN( KL+KU+1, KU+1+( N-J ) ) 00304 A( IOFF+I ) = ZERO 00305 30 CONTINUE 00306 IOFF = IOFF + LDA 00307 40 CONTINUE 00308 END IF 00309 END IF 00310 * 00311 * Save a copy of the matrix A in ASAV. 00312 * 00313 CALL DLACPY( 'Full', KL+KU+1, N, A, LDA, ASAV, LDA ) 00314 * 00315 DO 110 IEQUED = 1, 4 00316 EQUED = EQUEDS( IEQUED ) 00317 IF( IEQUED.EQ.1 ) THEN 00318 NFACT = 3 00319 ELSE 00320 NFACT = 1 00321 END IF 00322 * 00323 DO 100 IFACT = 1, NFACT 00324 FACT = FACTS( IFACT ) 00325 PREFAC = LSAME( FACT, 'F' ) 00326 NOFACT = LSAME( FACT, 'N' ) 00327 EQUIL = LSAME( FACT, 'E' ) 00328 * 00329 IF( ZEROT ) THEN 00330 IF( PREFAC ) 00331 $ GO TO 100 00332 RCONDO = ZERO 00333 RCONDI = ZERO 00334 * 00335 ELSE IF( .NOT.NOFACT ) THEN 00336 * 00337 * Compute the condition number for comparison 00338 * with the value returned by DGESVX (FACT = 00339 * 'N' reuses the condition number from the 00340 * previous iteration with FACT = 'F'). 00341 * 00342 CALL DLACPY( 'Full', KL+KU+1, N, ASAV, LDA, 00343 $ AFB( KL+1 ), LDAFB ) 00344 IF( EQUIL .OR. IEQUED.GT.1 ) THEN 00345 * 00346 * Compute row and column scale factors to 00347 * equilibrate the matrix A. 00348 * 00349 CALL DGBEQU( N, N, KL, KU, AFB( KL+1 ), 00350 $ LDAFB, S, S( N+1 ), ROWCND, 00351 $ COLCND, AMAX, INFO ) 00352 IF( INFO.EQ.0 .AND. N.GT.0 ) THEN 00353 IF( LSAME( EQUED, 'R' ) ) THEN 00354 ROWCND = ZERO 00355 COLCND = ONE 00356 ELSE IF( LSAME( EQUED, 'C' ) ) THEN 00357 ROWCND = ONE 00358 COLCND = ZERO 00359 ELSE IF( LSAME( EQUED, 'B' ) ) THEN 00360 ROWCND = ZERO 00361 COLCND = ZERO 00362 END IF 00363 * 00364 * Equilibrate the matrix. 00365 * 00366 CALL DLAQGB( N, N, KL, KU, AFB( KL+1 ), 00367 $ LDAFB, S, S( N+1 ), 00368 $ ROWCND, COLCND, AMAX, 00369 $ EQUED ) 00370 END IF 00371 END IF 00372 * 00373 * Save the condition number of the 00374 * non-equilibrated system for use in DGET04. 00375 * 00376 IF( EQUIL ) THEN 00377 ROLDO = RCONDO 00378 ROLDI = RCONDI 00379 END IF 00380 * 00381 * Compute the 1-norm and infinity-norm of A. 00382 * 00383 ANORMO = DLANGB( '1', N, KL, KU, AFB( KL+1 ), 00384 $ LDAFB, RWORK ) 00385 ANORMI = DLANGB( 'I', N, KL, KU, AFB( KL+1 ), 00386 $ LDAFB, RWORK ) 00387 * 00388 * Factor the matrix A. 00389 * 00390 CALL DGBTRF( N, N, KL, KU, AFB, LDAFB, IWORK, 00391 $ INFO ) 00392 * 00393 * Form the inverse of A. 00394 * 00395 CALL DLASET( 'Full', N, N, ZERO, ONE, WORK, 00396 $ LDB ) 00397 SRNAMT = 'DGBTRS' 00398 CALL DGBTRS( 'No transpose', N, KL, KU, N, 00399 $ AFB, LDAFB, IWORK, WORK, LDB, 00400 $ INFO ) 00401 * 00402 * Compute the 1-norm condition number of A. 00403 * 00404 AINVNM = DLANGE( '1', N, N, WORK, LDB, 00405 $ RWORK ) 00406 IF( ANORMO.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN 00407 RCONDO = ONE 00408 ELSE 00409 RCONDO = ( ONE / ANORMO ) / AINVNM 00410 END IF 00411 * 00412 * Compute the infinity-norm condition number 00413 * of A. 00414 * 00415 AINVNM = DLANGE( 'I', N, N, WORK, LDB, 00416 $ RWORK ) 00417 IF( ANORMI.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN 00418 RCONDI = ONE 00419 ELSE 00420 RCONDI = ( ONE / ANORMI ) / AINVNM 00421 END IF 00422 END IF 00423 * 00424 DO 90 ITRAN = 1, NTRAN 00425 * 00426 * Do for each value of TRANS. 00427 * 00428 TRANS = TRANSS( ITRAN ) 00429 IF( ITRAN.EQ.1 ) THEN 00430 RCONDC = RCONDO 00431 ELSE 00432 RCONDC = RCONDI 00433 END IF 00434 * 00435 * Restore the matrix A. 00436 * 00437 CALL DLACPY( 'Full', KL+KU+1, N, ASAV, LDA, 00438 $ A, LDA ) 00439 * 00440 * Form an exact solution and set the right hand 00441 * side. 00442 * 00443 SRNAMT = 'DLARHS' 00444 CALL DLARHS( PATH, XTYPE, 'Full', TRANS, N, 00445 $ N, KL, KU, NRHS, A, LDA, XACT, 00446 $ LDB, B, LDB, ISEED, INFO ) 00447 XTYPE = 'C' 00448 CALL DLACPY( 'Full', N, NRHS, B, LDB, BSAV, 00449 $ LDB ) 00450 * 00451 IF( NOFACT .AND. ITRAN.EQ.1 ) THEN 00452 * 00453 * --- Test DGBSV --- 00454 * 00455 * Compute the LU factorization of the matrix 00456 * and solve the system. 00457 * 00458 CALL DLACPY( 'Full', KL+KU+1, N, A, LDA, 00459 $ AFB( KL+1 ), LDAFB ) 00460 CALL DLACPY( 'Full', N, NRHS, B, LDB, X, 00461 $ LDB ) 00462 * 00463 SRNAMT = 'DGBSV ' 00464 CALL DGBSV( N, KL, KU, NRHS, AFB, LDAFB, 00465 $ IWORK, X, LDB, INFO ) 00466 * 00467 * Check error code from DGBSV . 00468 * 00469 IF( INFO.NE.IZERO ) 00470 $ CALL ALAERH( PATH, 'DGBSV ', INFO, 00471 $ IZERO, ' ', N, N, KL, KU, 00472 $ NRHS, IMAT, NFAIL, NERRS, 00473 $ NOUT ) 00474 * 00475 * Reconstruct matrix from factors and 00476 * compute residual. 00477 * 00478 CALL DGBT01( N, N, KL, KU, A, LDA, AFB, 00479 $ LDAFB, IWORK, WORK, 00480 $ RESULT( 1 ) ) 00481 NT = 1 00482 IF( IZERO.EQ.0 ) THEN 00483 * 00484 * Compute residual of the computed 00485 * solution. 00486 * 00487 CALL DLACPY( 'Full', N, NRHS, B, LDB, 00488 $ WORK, LDB ) 00489 CALL DGBT02( 'No transpose', N, N, KL, 00490 $ KU, NRHS, A, LDA, X, LDB, 00491 $ WORK, LDB, RESULT( 2 ) ) 00492 * 00493 * Check solution from generated exact 00494 * solution. 00495 * 00496 CALL DGET04( N, NRHS, X, LDB, XACT, 00497 $ LDB, RCONDC, RESULT( 3 ) ) 00498 NT = 3 00499 END IF 00500 * 00501 * Print information about the tests that did 00502 * not pass the threshold. 00503 * 00504 DO 50 K = 1, NT 00505 IF( RESULT( K ).GE.THRESH ) THEN 00506 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00507 $ CALL ALADHD( NOUT, PATH ) 00508 WRITE( NOUT, FMT = 9997 )'DGBSV ', 00509 $ N, KL, KU, IMAT, K, RESULT( K ) 00510 NFAIL = NFAIL + 1 00511 END IF 00512 50 CONTINUE 00513 NRUN = NRUN + NT 00514 END IF 00515 * 00516 * --- Test DGBSVX --- 00517 * 00518 IF( .NOT.PREFAC ) 00519 $ CALL DLASET( 'Full', 2*KL+KU+1, N, ZERO, 00520 $ ZERO, AFB, LDAFB ) 00521 CALL DLASET( 'Full', N, NRHS, ZERO, ZERO, X, 00522 $ LDB ) 00523 IF( IEQUED.GT.1 .AND. N.GT.0 ) THEN 00524 * 00525 * Equilibrate the matrix if FACT = 'F' and 00526 * EQUED = 'R', 'C', or 'B'. 00527 * 00528 CALL DLAQGB( N, N, KL, KU, A, LDA, S, 00529 $ S( N+1 ), ROWCND, COLCND, 00530 $ AMAX, EQUED ) 00531 END IF 00532 * 00533 * Solve the system and compute the condition 00534 * number and error bounds using DGBSVX. 00535 * 00536 SRNAMT = 'DGBSVX' 00537 CALL DGBSVX( FACT, TRANS, N, KL, KU, NRHS, A, 00538 $ LDA, AFB, LDAFB, IWORK, EQUED, 00539 $ S, S( N+1 ), B, LDB, X, LDB, 00540 $ RCOND, RWORK, RWORK( NRHS+1 ), 00541 $ WORK, IWORK( N+1 ), INFO ) 00542 * 00543 * Check the error code from DGBSVX. 00544 * 00545 IF( INFO.NE.IZERO ) 00546 $ CALL ALAERH( PATH, 'DGBSVX', INFO, IZERO, 00547 $ FACT // TRANS, N, N, KL, KU, 00548 $ NRHS, IMAT, NFAIL, NERRS, 00549 $ NOUT ) 00550 * 00551 * Compare WORK(1) from DGBSVX with the computed 00552 * reciprocal pivot growth factor RPVGRW 00553 * 00554 IF( INFO.NE.0 ) THEN 00555 ANRMPV = ZERO 00556 DO 70 J = 1, INFO 00557 DO 60 I = MAX( KU+2-J, 1 ), 00558 $ MIN( N+KU+1-J, KL+KU+1 ) 00559 ANRMPV = MAX( ANRMPV, 00560 $ ABS( A( I+( J-1 )*LDA ) ) ) 00561 60 CONTINUE 00562 70 CONTINUE 00563 RPVGRW = DLANTB( 'M', 'U', 'N', INFO, 00564 $ MIN( INFO-1, KL+KU ), 00565 $ AFB( MAX( 1, KL+KU+2-INFO ) ), 00566 $ LDAFB, WORK ) 00567 IF( RPVGRW.EQ.ZERO ) THEN 00568 RPVGRW = ONE 00569 ELSE 00570 RPVGRW = ANRMPV / RPVGRW 00571 END IF 00572 ELSE 00573 RPVGRW = DLANTB( 'M', 'U', 'N', N, KL+KU, 00574 $ AFB, LDAFB, WORK ) 00575 IF( RPVGRW.EQ.ZERO ) THEN 00576 RPVGRW = ONE 00577 ELSE 00578 RPVGRW = DLANGB( 'M', N, KL, KU, A, 00579 $ LDA, WORK ) / RPVGRW 00580 END IF 00581 END IF 00582 RESULT( 7 ) = ABS( RPVGRW-WORK( 1 ) ) / 00583 $ MAX( WORK( 1 ), RPVGRW ) / 00584 $ DLAMCH( 'E' ) 00585 * 00586 IF( .NOT.PREFAC ) THEN 00587 * 00588 * Reconstruct matrix from factors and 00589 * compute residual. 00590 * 00591 CALL DGBT01( N, N, KL, KU, A, LDA, AFB, 00592 $ LDAFB, IWORK, WORK, 00593 $ RESULT( 1 ) ) 00594 K1 = 1 00595 ELSE 00596 K1 = 2 00597 END IF 00598 * 00599 IF( INFO.EQ.0 ) THEN 00600 TRFCON = .FALSE. 00601 * 00602 * Compute residual of the computed solution. 00603 * 00604 CALL DLACPY( 'Full', N, NRHS, BSAV, LDB, 00605 $ WORK, LDB ) 00606 CALL DGBT02( TRANS, N, N, KL, KU, NRHS, 00607 $ ASAV, LDA, X, LDB, WORK, LDB, 00608 $ RESULT( 2 ) ) 00609 * 00610 * Check solution from generated exact 00611 * solution. 00612 * 00613 IF( NOFACT .OR. ( PREFAC .AND. 00614 $ LSAME( EQUED, 'N' ) ) ) THEN 00615 CALL DGET04( N, NRHS, X, LDB, XACT, 00616 $ LDB, RCONDC, RESULT( 3 ) ) 00617 ELSE 00618 IF( ITRAN.EQ.1 ) THEN 00619 ROLDC = ROLDO 00620 ELSE 00621 ROLDC = ROLDI 00622 END IF 00623 CALL DGET04( N, NRHS, X, LDB, XACT, 00624 $ LDB, ROLDC, RESULT( 3 ) ) 00625 END IF 00626 * 00627 * Check the error bounds from iterative 00628 * refinement. 00629 * 00630 CALL DGBT05( TRANS, N, KL, KU, NRHS, ASAV, 00631 $ LDA, B, LDB, X, LDB, XACT, 00632 $ LDB, RWORK, RWORK( NRHS+1 ), 00633 $ RESULT( 4 ) ) 00634 ELSE 00635 TRFCON = .TRUE. 00636 END IF 00637 * 00638 * Compare RCOND from DGBSVX with the computed 00639 * value in RCONDC. 00640 * 00641 RESULT( 6 ) = DGET06( RCOND, RCONDC ) 00642 * 00643 * Print information about the tests that did 00644 * not pass the threshold. 00645 * 00646 IF( .NOT.TRFCON ) THEN 00647 DO 80 K = K1, NTESTS 00648 IF( RESULT( K ).GE.THRESH ) THEN 00649 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00650 $ CALL ALADHD( NOUT, PATH ) 00651 IF( PREFAC ) THEN 00652 WRITE( NOUT, FMT = 9995 ) 00653 $ 'DGBSVX', FACT, TRANS, N, KL, 00654 $ KU, EQUED, IMAT, K, 00655 $ RESULT( K ) 00656 ELSE 00657 WRITE( NOUT, FMT = 9996 ) 00658 $ 'DGBSVX', FACT, TRANS, N, KL, 00659 $ KU, IMAT, K, RESULT( K ) 00660 END IF 00661 NFAIL = NFAIL + 1 00662 END IF 00663 80 CONTINUE 00664 NRUN = NRUN + 7 - K1 00665 ELSE 00666 IF( RESULT( 1 ).GE.THRESH .AND. .NOT. 00667 $ PREFAC ) THEN 00668 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00669 $ CALL ALADHD( NOUT, PATH ) 00670 IF( PREFAC ) THEN 00671 WRITE( NOUT, FMT = 9995 )'DGBSVX', 00672 $ FACT, TRANS, N, KL, KU, EQUED, 00673 $ IMAT, 1, RESULT( 1 ) 00674 ELSE 00675 WRITE( NOUT, FMT = 9996 )'DGBSVX', 00676 $ FACT, TRANS, N, KL, KU, IMAT, 1, 00677 $ RESULT( 1 ) 00678 END IF 00679 NFAIL = NFAIL + 1 00680 NRUN = NRUN + 1 00681 END IF 00682 IF( RESULT( 6 ).GE.THRESH ) THEN 00683 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00684 $ CALL ALADHD( NOUT, PATH ) 00685 IF( PREFAC ) THEN 00686 WRITE( NOUT, FMT = 9995 )'DGBSVX', 00687 $ FACT, TRANS, N, KL, KU, EQUED, 00688 $ IMAT, 6, RESULT( 6 ) 00689 ELSE 00690 WRITE( NOUT, FMT = 9996 )'DGBSVX', 00691 $ FACT, TRANS, N, KL, KU, IMAT, 6, 00692 $ RESULT( 6 ) 00693 END IF 00694 NFAIL = NFAIL + 1 00695 NRUN = NRUN + 1 00696 END IF 00697 IF( RESULT( 7 ).GE.THRESH ) THEN 00698 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00699 $ CALL ALADHD( NOUT, PATH ) 00700 IF( PREFAC ) THEN 00701 WRITE( NOUT, FMT = 9995 )'DGBSVX', 00702 $ FACT, TRANS, N, KL, KU, EQUED, 00703 $ IMAT, 7, RESULT( 7 ) 00704 ELSE 00705 WRITE( NOUT, FMT = 9996 )'DGBSVX', 00706 $ FACT, TRANS, N, KL, KU, IMAT, 7, 00707 $ RESULT( 7 ) 00708 END IF 00709 NFAIL = NFAIL + 1 00710 NRUN = NRUN + 1 00711 END IF 00712 * 00713 END IF 00714 * 00715 * --- Test DGBSVXX --- 00716 * 00717 * Restore the matrices A and B. 00718 * 00719 CALL DLACPY( 'Full', KL+KU+1, N, ASAV, LDA, A, 00720 $ LDA ) 00721 CALL DLACPY( 'Full', N, NRHS, BSAV, LDB, B, LDB ) 00722 00723 IF( .NOT.PREFAC ) 00724 $ CALL DLASET( 'Full', 2*KL+KU+1, N, ZERO, ZERO, 00725 $ AFB, LDAFB ) 00726 CALL DLASET( 'Full', N, NRHS, ZERO, ZERO, X, LDB ) 00727 IF( IEQUED.GT.1 .AND. N.GT.0 ) THEN 00728 * 00729 * Equilibrate the matrix if FACT = 'F' and 00730 * EQUED = 'R', 'C', or 'B'. 00731 * 00732 CALL DLAQGB( N, N, KL, KU, A, LDA, S, S( N+1 ), 00733 $ ROWCND, COLCND, AMAX, EQUED ) 00734 END IF 00735 * 00736 * Solve the system and compute the condition number 00737 * and error bounds using DGBSVXX. 00738 * 00739 SRNAMT = 'DGBSVXX' 00740 N_ERR_BNDS = 3 00741 CALL DGBSVXX( FACT, TRANS, N, KL, KU, NRHS, A, LDA, 00742 $ AFB, LDAFB, IWORK, EQUED, S, S( N+1 ), B, LDB, 00743 $ X, LDB, RCOND, RPVGRW_SVXX, BERR, N_ERR_BNDS, 00744 $ ERRBNDS_N, ERRBNDS_C, 0, ZERO, WORK, 00745 $ IWORK( N+1 ), INFO ) 00746 * 00747 * Check the error code from DGBSVXX. 00748 * 00749 IF( INFO.EQ.N+1 ) GOTO 90 00750 IF( INFO.NE.IZERO ) THEN 00751 CALL ALAERH( PATH, 'DGBSVXX', INFO, IZERO, 00752 $ FACT // TRANS, N, N, -1, -1, NRHS, 00753 $ IMAT, NFAIL, NERRS, NOUT ) 00754 GOTO 90 00755 END IF 00756 * 00757 * Compare rpvgrw_svxx from DGBSVXX with the computed 00758 * reciprocal pivot growth factor RPVGRW 00759 * 00760 00761 IF ( INFO .GT. 0 .AND. INFO .LT. N+1 ) THEN 00762 RPVGRW = DLA_GBRPVGRW(N, KL, KU, INFO, A, LDA, 00763 $ AFB, LDAFB) 00764 ELSE 00765 RPVGRW = DLA_GBRPVGRW(N, KL, KU, N, A, LDA, 00766 $ AFB, LDAFB) 00767 ENDIF 00768 00769 RESULT( 7 ) = ABS( RPVGRW-rpvgrw_svxx ) / 00770 $ MAX( rpvgrw_svxx, RPVGRW ) / 00771 $ DLAMCH( 'E' ) 00772 * 00773 IF( .NOT.PREFAC ) THEN 00774 * 00775 * Reconstruct matrix from factors and compute 00776 * residual. 00777 * 00778 CALL DGBT01( N, N, KL, KU, A, LDA, AFB, LDAFB, 00779 $ IWORK, WORK, RESULT( 1 ) ) 00780 K1 = 1 00781 ELSE 00782 K1 = 2 00783 END IF 00784 * 00785 IF( INFO.EQ.0 ) THEN 00786 TRFCON = .FALSE. 00787 * 00788 * Compute residual of the computed solution. 00789 * 00790 CALL DLACPY( 'Full', N, NRHS, BSAV, LDB, WORK, 00791 $ LDB ) 00792 CALL DGBT02( TRANS, N, N, KL, KU, NRHS, ASAV, 00793 $ LDA, X, LDB, WORK, LDB, 00794 $ WORK, RESULT( 2 ) ) 00795 * 00796 * Check solution from generated exact solution. 00797 * 00798 IF( NOFACT .OR. ( PREFAC .AND. LSAME( EQUED, 00799 $ 'N' ) ) ) THEN 00800 CALL DGET04( N, NRHS, X, LDB, XACT, LDB, 00801 $ RCONDC, RESULT( 3 ) ) 00802 ELSE 00803 IF( ITRAN.EQ.1 ) THEN 00804 ROLDC = ROLDO 00805 ELSE 00806 ROLDC = ROLDI 00807 END IF 00808 CALL DGET04( N, NRHS, X, LDB, XACT, LDB, 00809 $ ROLDC, RESULT( 3 ) ) 00810 END IF 00811 ELSE 00812 TRFCON = .TRUE. 00813 END IF 00814 * 00815 * Compare RCOND from DGBSVXX with the computed value 00816 * in RCONDC. 00817 * 00818 RESULT( 6 ) = DGET06( RCOND, RCONDC ) 00819 * 00820 * Print information about the tests that did not pass 00821 * the threshold. 00822 * 00823 IF( .NOT.TRFCON ) THEN 00824 DO 45 K = K1, NTESTS 00825 IF( RESULT( K ).GE.THRESH ) THEN 00826 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00827 $ CALL ALADHD( NOUT, PATH ) 00828 IF( PREFAC ) THEN 00829 WRITE( NOUT, FMT = 9995 )'DGBSVXX', 00830 $ FACT, TRANS, N, KL, KU, EQUED, 00831 $ IMAT, K, RESULT( K ) 00832 ELSE 00833 WRITE( NOUT, FMT = 9996 )'DGBSVXX', 00834 $ FACT, TRANS, N, KL, KU, IMAT, K, 00835 $ RESULT( K ) 00836 END IF 00837 NFAIL = NFAIL + 1 00838 END IF 00839 45 CONTINUE 00840 NRUN = NRUN + 7 - K1 00841 ELSE 00842 IF( RESULT( 1 ).GE.THRESH .AND. .NOT.PREFAC ) 00843 $ THEN 00844 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00845 $ CALL ALADHD( NOUT, PATH ) 00846 IF( PREFAC ) THEN 00847 WRITE( NOUT, FMT = 9995 )'DGBSVXX', FACT, 00848 $ TRANS, N, KL, KU, EQUED, IMAT, 1, 00849 $ RESULT( 1 ) 00850 ELSE 00851 WRITE( NOUT, FMT = 9996 )'DGBSVXX', FACT, 00852 $ TRANS, N, KL, KU, IMAT, 1, 00853 $ RESULT( 1 ) 00854 END IF 00855 NFAIL = NFAIL + 1 00856 NRUN = NRUN + 1 00857 END IF 00858 IF( RESULT( 6 ).GE.THRESH ) THEN 00859 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00860 $ CALL ALADHD( NOUT, PATH ) 00861 IF( PREFAC ) THEN 00862 WRITE( NOUT, FMT = 9995 )'DGBSVXX', FACT, 00863 $ TRANS, N, KL, KU, EQUED, IMAT, 6, 00864 $ RESULT( 6 ) 00865 ELSE 00866 WRITE( NOUT, FMT = 9996 )'DGBSVXX', FACT, 00867 $ TRANS, N, KL, KU, IMAT, 6, 00868 $ RESULT( 6 ) 00869 END IF 00870 NFAIL = NFAIL + 1 00871 NRUN = NRUN + 1 00872 END IF 00873 IF( RESULT( 7 ).GE.THRESH ) THEN 00874 IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 ) 00875 $ CALL ALADHD( NOUT, PATH ) 00876 IF( PREFAC ) THEN 00877 WRITE( NOUT, FMT = 9995 )'DGBSVXX', FACT, 00878 $ TRANS, N, KL, KU, EQUED, IMAT, 7, 00879 $ RESULT( 7 ) 00880 ELSE 00881 WRITE( NOUT, FMT = 9996 )'DGBSVXX', FACT, 00882 $ TRANS, N, KL, KU, IMAT, 7, 00883 $ RESULT( 7 ) 00884 END IF 00885 NFAIL = NFAIL + 1 00886 NRUN = NRUN + 1 00887 END IF 00888 * 00889 END IF 00890 90 CONTINUE 00891 100 CONTINUE 00892 110 CONTINUE 00893 120 CONTINUE 00894 130 CONTINUE 00895 140 CONTINUE 00896 150 CONTINUE 00897 * 00898 * Print a summary of the results. 00899 * 00900 CALL ALASVM( PATH, NOUT, NFAIL, NRUN, NERRS ) 00901 00902 * Test Error Bounds from DGBSVXX 00903 00904 CALL DEBCHVXX(THRESH, PATH) 00905 00906 9999 FORMAT( ' *** In DDRVGB, LA=', I5, ' is too small for N=', I5, 00907 $ ', KU=', I5, ', KL=', I5, / ' ==> Increase LA to at least ', 00908 $ I5 ) 00909 9998 FORMAT( ' *** In DDRVGB, LAFB=', I5, ' is too small for N=', I5, 00910 $ ', KU=', I5, ', KL=', I5, / 00911 $ ' ==> Increase LAFB to at least ', I5 ) 00912 9997 FORMAT( 1X, A, ', N=', I5, ', KL=', I5, ', KU=', I5, ', type ', 00913 $ I1, ', test(', I1, ')=', G12.5 ) 00914 9996 FORMAT( 1X, A, '( ''', A1, ''',''', A1, ''',', I5, ',', I5, ',', 00915 $ I5, ',...), type ', I1, ', test(', I1, ')=', G12.5 ) 00916 9995 FORMAT( 1X, A, '( ''', A1, ''',''', A1, ''',', I5, ',', I5, ',', 00917 $ I5, ',...), EQUED=''', A1, ''', type ', I1, ', test(', I1, 00918 $ ')=', G12.5 ) 00919 * 00920 RETURN 00921 * 00922 * End of DDRVGB 00923 * 00924 END