LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
sdrvpb.f
Go to the documentation of this file.
1 *> \brief \b SDRVPB
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE SDRVPB( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12 * A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
13 * RWORK, IWORK, NOUT )
14 *
15 * .. Scalar Arguments ..
16 * LOGICAL TSTERR
17 * INTEGER NMAX, NN, NOUT, NRHS
18 * REAL THRESH
19 * ..
20 * .. Array Arguments ..
21 * LOGICAL DOTYPE( * )
22 * INTEGER IWORK( * ), NVAL( * )
23 * REAL A( * ), AFAC( * ), ASAV( * ), B( * ),
24 * $ BSAV( * ), RWORK( * ), S( * ), WORK( * ),
25 * $ X( * ), XACT( * )
26 * ..
27 *
28 *
29 *> \par Purpose:
30 * =============
31 *>
32 *> \verbatim
33 *>
34 *> SDRVPB tests the driver routines SPBSV and -SVX.
35 *> \endverbatim
36 *
37 * Arguments:
38 * ==========
39 *
40 *> \param[in] DOTYPE
41 *> \verbatim
42 *> DOTYPE is LOGICAL array, dimension (NTYPES)
43 *> The matrix types to be used for testing. Matrices of type j
44 *> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
45 *> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
46 *> \endverbatim
47 *>
48 *> \param[in] NN
49 *> \verbatim
50 *> NN is INTEGER
51 *> The number of values of N contained in the vector NVAL.
52 *> \endverbatim
53 *>
54 *> \param[in] NVAL
55 *> \verbatim
56 *> NVAL is INTEGER array, dimension (NN)
57 *> The values of the matrix dimension N.
58 *> \endverbatim
59 *>
60 *> \param[in] NRHS
61 *> \verbatim
62 *> NRHS is INTEGER
63 *> The number of right hand side vectors to be generated for
64 *> each linear system.
65 *> \endverbatim
66 *>
67 *> \param[in] THRESH
68 *> \verbatim
69 *> THRESH is REAL
70 *> The threshold value for the test ratios. A result is
71 *> included in the output file if RESULT >= THRESH. To have
72 *> every test ratio printed, use THRESH = 0.
73 *> \endverbatim
74 *>
75 *> \param[in] TSTERR
76 *> \verbatim
77 *> TSTERR is LOGICAL
78 *> Flag that indicates whether error exits are to be tested.
79 *> \endverbatim
80 *>
81 *> \param[in] NMAX
82 *> \verbatim
83 *> NMAX is INTEGER
84 *> The maximum value permitted for N, used in dimensioning the
85 *> work arrays.
86 *> \endverbatim
87 *>
88 *> \param[out] A
89 *> \verbatim
90 *> A is REAL array, dimension (NMAX*NMAX)
91 *> \endverbatim
92 *>
93 *> \param[out] AFAC
94 *> \verbatim
95 *> AFAC is REAL array, dimension (NMAX*NMAX)
96 *> \endverbatim
97 *>
98 *> \param[out] ASAV
99 *> \verbatim
100 *> ASAV is REAL array, dimension (NMAX*NMAX)
101 *> \endverbatim
102 *>
103 *> \param[out] B
104 *> \verbatim
105 *> B is REAL array, dimension (NMAX*NRHS)
106 *> \endverbatim
107 *>
108 *> \param[out] BSAV
109 *> \verbatim
110 *> BSAV is REAL array, dimension (NMAX*NRHS)
111 *> \endverbatim
112 *>
113 *> \param[out] X
114 *> \verbatim
115 *> X is REAL array, dimension (NMAX*NRHS)
116 *> \endverbatim
117 *>
118 *> \param[out] XACT
119 *> \verbatim
120 *> XACT is REAL array, dimension (NMAX*NRHS)
121 *> \endverbatim
122 *>
123 *> \param[out] S
124 *> \verbatim
125 *> S is REAL array, dimension (NMAX)
126 *> \endverbatim
127 *>
128 *> \param[out] WORK
129 *> \verbatim
130 *> WORK is REAL array, dimension
131 *> (NMAX*max(3,NRHS))
132 *> \endverbatim
133 *>
134 *> \param[out] RWORK
135 *> \verbatim
136 *> RWORK is REAL array, dimension (NMAX+2*NRHS)
137 *> \endverbatim
138 *>
139 *> \param[out] IWORK
140 *> \verbatim
141 *> IWORK is INTEGER array, dimension (NMAX)
142 *> \endverbatim
143 *>
144 *> \param[in] NOUT
145 *> \verbatim
146 *> NOUT is INTEGER
147 *> The unit number for output.
148 *> \endverbatim
149 *
150 * Authors:
151 * ========
152 *
153 *> \author Univ. of Tennessee
154 *> \author Univ. of California Berkeley
155 *> \author Univ. of Colorado Denver
156 *> \author NAG Ltd.
157 *
158 *> \date November 2011
159 *
160 *> \ingroup single_lin
161 *
162 * =====================================================================
163  SUBROUTINE sdrvpb( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
164  $ a, afac, asav, b, bsav, x, xact, s, work,
165  $ rwork, iwork, nout )
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  REAL THRESH
176 * ..
177 * .. Array Arguments ..
178  LOGICAL DOTYPE( * )
179  INTEGER IWORK( * ), NVAL( * )
180  REAL A( * ), AFAC( * ), ASAV( * ), B( * ),
181  $ bsav( * ), rwork( * ), s( * ), work( * ),
182  $ x( * ), xact( * )
183 * ..
184 *
185 * =====================================================================
186 *
187 * .. Parameters ..
188  REAL ONE, ZERO
189  parameter ( one = 1.0e+0, zero = 0.0e+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  REAL 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  REAL RESULT( ntests )
210 * ..
211 * .. External Functions ..
212  LOGICAL LSAME
213  REAL SGET06, SLANGE, SLANSB
214  EXTERNAL lsame, sget06, slange, slansb
215 * ..
216 * .. External Subroutines ..
217  EXTERNAL aladhd, alaerh, alasvm, scopy, serrvx, sget04,
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 ) = 'Single 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 serrvx( 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 SLATB4 and generate a test
322 * matrix with SLATMS.
323 *
324  CALL slatb4( path, imat, n, n, TYPE, KL, KU, ANORM,
325  $ mode, cndnum, dist )
326 *
327  srnamt = 'SLATMS'
328  CALL slatms( 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 SLATMS.
333 *
334  IF( info.NE.0 ) THEN
335  CALL alaerh( path, 'SLATMS', 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 scopy( izero-i1, work( iw ), 1,
349  $ a( ioff-izero+i1 ), 1 )
350  iw = iw + izero - i1
351  CALL scopy( i2-izero+1, work( iw ), 1,
352  $ a( ioff ), max( ldab-1, 1 ) )
353  ELSE
354  ioff = ( i1-1 )*ldab + 1
355  CALL scopy( 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 scopy( 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 sswap( izero-i1, a( ioff-izero+i1 ), 1,
391  $ work( iw ), 1 )
392  iw = iw + izero - i1
393  CALL sswap( i2-izero+1, a( ioff ),
394  $ max( ldab-1, 1 ), work( iw ), 1 )
395  ELSE
396  ioff = ( i1-1 )*ldab + 1
397  CALL sswap( 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 sswap( 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 slacpy( '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 SPBSVX (FACT =
433 * 'N' reuses the condition number from the
434 * previous iteration with FACT = 'F').
435 *
436  CALL slacpy( '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 spbequ( 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 slaqsb( 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 SGET04.
458 *
459  IF( equil )
460  $ roldc = rcondc
461 *
462 * Compute the 1-norm of A.
463 *
464  anorm = slansb( '1', uplo, n, kd, afac, ldab,
465  $ rwork )
466 *
467 * Factor the matrix A.
468 *
469  CALL spbtrf( uplo, n, kd, afac, ldab, info )
470 *
471 * Form the inverse of A.
472 *
473  CALL slaset( 'Full', n, n, zero, one, a,
474  $ lda )
475  srnamt = 'SPBTRS'
476  CALL spbtrs( uplo, n, kd, n, afac, ldab, a,
477  $ lda, info )
478 *
479 * Compute the 1-norm condition number of A.
480 *
481  ainvnm = slange( '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 slacpy( '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 = 'SLARHS'
498  CALL slarhs( path, xtype, uplo, ' ', n, n, kd,
499  $ kd, nrhs, a, ldab, xact, lda, b,
500  $ lda, iseed, info )
501  xtype = 'C'
502  CALL slacpy( 'Full', n, nrhs, b, lda, bsav,
503  $ lda )
504 *
505  IF( nofact ) THEN
506 *
507 * --- Test SPBSV ---
508 *
509 * Compute the L*L' or U'*U factorization of the
510 * matrix and solve the system.
511 *
512  CALL slacpy( 'Full', kd+1, n, a, ldab, afac,
513  $ ldab )
514  CALL slacpy( 'Full', n, nrhs, b, lda, x,
515  $ lda )
516 *
517  srnamt = 'SPBSV '
518  CALL spbsv( uplo, n, kd, nrhs, afac, ldab, x,
519  $ lda, info )
520 *
521 * Check error code from SPBSV .
522 *
523  IF( info.NE.izero ) THEN
524  CALL alaerh( path, 'SPBSV ', 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 spbt01( uplo, n, kd, a, ldab, afac,
536  $ ldab, rwork, result( 1 ) )
537 *
538 * Compute residual of the computed solution.
539 *
540  CALL slacpy( 'Full', n, nrhs, b, lda, work,
541  $ lda )
542  CALL spbt02( 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 sget04( 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 )'SPBSV ',
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 SPBSVX ---
569 *
570  IF( .NOT.prefac )
571  $ CALL slaset( 'Full', kd+1, n, zero, zero,
572  $ afac, ldab )
573  CALL slaset( '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 slaqsb( 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 SPBSVX.
586 *
587  srnamt = 'SPBSVX'
588  CALL spbsvx( 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 SPBSVX.
594 *
595  IF( info.NE.izero ) THEN
596  CALL alaerh( path, 'SPBSVX', 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 spbt01( 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 slacpy( 'Full', n, nrhs, bsav, lda,
619  $ work, lda )
620  CALL spbt02( 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 sget04( n, nrhs, x, lda, xact, lda,
629  $ rcondc, result( 3 ) )
630  ELSE
631  CALL sget04( 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 spbt05( 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 SPBSVX with the computed
647 * value in RCONDC.
648 *
649  result( 6 ) = sget06( 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 )'SPBSVX',
660  $ fact, uplo, n, kd, equed, imat, k,
661  $ result( k )
662  ELSE
663  WRITE( nout, fmt = 9998 )'SPBSVX',
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 SDRVPB
692 *
693  END
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 spbequ(UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, INFO)
SPBEQU
Definition: spbequ.f:131
subroutine slatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
SLATB4
Definition: slatb4.f:122
subroutine spbt02(UPLO, N, KD, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
SPBT02
Definition: spbt02.f:138
subroutine slarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
SLARHS
Definition: slarhs.f:206
subroutine spbsvx(FACT, UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
SPBSVX computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: spbsvx.f:345
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine spbtrs(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
SPBTRS
Definition: spbtrs.f:123
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
subroutine sget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
SGET04
Definition: sget04.f:104
subroutine spbtrf(UPLO, N, KD, AB, LDAB, INFO)
SPBTRF
Definition: spbtrf.f:144
subroutine slaqsb(UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, EQUED)
SLAQSB scales a symmetric/Hermitian band matrix, using scaling factors computed by spbequ...
Definition: slaqsb.f:142
subroutine spbsv(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
SPBSV computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: spbsv.f:166
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine sdrvpb(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, IWORK, NOUT)
SDRVPB
Definition: sdrvpb.f:166
subroutine spbt05(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
SPBT05
Definition: spbt05.f:173
subroutine serrvx(PATH, NUNIT)
SERRVX
Definition: serrvx.f:57
subroutine spbt01(UPLO, N, KD, A, LDA, AFAC, LDAFAC, RWORK, RESID)
SPBT01
Definition: spbt01.f:121
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP
Definition: sswap.f:53
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53