LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cdrvpb ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
real  THRESH,
logical  TSTERR,
integer  NMAX,
complex, dimension( * )  A,
complex, dimension( * )  AFAC,
complex, dimension( * )  ASAV,
complex, dimension( * )  B,
complex, dimension( * )  BSAV,
complex, dimension( * )  X,
complex, dimension( * )  XACT,
real, dimension( * )  S,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  NOUT 
)

CDRVPB

Purpose:
 CDRVPB tests the driver routines CPBSV 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 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.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is COMPLEX array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX array, dimension (NMAX*NMAX)
[out]ASAV
          ASAV is COMPLEX array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is COMPLEX array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NRHS)
[out]S
          S is REAL array, dimension (NMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is REAL array, dimension (NMAX+2*NRHS)
[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 161 of file cdrvpb.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: