LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ddrvpb ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
double precision, dimension( * )  A,
double precision, dimension( * )  AFAC,
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 
)

DDRVPB

Purpose:
 DDRVPB tests the driver routines DPBSV 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 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.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]ASAV
          ASAV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[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 (NMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[out]IWORK
          IWORK is INTEGER array, dimension (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 166 of file ddrvpb.f.

166 *
167 * -- LAPACK test routine (version 3.4.0) --
168 * -- LAPACK is a software package provided by Univ. of Tennessee, --
169 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
170 * November 2011
171 *
172 * .. Scalar Arguments ..
173  LOGICAL tsterr
174  INTEGER nmax, nn, nout, nrhs
175  DOUBLE PRECISION thresh
176 * ..
177 * .. Array Arguments ..
178  LOGICAL dotype( * )
179  INTEGER iwork( * ), nval( * )
180  DOUBLE PRECISION a( * ), afac( * ), asav( * ), b( * ),
181  $ bsav( * ), rwork( * ), s( * ), work( * ),
182  $ x( * ), xact( * )
183 * ..
184 *
185 * =====================================================================
186 *
187 * .. Parameters ..
188  DOUBLE PRECISION one, zero
189  parameter ( one = 1.0d+0, zero = 0.0d+0 )
190  INTEGER ntypes, ntests
191  parameter ( ntypes = 8, ntests = 6 )
192  INTEGER nbw
193  parameter ( nbw = 4 )
194 * ..
195 * .. Local Scalars ..
196  LOGICAL equil, nofact, prefac, zerot
197  CHARACTER dist, equed, fact, packit, TYPE, uplo, xtype
198  CHARACTER*3 path
199  INTEGER i, i1, i2, iequed, ifact, ikd, imat, in, info,
200  $ ioff, iuplo, iw, izero, k, k1, kd, kl, koff,
201  $ ku, lda, ldab, mode, n, nb, nbmin, nerrs,
202  $ nfact, nfail, nimat, nkd, nrun, nt
203  DOUBLE PRECISION ainvnm, amax, anorm, cndnum, rcond, rcondc,
204  $ roldc, scond
205 * ..
206 * .. Local Arrays ..
207  CHARACTER equeds( 2 ), facts( 3 )
208  INTEGER iseed( 4 ), iseedy( 4 ), kdval( nbw )
209  DOUBLE PRECISION result( ntests )
210 * ..
211 * .. External Functions ..
212  LOGICAL lsame
213  DOUBLE PRECISION dget06, dlange, dlansb
214  EXTERNAL lsame, dget06, dlange, dlansb
215 * ..
216 * .. External Subroutines ..
217  EXTERNAL aladhd, alaerh, alasvm, dcopy, derrvx, dget04,
221 * ..
222 * .. Intrinsic Functions ..
223  INTRINSIC max, min
224 * ..
225 * .. Scalars in Common ..
226  LOGICAL lerr, ok
227  CHARACTER*32 srnamt
228  INTEGER infot, nunit
229 * ..
230 * .. Common blocks ..
231  COMMON / infoc / infot, nunit, ok, lerr
232  COMMON / srnamc / srnamt
233 * ..
234 * .. Data statements ..
235  DATA iseedy / 1988, 1989, 1990, 1991 /
236  DATA facts / 'F', 'N', 'E' /
237  DATA equeds / 'N', 'Y' /
238 * ..
239 * .. Executable Statements ..
240 *
241 * Initialize constants and the random number seed.
242 *
243  path( 1: 1 ) = 'Double precision'
244  path( 2: 3 ) = 'PB'
245  nrun = 0
246  nfail = 0
247  nerrs = 0
248  DO 10 i = 1, 4
249  iseed( i ) = iseedy( i )
250  10 CONTINUE
251 *
252 * Test the error exits
253 *
254  IF( tsterr )
255  $ CALL derrvx( path, nout )
256  infot = 0
257  kdval( 1 ) = 0
258 *
259 * Set the block size and minimum block size for testing.
260 *
261  nb = 1
262  nbmin = 2
263  CALL xlaenv( 1, nb )
264  CALL xlaenv( 2, nbmin )
265 *
266 * Do for each value of N in NVAL
267 *
268  DO 110 in = 1, nn
269  n = nval( in )
270  lda = max( n, 1 )
271  xtype = 'N'
272 *
273 * Set limits on the number of loop iterations.
274 *
275  nkd = max( 1, min( n, 4 ) )
276  nimat = ntypes
277  IF( n.EQ.0 )
278  $ nimat = 1
279 *
280  kdval( 2 ) = n + ( n+1 ) / 4
281  kdval( 3 ) = ( 3*n-1 ) / 4
282  kdval( 4 ) = ( n+1 ) / 4
283 *
284  DO 100 ikd = 1, nkd
285 *
286 * Do for KD = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This order
287 * makes it easier to skip redundant values for small values
288 * of N.
289 *
290  kd = kdval( ikd )
291  ldab = kd + 1
292 *
293 * Do first for UPLO = 'U', then for UPLO = 'L'
294 *
295  DO 90 iuplo = 1, 2
296  koff = 1
297  IF( iuplo.EQ.1 ) THEN
298  uplo = 'U'
299  packit = 'Q'
300  koff = max( 1, kd+2-n )
301  ELSE
302  uplo = 'L'
303  packit = 'B'
304  END IF
305 *
306  DO 80 imat = 1, nimat
307 *
308 * Do the tests only if DOTYPE( IMAT ) is true.
309 *
310  IF( .NOT.dotype( imat ) )
311  $ GO TO 80
312 *
313 * Skip types 2, 3, or 4 if the matrix size is too small.
314 *
315  zerot = imat.GE.2 .AND. imat.LE.4
316  IF( zerot .AND. n.LT.imat-1 )
317  $ GO TO 80
318 *
319  IF( .NOT.zerot .OR. .NOT.dotype( 1 ) ) THEN
320 *
321 * Set up parameters with DLATB4 and generate a test
322 * matrix with DLATMS.
323 *
324  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm,
325  $ mode, cndnum, dist )
326 *
327  srnamt = 'DLATMS'
328  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
329  $ cndnum, anorm, kd, kd, packit,
330  $ a( koff ), ldab, work, info )
331 *
332 * Check error code from DLATMS.
333 *
334  IF( info.NE.0 ) THEN
335  CALL alaerh( path, 'DLATMS', info, 0, uplo, n,
336  $ n, -1, -1, -1, imat, nfail, nerrs,
337  $ nout )
338  GO TO 80
339  END IF
340  ELSE IF( izero.GT.0 ) THEN
341 *
342 * Use the same matrix for types 3 and 4 as for type
343 * 2 by copying back the zeroed out column,
344 *
345  iw = 2*lda + 1
346  IF( iuplo.EQ.1 ) THEN
347  ioff = ( izero-1 )*ldab + kd + 1
348  CALL dcopy( izero-i1, work( iw ), 1,
349  $ a( ioff-izero+i1 ), 1 )
350  iw = iw + izero - i1
351  CALL dcopy( i2-izero+1, work( iw ), 1,
352  $ a( ioff ), max( ldab-1, 1 ) )
353  ELSE
354  ioff = ( i1-1 )*ldab + 1
355  CALL dcopy( izero-i1, work( iw ), 1,
356  $ a( ioff+izero-i1 ),
357  $ max( ldab-1, 1 ) )
358  ioff = ( izero-1 )*ldab + 1
359  iw = iw + izero - i1
360  CALL dcopy( i2-izero+1, work( iw ), 1,
361  $ a( ioff ), 1 )
362  END IF
363  END IF
364 *
365 * For types 2-4, zero one row and column of the matrix
366 * to test that INFO is returned correctly.
367 *
368  izero = 0
369  IF( zerot ) THEN
370  IF( imat.EQ.2 ) THEN
371  izero = 1
372  ELSE IF( imat.EQ.3 ) THEN
373  izero = n
374  ELSE
375  izero = n / 2 + 1
376  END IF
377 *
378 * Save the zeroed out row and column in WORK(*,3)
379 *
380  iw = 2*lda
381  DO 20 i = 1, min( 2*kd+1, n )
382  work( iw+i ) = zero
383  20 CONTINUE
384  iw = iw + 1
385  i1 = max( izero-kd, 1 )
386  i2 = min( izero+kd, n )
387 *
388  IF( iuplo.EQ.1 ) THEN
389  ioff = ( izero-1 )*ldab + kd + 1
390  CALL dswap( izero-i1, a( ioff-izero+i1 ), 1,
391  $ work( iw ), 1 )
392  iw = iw + izero - i1
393  CALL dswap( i2-izero+1, a( ioff ),
394  $ max( ldab-1, 1 ), work( iw ), 1 )
395  ELSE
396  ioff = ( i1-1 )*ldab + 1
397  CALL dswap( izero-i1, a( ioff+izero-i1 ),
398  $ max( ldab-1, 1 ), work( iw ), 1 )
399  ioff = ( izero-1 )*ldab + 1
400  iw = iw + izero - i1
401  CALL dswap( i2-izero+1, a( ioff ), 1,
402  $ work( iw ), 1 )
403  END IF
404  END IF
405 *
406 * Save a copy of the matrix A in ASAV.
407 *
408  CALL dlacpy( 'Full', kd+1, n, a, ldab, asav, ldab )
409 *
410  DO 70 iequed = 1, 2
411  equed = equeds( iequed )
412  IF( iequed.EQ.1 ) THEN
413  nfact = 3
414  ELSE
415  nfact = 1
416  END IF
417 *
418  DO 60 ifact = 1, nfact
419  fact = facts( ifact )
420  prefac = lsame( fact, 'F' )
421  nofact = lsame( fact, 'N' )
422  equil = lsame( fact, 'E' )
423 *
424  IF( zerot ) THEN
425  IF( prefac )
426  $ GO TO 60
427  rcondc = zero
428 *
429  ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
430 *
431 * Compute the condition number for comparison
432 * with the value returned by DPBSVX (FACT =
433 * 'N' reuses the condition number from the
434 * previous iteration with FACT = 'F').
435 *
436  CALL dlacpy( 'Full', kd+1, n, asav, ldab,
437  $ afac, ldab )
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 dpbequ( uplo, n, kd, afac, ldab, s,
444  $ scond, amax, info )
445  IF( info.EQ.0 .AND. n.GT.0 ) THEN
446  IF( iequed.GT.1 )
447  $ scond = zero
448 *
449 * Equilibrate the matrix.
450 *
451  CALL dlaqsb( uplo, n, kd, afac, ldab,
452  $ s, scond, amax, equed )
453  END IF
454  END IF
455 *
456 * Save the condition number of the
457 * non-equilibrated system for use in DGET04.
458 *
459  IF( equil )
460  $ roldc = rcondc
461 *
462 * Compute the 1-norm of A.
463 *
464  anorm = dlansb( '1', uplo, n, kd, afac, ldab,
465  $ rwork )
466 *
467 * Factor the matrix A.
468 *
469  CALL dpbtrf( uplo, n, kd, afac, ldab, info )
470 *
471 * Form the inverse of A.
472 *
473  CALL dlaset( 'Full', n, n, zero, one, a,
474  $ lda )
475  srnamt = 'DPBTRS'
476  CALL dpbtrs( uplo, n, kd, n, afac, ldab, a,
477  $ lda, info )
478 *
479 * Compute the 1-norm condition number of A.
480 *
481  ainvnm = dlange( '1', n, n, a, lda, rwork )
482  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
483  rcondc = one
484  ELSE
485  rcondc = ( one / anorm ) / ainvnm
486  END IF
487  END IF
488 *
489 * Restore the matrix A.
490 *
491  CALL dlacpy( 'Full', kd+1, n, asav, ldab, a,
492  $ ldab )
493 *
494 * Form an exact solution and set the right hand
495 * side.
496 *
497  srnamt = 'DLARHS'
498  CALL dlarhs( path, xtype, uplo, ' ', n, n, kd,
499  $ kd, nrhs, a, ldab, xact, lda, b,
500  $ lda, iseed, info )
501  xtype = 'C'
502  CALL dlacpy( 'Full', n, nrhs, b, lda, bsav,
503  $ lda )
504 *
505  IF( nofact ) THEN
506 *
507 * --- Test DPBSV ---
508 *
509 * Compute the L*L' or U'*U factorization of the
510 * matrix and solve the system.
511 *
512  CALL dlacpy( 'Full', kd+1, n, a, ldab, afac,
513  $ ldab )
514  CALL dlacpy( 'Full', n, nrhs, b, lda, x,
515  $ lda )
516 *
517  srnamt = 'DPBSV '
518  CALL dpbsv( uplo, n, kd, nrhs, afac, ldab, x,
519  $ lda, info )
520 *
521 * Check error code from DPBSV .
522 *
523  IF( info.NE.izero ) THEN
524  CALL alaerh( path, 'DPBSV ', info, izero,
525  $ uplo, n, n, kd, kd, nrhs,
526  $ imat, nfail, nerrs, nout )
527  GO TO 40
528  ELSE IF( info.NE.0 ) THEN
529  GO TO 40
530  END IF
531 *
532 * Reconstruct matrix from factors and compute
533 * residual.
534 *
535  CALL dpbt01( uplo, n, kd, a, ldab, afac,
536  $ ldab, rwork, result( 1 ) )
537 *
538 * Compute residual of the computed solution.
539 *
540  CALL dlacpy( 'Full', n, nrhs, b, lda, work,
541  $ lda )
542  CALL dpbt02( uplo, n, kd, nrhs, a, ldab, x,
543  $ lda, work, lda, rwork,
544  $ result( 2 ) )
545 *
546 * Check solution from generated exact solution.
547 *
548  CALL dget04( n, nrhs, x, lda, xact, lda,
549  $ rcondc, result( 3 ) )
550  nt = 3
551 *
552 * Print information about the tests that did
553 * not pass the threshold.
554 *
555  DO 30 k = 1, nt
556  IF( result( k ).GE.thresh ) THEN
557  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
558  $ CALL aladhd( nout, path )
559  WRITE( nout, fmt = 9999 )'DPBSV ',
560  $ uplo, n, kd, imat, k, result( k )
561  nfail = nfail + 1
562  END IF
563  30 CONTINUE
564  nrun = nrun + nt
565  40 CONTINUE
566  END IF
567 *
568 * --- Test DPBSVX ---
569 *
570  IF( .NOT.prefac )
571  $ CALL dlaset( 'Full', kd+1, n, zero, zero,
572  $ afac, ldab )
573  CALL dlaset( 'Full', n, nrhs, zero, zero, x,
574  $ lda )
575  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
576 *
577 * Equilibrate the matrix if FACT='F' and
578 * EQUED='Y'
579 *
580  CALL dlaqsb( uplo, n, kd, a, ldab, s, scond,
581  $ amax, equed )
582  END IF
583 *
584 * Solve the system and compute the condition
585 * number and error bounds using DPBSVX.
586 *
587  srnamt = 'DPBSVX'
588  CALL dpbsvx( fact, uplo, n, kd, nrhs, a, ldab,
589  $ afac, ldab, equed, s, b, lda, x,
590  $ lda, rcond, rwork, rwork( nrhs+1 ),
591  $ work, iwork, info )
592 *
593 * Check the error code from DPBSVX.
594 *
595  IF( info.NE.izero ) THEN
596  CALL alaerh( path, 'DPBSVX', info, izero,
597  $ fact // uplo, n, n, kd, kd,
598  $ nrhs, imat, nfail, nerrs, nout )
599  GO TO 60
600  END IF
601 *
602  IF( info.EQ.0 ) THEN
603  IF( .NOT.prefac ) THEN
604 *
605 * Reconstruct matrix from factors and
606 * compute residual.
607 *
608  CALL dpbt01( uplo, n, kd, a, ldab, afac,
609  $ ldab, rwork( 2*nrhs+1 ),
610  $ result( 1 ) )
611  k1 = 1
612  ELSE
613  k1 = 2
614  END IF
615 *
616 * Compute residual of the computed solution.
617 *
618  CALL dlacpy( 'Full', n, nrhs, bsav, lda,
619  $ work, lda )
620  CALL dpbt02( uplo, n, kd, nrhs, asav, ldab,
621  $ x, lda, work, lda,
622  $ rwork( 2*nrhs+1 ), result( 2 ) )
623 *
624 * Check solution from generated exact solution.
625 *
626  IF( nofact .OR. ( prefac .AND. lsame( equed,
627  $ 'N' ) ) ) THEN
628  CALL dget04( n, nrhs, x, lda, xact, lda,
629  $ rcondc, result( 3 ) )
630  ELSE
631  CALL dget04( n, nrhs, x, lda, xact, lda,
632  $ roldc, result( 3 ) )
633  END IF
634 *
635 * Check the error bounds from iterative
636 * refinement.
637 *
638  CALL dpbt05( uplo, n, kd, nrhs, asav, ldab,
639  $ b, lda, x, lda, xact, lda,
640  $ rwork, rwork( nrhs+1 ),
641  $ result( 4 ) )
642  ELSE
643  k1 = 6
644  END IF
645 *
646 * Compare RCOND from DPBSVX with the computed
647 * value in RCONDC.
648 *
649  result( 6 ) = dget06( rcond, rcondc )
650 *
651 * Print information about the tests that did not
652 * pass the threshold.
653 *
654  DO 50 k = k1, 6
655  IF( result( k ).GE.thresh ) THEN
656  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
657  $ CALL aladhd( nout, path )
658  IF( prefac ) THEN
659  WRITE( nout, fmt = 9997 )'DPBSVX',
660  $ fact, uplo, n, kd, equed, imat, k,
661  $ result( k )
662  ELSE
663  WRITE( nout, fmt = 9998 )'DPBSVX',
664  $ fact, uplo, n, kd, imat, k,
665  $ result( k )
666  END IF
667  nfail = nfail + 1
668  END IF
669  50 CONTINUE
670  nrun = nrun + 7 - k1
671  60 CONTINUE
672  70 CONTINUE
673  80 CONTINUE
674  90 CONTINUE
675  100 CONTINUE
676  110 CONTINUE
677 *
678 * Print a summary of the results.
679 *
680  CALL alasvm( path, nout, nfail, nrun, nerrs )
681 *
682  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', KD =', i5,
683  $ ', type ', i1, ', test(', i1, ')=', g12.5 )
684  9998 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ', i5, ', ', i5,
685  $ ', ... ), type ', i1, ', test(', i1, ')=', g12.5 )
686  9997 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ', i5, ', ', i5,
687  $ ', ... ), EQUED=''', a1, ''', type ', i1, ', test(', i1,
688  $ ')=', g12.5 )
689  RETURN
690 *
691 * End of DDRVPB
692 *
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
double precision function dlansb(NORM, UPLO, N, K, AB, LDAB, WORK)
DLANSB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric band matrix.
Definition: dlansb.f:131
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 dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dlaqsb(UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, EQUED)
DLAQSB scales a symmetric/Hermitian band matrix, using scaling factors computed by spbequ...
Definition: dlaqsb.f:142
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 dpbequ(UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, INFO)
DPBEQU
Definition: dpbequ.f:131
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine dpbt01(UPLO, N, KD, A, LDA, AFAC, LDAFAC, RWORK, RESID)
DPBT01
Definition: dpbt01.f:121
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 dpbtrf(UPLO, N, KD, AB, LDAB, INFO)
DPBTRF
Definition: dpbtrf.f:144
subroutine dpbsv(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
DPBSV computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: dpbsv.f:166
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 dpbsvx(FACT, UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
DPBSVX computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: dpbsvx.f:345
subroutine derrvx(PATH, NUNIT)
DERRVX
Definition: derrvx.f:57
subroutine dpbtrs(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
DPBTRS
Definition: dpbtrs.f:123
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dpbt05(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DPBT05
Definition: dpbt05.f:173
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 dpbt02(UPLO, N, KD, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DPBT02
Definition: dpbt02.f:138

Here is the call graph for this function:

Here is the caller graph for this function: