LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cchkgb ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
real  THRESH,
logical  TSTERR,
complex, dimension( * )  A,
integer  LA,
complex, dimension( * )  AFAC,
integer  LAFAC,
complex, dimension( * )  B,
complex, dimension( * )  X,
complex, dimension( * )  XACT,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

CCHKGB

Purpose:
 CCHKGB tests CGBTRF, -TRS, -RFS, and -CON
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]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[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]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the blocksize NB.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[in]THRESH
          THRESH is REAL
          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 COMPLEX array, dimension (LA)
[in]LA
          LA is INTEGER
          The length of the array A.  LA >= (KLMAX+KUMAX+1)*NMAX
          where KLMAX is the largest entry in the local array KLVAL,
                KUMAX is the largest entry in the local array KUVAL and
                NMAX is the largest entry in the input array NVAL.
[out]AFAC
          AFAC is COMPLEX array, dimension (LAFAC)
[in]LAFAC
          LAFAC is INTEGER
          The length of the array AFAC. LAFAC >= (2*KLMAX+KUMAX+1)*NMAX
          where KLMAX is the largest entry in the local array KLVAL,
                KUMAX is the largest entry in the local array KUVAL and
                NMAX is the largest entry in the input array NVAL.
[out]B
          B is COMPLEX array, dimension (NMAX*NSMAX)
[out]X
          X is COMPLEX array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NSMAX,NMAX))
[out]RWORK
          RWORK is REAL array, dimension
                      (max(NMAX,2*NSMAX))
[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 193 of file cchkgb.f.

193 *
194 * -- LAPACK test routine (version 3.4.0) --
195 * -- LAPACK is a software package provided by Univ. of Tennessee, --
196 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
197 * November 2011
198 *
199 * .. Scalar Arguments ..
200  LOGICAL tsterr
201  INTEGER la, lafac, nm, nn, nnb, nns, nout
202  REAL thresh
203 * ..
204 * .. Array Arguments ..
205  LOGICAL dotype( * )
206  INTEGER iwork( * ), mval( * ), nbval( * ), nsval( * ),
207  $ nval( * )
208  REAL rwork( * )
209  COMPLEX a( * ), afac( * ), b( * ), work( * ), x( * ),
210  $ xact( * )
211 * ..
212 *
213 * =====================================================================
214 *
215 * .. Parameters ..
216  REAL one, zero
217  parameter ( one = 1.0e+0, zero = 0.0e+0 )
218  INTEGER ntypes, ntests
219  parameter ( ntypes = 8, ntests = 7 )
220  INTEGER nbw, ntran
221  parameter ( nbw = 4, ntran = 3 )
222 * ..
223 * .. Local Scalars ..
224  LOGICAL trfcon, zerot
225  CHARACTER dist, norm, trans, TYPE, xtype
226  CHARACTER*3 path
227  INTEGER i, i1, i2, ikl, iku, im, imat, in, inb, info,
228  $ ioff, irhs, itran, izero, j, k, kl, koff, ku,
229  $ lda, ldafac, ldb, m, mode, n, nb, nerrs, nfail,
230  $ nimat, nkl, nku, nrhs, nrun
231  REAL ainvnm, anorm, anormi, anormo, cndnum, rcond,
232  $ rcondc, rcondi, rcondo
233 * ..
234 * .. Local Arrays ..
235  CHARACTER transs( ntran )
236  INTEGER iseed( 4 ), iseedy( 4 ), klval( nbw ),
237  $ kuval( nbw )
238  REAL result( ntests )
239 * ..
240 * .. External Functions ..
241  REAL clangb, clange, sget06
242  EXTERNAL clangb, clange, sget06
243 * ..
244 * .. External Subroutines ..
245  EXTERNAL alaerh, alahd, alasum, ccopy, cerrge, cgbcon,
248  $ xlaenv
249 * ..
250 * .. Intrinsic Functions ..
251  INTRINSIC cmplx, max, min
252 * ..
253 * .. Scalars in Common ..
254  LOGICAL lerr, ok
255  CHARACTER*32 srnamt
256  INTEGER infot, nunit
257 * ..
258 * .. Common blocks ..
259  COMMON / infoc / infot, nunit, ok, lerr
260  COMMON / srnamc / srnamt
261 * ..
262 * .. Data statements ..
263  DATA iseedy / 1988, 1989, 1990, 1991 / ,
264  $ transs / 'N', 'T', 'C' /
265 * ..
266 * .. Executable Statements ..
267 *
268 * Initialize constants and the random number seed.
269 *
270  path( 1: 1 ) = 'Complex precision'
271  path( 2: 3 ) = 'GB'
272  nrun = 0
273  nfail = 0
274  nerrs = 0
275  DO 10 i = 1, 4
276  iseed( i ) = iseedy( i )
277  10 CONTINUE
278 *
279 * Test the error exits
280 *
281  IF( tsterr )
282  $ CALL cerrge( path, nout )
283  infot = 0
284 *
285 * Initialize the first value for the lower and upper bandwidths.
286 *
287  klval( 1 ) = 0
288  kuval( 1 ) = 0
289 *
290 * Do for each value of M in MVAL
291 *
292  DO 160 im = 1, nm
293  m = mval( im )
294 *
295 * Set values to use for the lower bandwidth.
296 *
297  klval( 2 ) = m + ( m+1 ) / 4
298 *
299 * KLVAL( 2 ) = MAX( M-1, 0 )
300 *
301  klval( 3 ) = ( 3*m-1 ) / 4
302  klval( 4 ) = ( m+1 ) / 4
303 *
304 * Do for each value of N in NVAL
305 *
306  DO 150 in = 1, nn
307  n = nval( in )
308  xtype = 'N'
309 *
310 * Set values to use for the upper bandwidth.
311 *
312  kuval( 2 ) = n + ( n+1 ) / 4
313 *
314 * KUVAL( 2 ) = MAX( N-1, 0 )
315 *
316  kuval( 3 ) = ( 3*n-1 ) / 4
317  kuval( 4 ) = ( n+1 ) / 4
318 *
319 * Set limits on the number of loop iterations.
320 *
321  nkl = min( m+1, 4 )
322  IF( n.EQ.0 )
323  $ nkl = 2
324  nku = min( n+1, 4 )
325  IF( m.EQ.0 )
326  $ nku = 2
327  nimat = ntypes
328  IF( m.LE.0 .OR. n.LE.0 )
329  $ nimat = 1
330 *
331  DO 140 ikl = 1, nkl
332 *
333 * Do for KL = 0, (5*M+1)/4, (3M-1)/4, and (M+1)/4. This
334 * order makes it easier to skip redundant values for small
335 * values of M.
336 *
337  kl = klval( ikl )
338  DO 130 iku = 1, nku
339 *
340 * Do for KU = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This
341 * order makes it easier to skip redundant values for
342 * small values of N.
343 *
344  ku = kuval( iku )
345 *
346 * Check that A and AFAC are big enough to generate this
347 * matrix.
348 *
349  lda = kl + ku + 1
350  ldafac = 2*kl + ku + 1
351  IF( ( lda*n ).GT.la .OR. ( ldafac*n ).GT.lafac ) THEN
352  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
353  $ CALL alahd( nout, path )
354  IF( n*( kl+ku+1 ).GT.la ) THEN
355  WRITE( nout, fmt = 9999 )la, m, n, kl, ku,
356  $ n*( kl+ku+1 )
357  nerrs = nerrs + 1
358  END IF
359  IF( n*( 2*kl+ku+1 ).GT.lafac ) THEN
360  WRITE( nout, fmt = 9998 )lafac, m, n, kl, ku,
361  $ n*( 2*kl+ku+1 )
362  nerrs = nerrs + 1
363  END IF
364  GO TO 130
365  END IF
366 *
367  DO 120 imat = 1, nimat
368 *
369 * Do the tests only if DOTYPE( IMAT ) is true.
370 *
371  IF( .NOT.dotype( imat ) )
372  $ GO TO 120
373 *
374 * Skip types 2, 3, or 4 if the matrix size is too
375 * small.
376 *
377  zerot = imat.GE.2 .AND. imat.LE.4
378  IF( zerot .AND. n.LT.imat-1 )
379  $ GO TO 120
380 *
381  IF( .NOT.zerot .OR. .NOT.dotype( 1 ) ) THEN
382 *
383 * Set up parameters with CLATB4 and generate a
384 * test matrix with CLATMS.
385 *
386  CALL clatb4( path, imat, m, n, TYPE, kl, ku,
387  $ anorm, mode, cndnum, dist )
388 *
389  koff = max( 1, ku+2-n )
390  DO 20 i = 1, koff - 1
391  a( i ) = zero
392  20 CONTINUE
393  srnamt = 'CLATMS'
394  CALL clatms( m, n, dist, iseed, TYPE, rwork,
395  $ mode, cndnum, anorm, kl, ku, 'Z',
396  $ a( koff ), lda, work, info )
397 *
398 * Check the error code from CLATMS.
399 *
400  IF( info.NE.0 ) THEN
401  CALL alaerh( path, 'CLATMS', info, 0, ' ', m,
402  $ n, kl, ku, -1, imat, nfail,
403  $ nerrs, nout )
404  GO TO 120
405  END IF
406  ELSE IF( izero.GT.0 ) THEN
407 *
408 * Use the same matrix for types 3 and 4 as for
409 * type 2 by copying back the zeroed out column.
410 *
411  CALL ccopy( i2-i1+1, b, 1, a( ioff+i1 ), 1 )
412  END IF
413 *
414 * For types 2, 3, and 4, zero one or more columns of
415 * the matrix to test that INFO is returned correctly.
416 *
417  izero = 0
418  IF( zerot ) THEN
419  IF( imat.EQ.2 ) THEN
420  izero = 1
421  ELSE IF( imat.EQ.3 ) THEN
422  izero = min( m, n )
423  ELSE
424  izero = min( m, n ) / 2 + 1
425  END IF
426  ioff = ( izero-1 )*lda
427  IF( imat.LT.4 ) THEN
428 *
429 * Store the column to be zeroed out in B.
430 *
431  i1 = max( 1, ku+2-izero )
432  i2 = min( kl+ku+1, ku+1+( m-izero ) )
433  CALL ccopy( i2-i1+1, a( ioff+i1 ), 1, b, 1 )
434 *
435  DO 30 i = i1, i2
436  a( ioff+i ) = zero
437  30 CONTINUE
438  ELSE
439  DO 50 j = izero, n
440  DO 40 i = max( 1, ku+2-j ),
441  $ min( kl+ku+1, ku+1+( m-j ) )
442  a( ioff+i ) = zero
443  40 CONTINUE
444  ioff = ioff + lda
445  50 CONTINUE
446  END IF
447  END IF
448 *
449 * These lines, if used in place of the calls in the
450 * loop over INB, cause the code to bomb on a Sun
451 * SPARCstation.
452 *
453 * ANORMO = CLANGB( 'O', N, KL, KU, A, LDA, RWORK )
454 * ANORMI = CLANGB( 'I', N, KL, KU, A, LDA, RWORK )
455 *
456 * Do for each blocksize in NBVAL
457 *
458  DO 110 inb = 1, nnb
459  nb = nbval( inb )
460  CALL xlaenv( 1, nb )
461 *
462 * Compute the LU factorization of the band matrix.
463 *
464  IF( m.GT.0 .AND. n.GT.0 )
465  $ CALL clacpy( 'Full', kl+ku+1, n, a, lda,
466  $ afac( kl+1 ), ldafac )
467  srnamt = 'CGBTRF'
468  CALL cgbtrf( m, n, kl, ku, afac, ldafac, iwork,
469  $ info )
470 *
471 * Check error code from CGBTRF.
472 *
473  IF( info.NE.izero )
474  $ CALL alaerh( path, 'CGBTRF', info, izero,
475  $ ' ', m, n, kl, ku, nb, imat,
476  $ nfail, nerrs, nout )
477  trfcon = .false.
478 *
479 *+ TEST 1
480 * Reconstruct matrix from factors and compute
481 * residual.
482 *
483  CALL cgbt01( m, n, kl, ku, a, lda, afac, ldafac,
484  $ iwork, work, result( 1 ) )
485 *
486 * Print information about the tests so far that
487 * did not pass the threshold.
488 *
489  IF( result( 1 ).GE.thresh ) THEN
490  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
491  $ CALL alahd( nout, path )
492  WRITE( nout, fmt = 9997 )m, n, kl, ku, nb,
493  $ imat, 1, result( 1 )
494  nfail = nfail + 1
495  END IF
496  nrun = nrun + 1
497 *
498 * Skip the remaining tests if this is not the
499 * first block size or if M .ne. N.
500 *
501  IF( inb.GT.1 .OR. m.NE.n )
502  $ GO TO 110
503 *
504  anormo = clangb( 'O', n, kl, ku, a, lda, rwork )
505  anormi = clangb( 'I', n, kl, ku, a, lda, rwork )
506 *
507  IF( info.EQ.0 ) THEN
508 *
509 * Form the inverse of A so we can get a good
510 * estimate of CNDNUM = norm(A) * norm(inv(A)).
511 *
512  ldb = max( 1, n )
513  CALL claset( 'Full', n, n, cmplx( zero ),
514  $ cmplx( one ), work, ldb )
515  srnamt = 'CGBTRS'
516  CALL cgbtrs( 'No transpose', n, kl, ku, n,
517  $ afac, ldafac, iwork, work, ldb,
518  $ info )
519 *
520 * Compute the 1-norm condition number of A.
521 *
522  ainvnm = clange( 'O', n, n, work, ldb,
523  $ rwork )
524  IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
525  rcondo = one
526  ELSE
527  rcondo = ( one / anormo ) / ainvnm
528  END IF
529 *
530 * Compute the infinity-norm condition number of
531 * A.
532 *
533  ainvnm = clange( 'I', n, n, work, ldb,
534  $ rwork )
535  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
536  rcondi = one
537  ELSE
538  rcondi = ( one / anormi ) / ainvnm
539  END IF
540  ELSE
541 *
542 * Do only the condition estimate if INFO.NE.0.
543 *
544  trfcon = .true.
545  rcondo = zero
546  rcondi = zero
547  END IF
548 *
549 * Skip the solve tests if the matrix is singular.
550 *
551  IF( trfcon )
552  $ GO TO 90
553 *
554  DO 80 irhs = 1, nns
555  nrhs = nsval( irhs )
556  xtype = 'N'
557 *
558  DO 70 itran = 1, ntran
559  trans = transs( itran )
560  IF( itran.EQ.1 ) THEN
561  rcondc = rcondo
562  norm = 'O'
563  ELSE
564  rcondc = rcondi
565  norm = 'I'
566  END IF
567 *
568 *+ TEST 2:
569 * Solve and compute residual for A * X = B.
570 *
571  srnamt = 'CLARHS'
572  CALL clarhs( path, xtype, ' ', trans, n,
573  $ n, kl, ku, nrhs, a, lda,
574  $ xact, ldb, b, ldb, iseed,
575  $ info )
576  xtype = 'C'
577  CALL clacpy( 'Full', n, nrhs, b, ldb, x,
578  $ ldb )
579 *
580  srnamt = 'CGBTRS'
581  CALL cgbtrs( trans, n, kl, ku, nrhs, afac,
582  $ ldafac, iwork, x, ldb, info )
583 *
584 * Check error code from CGBTRS.
585 *
586  IF( info.NE.0 )
587  $ CALL alaerh( path, 'CGBTRS', info, 0,
588  $ trans, n, n, kl, ku, -1,
589  $ imat, nfail, nerrs, nout )
590 *
591  CALL clacpy( 'Full', n, nrhs, b, ldb,
592  $ work, ldb )
593  CALL cgbt02( trans, m, n, kl, ku, nrhs, a,
594  $ lda, x, ldb, work, ldb,
595  $ result( 2 ) )
596 *
597 *+ TEST 3:
598 * Check solution from generated exact
599 * solution.
600 *
601  CALL cget04( n, nrhs, x, ldb, xact, ldb,
602  $ rcondc, result( 3 ) )
603 *
604 *+ TESTS 4, 5, 6:
605 * Use iterative refinement to improve the
606 * solution.
607 *
608  srnamt = 'CGBRFS'
609  CALL cgbrfs( trans, n, kl, ku, nrhs, a,
610  $ lda, afac, ldafac, iwork, b,
611  $ ldb, x, ldb, rwork,
612  $ rwork( nrhs+1 ), work,
613  $ rwork( 2*nrhs+1 ), info )
614 *
615 * Check error code from CGBRFS.
616 *
617  IF( info.NE.0 )
618  $ CALL alaerh( path, 'CGBRFS', info, 0,
619  $ trans, n, n, kl, ku, nrhs,
620  $ imat, nfail, nerrs, nout )
621 *
622  CALL cget04( n, nrhs, x, ldb, xact, ldb,
623  $ rcondc, result( 4 ) )
624  CALL cgbt05( trans, n, kl, ku, nrhs, a,
625  $ lda, b, ldb, x, ldb, xact,
626  $ ldb, rwork, rwork( nrhs+1 ),
627  $ result( 5 ) )
628 *
629 * Print information about the tests that did
630 * not pass the threshold.
631 *
632  DO 60 k = 2, 6
633  IF( result( k ).GE.thresh ) THEN
634  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
635  $ CALL alahd( nout, path )
636  WRITE( nout, fmt = 9996 )trans, n,
637  $ kl, ku, nrhs, imat, k,
638  $ result( k )
639  nfail = nfail + 1
640  END IF
641  60 CONTINUE
642  nrun = nrun + 5
643  70 CONTINUE
644  80 CONTINUE
645 *
646 *+ TEST 7:
647 * Get an estimate of RCOND = 1/CNDNUM.
648 *
649  90 CONTINUE
650  DO 100 itran = 1, 2
651  IF( itran.EQ.1 ) THEN
652  anorm = anormo
653  rcondc = rcondo
654  norm = 'O'
655  ELSE
656  anorm = anormi
657  rcondc = rcondi
658  norm = 'I'
659  END IF
660  srnamt = 'CGBCON'
661  CALL cgbcon( norm, n, kl, ku, afac, ldafac,
662  $ iwork, anorm, rcond, work,
663  $ rwork, info )
664 *
665 * Check error code from CGBCON.
666 *
667  IF( info.NE.0 )
668  $ CALL alaerh( path, 'CGBCON', info, 0,
669  $ norm, n, n, kl, ku, -1, imat,
670  $ nfail, nerrs, nout )
671 *
672  result( 7 ) = sget06( rcond, rcondc )
673 *
674 * Print information about the tests that did
675 * not pass the threshold.
676 *
677  IF( result( 7 ).GE.thresh ) THEN
678  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
679  $ CALL alahd( nout, path )
680  WRITE( nout, fmt = 9995 )norm, n, kl, ku,
681  $ imat, 7, result( 7 )
682  nfail = nfail + 1
683  END IF
684  nrun = nrun + 1
685  100 CONTINUE
686  110 CONTINUE
687  120 CONTINUE
688  130 CONTINUE
689  140 CONTINUE
690  150 CONTINUE
691  160 CONTINUE
692 *
693 * Print a summary of the results.
694 *
695  CALL alasum( path, nout, nfail, nrun, nerrs )
696 *
697  9999 FORMAT( ' *** In CCHKGB, LA=', i5, ' is too small for M=', i5,
698  $ ', N=', i5, ', KL=', i4, ', KU=', i4,
699  $ / ' ==> Increase LA to at least ', i5 )
700  9998 FORMAT( ' *** In CCHKGB, LAFAC=', i5, ' is too small for M=', i5,
701  $ ', N=', i5, ', KL=', i4, ', KU=', i4,
702  $ / ' ==> Increase LAFAC to at least ', i5 )
703  9997 FORMAT( ' M =', i5, ', N =', i5, ', KL=', i5, ', KU=', i5,
704  $ ', NB =', i4, ', type ', i1, ', test(', i1, ')=', g12.5 )
705  9996 FORMAT( ' TRANS=''', a1, ''', N=', i5, ', KL=', i5, ', KU=', i5,
706  $ ', NRHS=', i3, ', type ', i1, ', test(', i1, ')=', g12.5 )
707  9995 FORMAT( ' NORM =''', a1, ''', N=', i5, ', KL=', i5, ', KU=', i5,
708  $ ',', 10x, ' type ', i1, ', test(', i1, ')=', g12.5 )
709 *
710  RETURN
711 *
712 * End of CCHKGB
713 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine clarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
CLARHS
Definition: clarhs.f:211
subroutine cgbtrf(M, N, KL, KU, AB, LDAB, IPIV, INFO)
CGBTRF
Definition: cgbtrf.f:146
subroutine cgbt02(TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, RESID)
CGBT02
Definition: cgbt02.f:141
subroutine cgbt01(M, N, KL, KU, A, LDA, AFAC, LDAFAC, IPIV, WORK, RESID)
CGBT01
Definition: cgbt01.f:128
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
real function clange(NORM, M, N, A, LDA, WORK)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: clange.f:117
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine claset(UPLO, M, N, ALPHA, BETA, A, LDA)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: claset.f:108
subroutine clatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
CLATMS
Definition: clatms.f:334
subroutine cgbt05(TRANS, N, KL, KU, NRHS, AB, LDAB, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
CGBT05
Definition: cgbt05.f:178
subroutine cgbcon(NORM, N, KL, KU, AB, LDAB, IPIV, ANORM, RCOND, WORK, RWORK, INFO)
CGBCON
Definition: cgbcon.f:149
real function clangb(NORM, N, KL, KU, AB, LDAB, WORK)
CLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: clangb.f:127
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52
subroutine cgbrfs(TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
CGBRFS
Definition: cgbrfs.f:208
subroutine clatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
CLATB4
Definition: clatb4.f:123
subroutine cerrge(PATH, NUNIT)
CERRGE
Definition: cerrge.f:57
subroutine cget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
CGET04
Definition: cget04.f:104
subroutine cgbtrs(TRANS, N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
CGBTRS
Definition: cgbtrs.f:140
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75

Here is the call graph for this function:

Here is the caller graph for this function: