LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ddrvgb ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
double precision  THRESH,
logical  TSTERR,
double precision, dimension( * )  A,
integer  LA,
double precision, dimension( * )  AFB,
integer  LAFB,
double precision, dimension( * )  ASAV,
double precision, dimension( * )  B,
double precision, dimension( * )  BSAV,
double precision, dimension( * )  X,
double precision, dimension( * )  XACT,
double precision, dimension( * )  S,
double precision, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

DDRVGB

DDRVGBX

Purpose:
 DDRVGB tests the driver routines DGBSV and -SVX.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix column dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[out]A
          A is DOUBLE PRECISION array, dimension (LA)
[in]LA
          LA is INTEGER
          The length of the array A.  LA >= (2*NMAX-1)*NMAX
          where NMAX is the largest entry in NVAL.
[out]AFB
          AFB is DOUBLE PRECISION array, dimension (LAFB)
[in]LAFB
          LAFB is INTEGER
          The length of the array AFB.  LAFB >= (3*NMAX-2)*NMAX
          where NMAX is the largest entry in NVAL.
[out]ASAV
          ASAV is DOUBLE PRECISION array, dimension (LA)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]S
          S is DOUBLE PRECISION array, dimension (2*NMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NRHS,NMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(NMAX,2*NRHS))
[out]IWORK
          IWORK is INTEGER array, dimension (2*NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015
Purpose:
 DDRVGB tests the driver routines DGBSV, -SVX, and -SVXX.

 Note that this file is used only when the XBLAS are available,
 otherwise ddrvgb.f defines this subroutine.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix column dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[out]A
          A is DOUBLE PRECISION array, dimension (LA)
[in]LA
          LA is INTEGER
          The length of the array A.  LA >= (2*NMAX-1)*NMAX
          where NMAX is the largest entry in NVAL.
[out]AFB
          AFB is DOUBLE PRECISION array, dimension (LAFB)
[in]LAFB
          LAFB is INTEGER
          The length of the array AFB.  LAFB >= (3*NMAX-2)*NMAX
          where NMAX is the largest entry in NVAL.
[out]ASAV
          ASAV is DOUBLE PRECISION array, dimension (LA)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]S
          S is DOUBLE PRECISION array, dimension (2*NMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NRHS,NMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(NMAX,2*NRHS))
[out]IWORK
          IWORK is INTEGER array, dimension (2*NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 174 of file ddrvgb.f.

174 *
175 * -- LAPACK test routine (version 3.6.0) --
176 * -- LAPACK is a software package provided by Univ. of Tennessee, --
177 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178 * November 2015
179 *
180 * .. Scalar Arguments ..
181  LOGICAL tsterr
182  INTEGER la, lafb, nn, nout, nrhs
183  DOUBLE PRECISION thresh
184 * ..
185 * .. Array Arguments ..
186  LOGICAL dotype( * )
187  INTEGER iwork( * ), nval( * )
188  DOUBLE PRECISION a( * ), afb( * ), asav( * ), b( * ), bsav( * ),
189  $ rwork( * ), s( * ), work( * ), x( * ),
190  $ xact( * )
191 * ..
192 *
193 * =====================================================================
194 *
195 * .. Parameters ..
196  DOUBLE PRECISION one, zero
197  parameter ( one = 1.0d+0, zero = 0.0d+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  DOUBLE PRECISION ainvnm, amax, anorm, anormi, anormo, anrmpv,
214  $ cndnum, colcnd, rcond, rcondc, rcondi, rcondo,
215  $ roldc, roldi, roldo, rowcnd, rpvgrw
216 * ..
217 * .. Local Arrays ..
218  CHARACTER equeds( 4 ), facts( 3 ), transs( ntran )
219  INTEGER iseed( 4 ), iseedy( 4 )
220  DOUBLE PRECISION result( ntests )
221 * ..
222 * .. External Functions ..
223  LOGICAL lsame
224  DOUBLE PRECISION dget06, dlamch, dlangb, dlange, dlantb
225  EXTERNAL lsame, dget06, dlamch, dlangb, dlange, dlantb
226 * ..
227 * .. External Subroutines ..
228  EXTERNAL aladhd, alaerh, alasvm, derrvx, dgbequ, dgbsv,
231  $ dlatms, xlaenv
232 * ..
233 * .. Intrinsic Functions ..
234  INTRINSIC abs, max, min
235 * ..
236 * .. Scalars in Common ..
237  LOGICAL lerr, ok
238  CHARACTER*32 srnamt
239  INTEGER infot, nunit
240 * ..
241 * .. Common blocks ..
242  COMMON / infoc / infot, nunit, ok, lerr
243  COMMON / srnamc / srnamt
244 * ..
245 * .. Data statements ..
246  DATA iseedy / 1988, 1989, 1990, 1991 /
247  DATA transs / 'N', 'T', 'C' /
248  DATA facts / 'F', 'N', 'E' /
249  DATA equeds / 'N', 'R', 'C', 'B' /
250 * ..
251 * .. Executable Statements ..
252 *
253 * Initialize constants and the random number seed.
254 *
255  path( 1: 1 ) = 'Double precision'
256  path( 2: 3 ) = 'GB'
257  nrun = 0
258  nfail = 0
259  nerrs = 0
260  DO 10 i = 1, 4
261  iseed( i ) = iseedy( i )
262  10 CONTINUE
263 *
264 * Test the error exits
265 *
266  IF( tsterr )
267  $ CALL derrvx( path, nout )
268  infot = 0
269 *
270 * Set the block size and minimum block size for testing.
271 *
272  nb = 1
273  nbmin = 2
274  CALL xlaenv( 1, nb )
275  CALL xlaenv( 2, nbmin )
276 *
277 * Do for each value of N in NVAL
278 *
279  DO 150 in = 1, nn
280  n = nval( in )
281  ldb = max( n, 1 )
282  xtype = 'N'
283 *
284 * Set limits on the number of loop iterations.
285 *
286  nkl = max( 1, min( n, 4 ) )
287  IF( n.EQ.0 )
288  $ nkl = 1
289  nku = nkl
290  nimat = ntypes
291  IF( n.LE.0 )
292  $ nimat = 1
293 *
294  DO 140 ikl = 1, nkl
295 *
296 * Do for KL = 0, N-1, (3N-1)/4, and (N+1)/4. This order makes
297 * it easier to skip redundant values for small values of N.
298 *
299  IF( ikl.EQ.1 ) THEN
300  kl = 0
301  ELSE IF( ikl.EQ.2 ) THEN
302  kl = max( n-1, 0 )
303  ELSE IF( ikl.EQ.3 ) THEN
304  kl = ( 3*n-1 ) / 4
305  ELSE IF( ikl.EQ.4 ) THEN
306  kl = ( n+1 ) / 4
307  END IF
308  DO 130 iku = 1, nku
309 *
310 * Do for KU = 0, N-1, (3N-1)/4, and (N+1)/4. This order
311 * makes it easier to skip redundant values for small
312 * values of N.
313 *
314  IF( iku.EQ.1 ) THEN
315  ku = 0
316  ELSE IF( iku.EQ.2 ) THEN
317  ku = max( n-1, 0 )
318  ELSE IF( iku.EQ.3 ) THEN
319  ku = ( 3*n-1 ) / 4
320  ELSE IF( iku.EQ.4 ) THEN
321  ku = ( n+1 ) / 4
322  END IF
323 *
324 * Check that A and AFB are big enough to generate this
325 * matrix.
326 *
327  lda = kl + ku + 1
328  ldafb = 2*kl + ku + 1
329  IF( lda*n.GT.la .OR. ldafb*n.GT.lafb ) THEN
330  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
331  $ CALL aladhd( nout, path )
332  IF( lda*n.GT.la ) THEN
333  WRITE( nout, fmt = 9999 )la, n, kl, ku,
334  $ n*( kl+ku+1 )
335  nerrs = nerrs + 1
336  END IF
337  IF( ldafb*n.GT.lafb ) THEN
338  WRITE( nout, fmt = 9998 )lafb, n, kl, ku,
339  $ n*( 2*kl+ku+1 )
340  nerrs = nerrs + 1
341  END IF
342  GO TO 130
343  END IF
344 *
345  DO 120 imat = 1, nimat
346 *
347 * Do the tests only if DOTYPE( IMAT ) is true.
348 *
349  IF( .NOT.dotype( imat ) )
350  $ GO TO 120
351 *
352 * Skip types 2, 3, or 4 if the matrix is too small.
353 *
354  zerot = imat.GE.2 .AND. imat.LE.4
355  IF( zerot .AND. n.LT.imat-1 )
356  $ GO TO 120
357 *
358 * Set up parameters with DLATB4 and generate a
359 * test matrix with DLATMS.
360 *
361  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm,
362  $ mode, cndnum, dist )
363  rcondc = one / cndnum
364 *
365  srnamt = 'DLATMS'
366  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
367  $ cndnum, anorm, kl, ku, 'Z', a, lda, work,
368  $ info )
369 *
370 * Check the error code from DLATMS.
371 *
372  IF( info.NE.0 ) THEN
373  CALL alaerh( path, 'DLATMS', info, 0, ' ', n, n,
374  $ kl, ku, -1, imat, nfail, nerrs, nout )
375  GO TO 120
376  END IF
377 *
378 * For types 2, 3, and 4, zero one or more columns of
379 * the matrix to test that INFO is returned correctly.
380 *
381  izero = 0
382  IF( zerot ) THEN
383  IF( imat.EQ.2 ) THEN
384  izero = 1
385  ELSE IF( imat.EQ.3 ) THEN
386  izero = n
387  ELSE
388  izero = n / 2 + 1
389  END IF
390  ioff = ( izero-1 )*lda
391  IF( imat.LT.4 ) THEN
392  i1 = max( 1, ku+2-izero )
393  i2 = min( kl+ku+1, ku+1+( n-izero ) )
394  DO 20 i = i1, i2
395  a( ioff+i ) = zero
396  20 CONTINUE
397  ELSE
398  DO 40 j = izero, n
399  DO 30 i = max( 1, ku+2-j ),
400  $ min( kl+ku+1, ku+1+( n-j ) )
401  a( ioff+i ) = zero
402  30 CONTINUE
403  ioff = ioff + lda
404  40 CONTINUE
405  END IF
406  END IF
407 *
408 * Save a copy of the matrix A in ASAV.
409 *
410  CALL dlacpy( 'Full', kl+ku+1, n, a, lda, asav, lda )
411 *
412  DO 110 iequed = 1, 4
413  equed = equeds( iequed )
414  IF( iequed.EQ.1 ) THEN
415  nfact = 3
416  ELSE
417  nfact = 1
418  END IF
419 *
420  DO 100 ifact = 1, nfact
421  fact = facts( ifact )
422  prefac = lsame( fact, 'F' )
423  nofact = lsame( fact, 'N' )
424  equil = lsame( fact, 'E' )
425 *
426  IF( zerot ) THEN
427  IF( prefac )
428  $ GO TO 100
429  rcondo = zero
430  rcondi = zero
431 *
432  ELSE IF( .NOT.nofact ) THEN
433 *
434 * Compute the condition number for comparison
435 * with the value returned by DGESVX (FACT =
436 * 'N' reuses the condition number from the
437 * previous iteration with FACT = 'F').
438 *
439  CALL dlacpy( 'Full', kl+ku+1, n, asav, lda,
440  $ afb( kl+1 ), ldafb )
441  IF( equil .OR. iequed.GT.1 ) THEN
442 *
443 * Compute row and column scale factors to
444 * equilibrate the matrix A.
445 *
446  CALL dgbequ( n, n, kl, ku, afb( kl+1 ),
447  $ ldafb, s, s( n+1 ), rowcnd,
448  $ colcnd, amax, info )
449  IF( info.EQ.0 .AND. n.GT.0 ) THEN
450  IF( lsame( equed, 'R' ) ) THEN
451  rowcnd = zero
452  colcnd = one
453  ELSE IF( lsame( equed, 'C' ) ) THEN
454  rowcnd = one
455  colcnd = zero
456  ELSE IF( lsame( equed, 'B' ) ) THEN
457  rowcnd = zero
458  colcnd = zero
459  END IF
460 *
461 * Equilibrate the matrix.
462 *
463  CALL dlaqgb( n, n, kl, ku, afb( kl+1 ),
464  $ ldafb, s, s( n+1 ),
465  $ rowcnd, colcnd, amax,
466  $ equed )
467  END IF
468  END IF
469 *
470 * Save the condition number of the
471 * non-equilibrated system for use in DGET04.
472 *
473  IF( equil ) THEN
474  roldo = rcondo
475  roldi = rcondi
476  END IF
477 *
478 * Compute the 1-norm and infinity-norm of A.
479 *
480  anormo = dlangb( '1', n, kl, ku, afb( kl+1 ),
481  $ ldafb, rwork )
482  anormi = dlangb( 'I', n, kl, ku, afb( kl+1 ),
483  $ ldafb, rwork )
484 *
485 * Factor the matrix A.
486 *
487  CALL dgbtrf( n, n, kl, ku, afb, ldafb, iwork,
488  $ info )
489 *
490 * Form the inverse of A.
491 *
492  CALL dlaset( 'Full', n, n, zero, one, work,
493  $ ldb )
494  srnamt = 'DGBTRS'
495  CALL dgbtrs( 'No transpose', n, kl, ku, n,
496  $ afb, ldafb, iwork, work, ldb,
497  $ info )
498 *
499 * Compute the 1-norm condition number of A.
500 *
501  ainvnm = dlange( '1', n, n, work, ldb,
502  $ rwork )
503  IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
504  rcondo = one
505  ELSE
506  rcondo = ( one / anormo ) / ainvnm
507  END IF
508 *
509 * Compute the infinity-norm condition number
510 * of A.
511 *
512  ainvnm = dlange( 'I', n, n, work, ldb,
513  $ rwork )
514  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
515  rcondi = one
516  ELSE
517  rcondi = ( one / anormi ) / ainvnm
518  END IF
519  END IF
520 *
521  DO 90 itran = 1, ntran
522 *
523 * Do for each value of TRANS.
524 *
525  trans = transs( itran )
526  IF( itran.EQ.1 ) THEN
527  rcondc = rcondo
528  ELSE
529  rcondc = rcondi
530  END IF
531 *
532 * Restore the matrix A.
533 *
534  CALL dlacpy( 'Full', kl+ku+1, n, asav, lda,
535  $ a, lda )
536 *
537 * Form an exact solution and set the right hand
538 * side.
539 *
540  srnamt = 'DLARHS'
541  CALL dlarhs( path, xtype, 'Full', trans, n,
542  $ n, kl, ku, nrhs, a, lda, xact,
543  $ ldb, b, ldb, iseed, info )
544  xtype = 'C'
545  CALL dlacpy( 'Full', n, nrhs, b, ldb, bsav,
546  $ ldb )
547 *
548  IF( nofact .AND. itran.EQ.1 ) THEN
549 *
550 * --- Test DGBSV ---
551 *
552 * Compute the LU factorization of the matrix
553 * and solve the system.
554 *
555  CALL dlacpy( 'Full', kl+ku+1, n, a, lda,
556  $ afb( kl+1 ), ldafb )
557  CALL dlacpy( 'Full', n, nrhs, b, ldb, x,
558  $ ldb )
559 *
560  srnamt = 'DGBSV '
561  CALL dgbsv( n, kl, ku, nrhs, afb, ldafb,
562  $ iwork, x, ldb, info )
563 *
564 * Check error code from DGBSV .
565 *
566  IF( info.NE.izero )
567  $ CALL alaerh( path, 'DGBSV ', info,
568  $ izero, ' ', n, n, kl, ku,
569  $ nrhs, imat, nfail, nerrs,
570  $ nout )
571 *
572 * Reconstruct matrix from factors and
573 * compute residual.
574 *
575  CALL dgbt01( n, n, kl, ku, a, lda, afb,
576  $ ldafb, iwork, work,
577  $ result( 1 ) )
578  nt = 1
579  IF( izero.EQ.0 ) THEN
580 *
581 * Compute residual of the computed
582 * solution.
583 *
584  CALL dlacpy( 'Full', n, nrhs, b, ldb,
585  $ work, ldb )
586  CALL dgbt02( 'No transpose', n, n, kl,
587  $ ku, nrhs, a, lda, x, ldb,
588  $ work, ldb, result( 2 ) )
589 *
590 * Check solution from generated exact
591 * solution.
592 *
593  CALL dget04( n, nrhs, x, ldb, xact,
594  $ ldb, rcondc, result( 3 ) )
595  nt = 3
596  END IF
597 *
598 * Print information about the tests that did
599 * not pass the threshold.
600 *
601  DO 50 k = 1, nt
602  IF( result( k ).GE.thresh ) THEN
603  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
604  $ CALL aladhd( nout, path )
605  WRITE( nout, fmt = 9997 )'DGBSV ',
606  $ n, kl, ku, imat, k, result( k )
607  nfail = nfail + 1
608  END IF
609  50 CONTINUE
610  nrun = nrun + nt
611  END IF
612 *
613 * --- Test DGBSVX ---
614 *
615  IF( .NOT.prefac )
616  $ CALL dlaset( 'Full', 2*kl+ku+1, n, zero,
617  $ zero, afb, ldafb )
618  CALL dlaset( 'Full', n, nrhs, zero, zero, x,
619  $ ldb )
620  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
621 *
622 * Equilibrate the matrix if FACT = 'F' and
623 * EQUED = 'R', 'C', or 'B'.
624 *
625  CALL dlaqgb( n, n, kl, ku, a, lda, s,
626  $ s( n+1 ), rowcnd, colcnd,
627  $ amax, equed )
628  END IF
629 *
630 * Solve the system and compute the condition
631 * number and error bounds using DGBSVX.
632 *
633  srnamt = 'DGBSVX'
634  CALL dgbsvx( fact, trans, n, kl, ku, nrhs, a,
635  $ lda, afb, ldafb, iwork, equed,
636  $ s, s( n+1 ), b, ldb, x, ldb,
637  $ rcond, rwork, rwork( nrhs+1 ),
638  $ work, iwork( n+1 ), info )
639 *
640 * Check the error code from DGBSVX.
641 *
642  IF( info.NE.izero )
643  $ CALL alaerh( path, 'DGBSVX', info, izero,
644  $ fact // trans, n, n, kl, ku,
645  $ nrhs, imat, nfail, nerrs,
646  $ nout )
647 *
648 * Compare WORK(1) from DGBSVX with the computed
649 * reciprocal pivot growth factor RPVGRW
650 *
651  IF( info.NE.0 .AND. info.LE.n) THEN
652  anrmpv = zero
653  DO 70 j = 1, info
654  DO 60 i = max( ku+2-j, 1 ),
655  $ min( n+ku+1-j, kl+ku+1 )
656  anrmpv = max( anrmpv,
657  $ abs( a( i+( j-1 )*lda ) ) )
658  60 CONTINUE
659  70 CONTINUE
660  rpvgrw = dlantb( 'M', 'U', 'N', info,
661  $ min( info-1, kl+ku ),
662  $ afb( max( 1, kl+ku+2-info ) ),
663  $ ldafb, work )
664  IF( rpvgrw.EQ.zero ) THEN
665  rpvgrw = one
666  ELSE
667  rpvgrw = anrmpv / rpvgrw
668  END IF
669  ELSE
670  rpvgrw = dlantb( 'M', 'U', 'N', n, kl+ku,
671  $ afb, ldafb, work )
672  IF( rpvgrw.EQ.zero ) THEN
673  rpvgrw = one
674  ELSE
675  rpvgrw = dlangb( 'M', n, kl, ku, a,
676  $ lda, work ) / rpvgrw
677  END IF
678  END IF
679  result( 7 ) = abs( rpvgrw-work( 1 ) ) /
680  $ max( work( 1 ), rpvgrw ) /
681  $ dlamch( 'E' )
682 *
683  IF( .NOT.prefac ) THEN
684 *
685 * Reconstruct matrix from factors and
686 * compute residual.
687 *
688  CALL dgbt01( n, n, kl, ku, a, lda, afb,
689  $ ldafb, iwork, work,
690  $ result( 1 ) )
691  k1 = 1
692  ELSE
693  k1 = 2
694  END IF
695 *
696  IF( info.EQ.0 ) THEN
697  trfcon = .false.
698 *
699 * Compute residual of the computed solution.
700 *
701  CALL dlacpy( 'Full', n, nrhs, bsav, ldb,
702  $ work, ldb )
703  CALL dgbt02( trans, n, n, kl, ku, nrhs,
704  $ asav, lda, x, ldb, work, ldb,
705  $ result( 2 ) )
706 *
707 * Check solution from generated exact
708 * solution.
709 *
710  IF( nofact .OR. ( prefac .AND.
711  $ lsame( equed, 'N' ) ) ) THEN
712  CALL dget04( n, nrhs, x, ldb, xact,
713  $ ldb, rcondc, result( 3 ) )
714  ELSE
715  IF( itran.EQ.1 ) THEN
716  roldc = roldo
717  ELSE
718  roldc = roldi
719  END IF
720  CALL dget04( n, nrhs, x, ldb, xact,
721  $ ldb, roldc, result( 3 ) )
722  END IF
723 *
724 * Check the error bounds from iterative
725 * refinement.
726 *
727  CALL dgbt05( trans, n, kl, ku, nrhs, asav,
728  $ lda, b, ldb, x, ldb, xact,
729  $ ldb, rwork, rwork( nrhs+1 ),
730  $ result( 4 ) )
731  ELSE
732  trfcon = .true.
733  END IF
734 *
735 * Compare RCOND from DGBSVX with the computed
736 * value in RCONDC.
737 *
738  result( 6 ) = dget06( rcond, rcondc )
739 *
740 * Print information about the tests that did
741 * not pass the threshold.
742 *
743  IF( .NOT.trfcon ) THEN
744  DO 80 k = k1, ntests
745  IF( result( k ).GE.thresh ) THEN
746  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
747  $ CALL aladhd( nout, path )
748  IF( prefac ) THEN
749  WRITE( nout, fmt = 9995 )
750  $ 'DGBSVX', fact, trans, n, kl,
751  $ ku, equed, imat, k,
752  $ result( k )
753  ELSE
754  WRITE( nout, fmt = 9996 )
755  $ 'DGBSVX', fact, trans, n, kl,
756  $ ku, imat, k, result( k )
757  END IF
758  nfail = nfail + 1
759  END IF
760  80 CONTINUE
761  nrun = nrun + ntests - k1 + 1
762  ELSE
763  IF( result( 1 ).GE.thresh .AND. .NOT.
764  $ prefac ) THEN
765  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
766  $ CALL aladhd( nout, path )
767  IF( prefac ) THEN
768  WRITE( nout, fmt = 9995 )'DGBSVX',
769  $ fact, trans, n, kl, ku, equed,
770  $ imat, 1, result( 1 )
771  ELSE
772  WRITE( nout, fmt = 9996 )'DGBSVX',
773  $ fact, trans, n, kl, ku, imat, 1,
774  $ result( 1 )
775  END IF
776  nfail = nfail + 1
777  nrun = nrun + 1
778  END IF
779  IF( result( 6 ).GE.thresh ) THEN
780  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
781  $ CALL aladhd( nout, path )
782  IF( prefac ) THEN
783  WRITE( nout, fmt = 9995 )'DGBSVX',
784  $ fact, trans, n, kl, ku, equed,
785  $ imat, 6, result( 6 )
786  ELSE
787  WRITE( nout, fmt = 9996 )'DGBSVX',
788  $ fact, trans, n, kl, ku, imat, 6,
789  $ result( 6 )
790  END IF
791  nfail = nfail + 1
792  nrun = nrun + 1
793  END IF
794  IF( result( 7 ).GE.thresh ) THEN
795  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
796  $ CALL aladhd( nout, path )
797  IF( prefac ) THEN
798  WRITE( nout, fmt = 9995 )'DGBSVX',
799  $ fact, trans, n, kl, ku, equed,
800  $ imat, 7, result( 7 )
801  ELSE
802  WRITE( nout, fmt = 9996 )'DGBSVX',
803  $ fact, trans, n, kl, ku, imat, 7,
804  $ result( 7 )
805  END IF
806  nfail = nfail + 1
807  nrun = nrun + 1
808  END IF
809 *
810  END IF
811  90 CONTINUE
812  100 CONTINUE
813  110 CONTINUE
814  120 CONTINUE
815  130 CONTINUE
816  140 CONTINUE
817  150 CONTINUE
818 *
819 * Print a summary of the results.
820 *
821  CALL alasvm( path, nout, nfail, nrun, nerrs )
822 *
823  9999 FORMAT( ' *** In DDRVGB, LA=', i5, ' is too small for N=', i5,
824  $ ', KU=', i5, ', KL=', i5, / ' ==> Increase LA to at least ',
825  $ i5 )
826  9998 FORMAT( ' *** In DDRVGB, LAFB=', i5, ' is too small for N=', i5,
827  $ ', KU=', i5, ', KL=', i5, /
828  $ ' ==> Increase LAFB to at least ', i5 )
829  9997 FORMAT( 1x, a, ', N=', i5, ', KL=', i5, ', KU=', i5, ', type ',
830  $ i1, ', test(', i1, ')=', g12.5 )
831  9996 FORMAT( 1x, a, '( ''', a1, ''',''', a1, ''',', i5, ',', i5, ',',
832  $ i5, ',...), type ', i1, ', test(', i1, ')=', g12.5 )
833  9995 FORMAT( 1x, a, '( ''', a1, ''',''', a1, ''',', i5, ',', i5, ',',
834  $ i5, ',...), EQUED=''', a1, ''', type ', i1, ', test(', i1,
835  $ ')=', g12.5 )
836 *
837  RETURN
838 *
839 * End of DDRVGB
840 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine dgbt05(TRANS, N, KL, KU, NRHS, AB, LDAB, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DGBT05
Definition: dgbt05.f:178
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dlaqgb(M, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND, AMAX, EQUED)
DLAQGB scales a general band matrix, using row and column scaling factors computed by sgbequ...
Definition: dlaqgb.f:161
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
double precision function dlange(NORM, M, N, A, LDA, WORK)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlange.f:116
subroutine dgbequ(M, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND, AMAX, INFO)
DGBEQU
Definition: dgbequ.f:155
double precision function dlangb(NORM, N, KL, KU, AB, LDAB, WORK)
DLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlangb.f:126
subroutine dgbtrf(M, N, KL, KU, AB, LDAB, IPIV, INFO)
DGBTRF
Definition: dgbtrf.f:146
subroutine dgbt02(TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, RESID)
DGBT02
Definition: dgbt02.f:141
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine derrvx(PATH, NUNIT)
DERRVX
Definition: derrvx.f:57
subroutine dgbsvx(FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
DGBSVX computes the solution to system of linear equations A * X = B for GB matrices ...
Definition: dgbsvx.f:371
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
double precision function dlantb(NORM, UPLO, DIAG, N, K, AB, LDAB, WORK)
DLANTB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a triangular band matrix.
Definition: dlantb.f:142
subroutine dgbsv(N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
DGBSV computes the solution to system of linear equations A * X = B for GB matrices (simple driver) ...
Definition: dgbsv.f:164
subroutine dgbt01(M, N, KL, KU, A, LDA, AFAC, LDAFAC, IPIV, WORK, RESID)
DGBT01
Definition: dgbt01.f:128
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dgbtrs(TRANS, N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
DGBTRS
Definition: dgbtrs.f:140

Here is the call graph for this function:

Here is the caller graph for this function: