LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvgbx.f
Go to the documentation of this file.
1*> \brief \b ZDRVGBX
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 ZDRVGB( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, LA,
12* AFB, LAFB, ASAV, B, BSAV, X, XACT, S, WORK,
13* RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER LA, LAFB, NN, NOUT, NRHS
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NVAL( * )
23* DOUBLE PRECISION RWORK( * ), S( * )
24* COMPLEX*16 A( * ), AFB( * ), ASAV( * ), B( * ), BSAV( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZDRVGB tests the driver routines ZGBSV, -SVX, and -SVXX.
35*>
36*> Note that this file is used only when the XBLAS are available,
37*> otherwise zdrvgb.f defines this subroutine.
38*> \endverbatim
39*
40* Arguments:
41* ==========
42*
43*> \param[in] DOTYPE
44*> \verbatim
45*> DOTYPE is LOGICAL array, dimension (NTYPES)
46*> The matrix types to be used for testing. Matrices of type j
47*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
48*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
49*> \endverbatim
50*>
51*> \param[in] NN
52*> \verbatim
53*> NN is INTEGER
54*> The number of values of N contained in the vector NVAL.
55*> \endverbatim
56*>
57*> \param[in] NVAL
58*> \verbatim
59*> NVAL is INTEGER array, dimension (NN)
60*> The values of the matrix column dimension N.
61*> \endverbatim
62*>
63*> \param[in] NRHS
64*> \verbatim
65*> NRHS is INTEGER
66*> The number of right hand side vectors to be generated for
67*> each linear system.
68*> \endverbatim
69*>
70*> \param[in] THRESH
71*> \verbatim
72*> THRESH is DOUBLE PRECISION
73*> The threshold value for the test ratios. A result is
74*> included in the output file if RESULT >= THRESH. To have
75*> every test ratio printed, use THRESH = 0.
76*> \endverbatim
77*>
78*> \param[in] TSTERR
79*> \verbatim
80*> TSTERR is LOGICAL
81*> Flag that indicates whether error exits are to be tested.
82*> \endverbatim
83*>
84*> \param[out] A
85*> \verbatim
86*> A is COMPLEX*16 array, dimension (LA)
87*> \endverbatim
88*>
89*> \param[in] LA
90*> \verbatim
91*> LA is INTEGER
92*> The length of the array A. LA >= (2*NMAX-1)*NMAX
93*> where NMAX is the largest entry in NVAL.
94*> \endverbatim
95*>
96*> \param[out] AFB
97*> \verbatim
98*> AFB is COMPLEX*16 array, dimension (LAFB)
99*> \endverbatim
100*>
101*> \param[in] LAFB
102*> \verbatim
103*> LAFB is INTEGER
104*> The length of the array AFB. LAFB >= (3*NMAX-2)*NMAX
105*> where NMAX is the largest entry in NVAL.
106*> \endverbatim
107*>
108*> \param[out] ASAV
109*> \verbatim
110*> ASAV is COMPLEX*16 array, dimension (LA)
111*> \endverbatim
112*>
113*> \param[out] B
114*> \verbatim
115*> B is COMPLEX*16 array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] BSAV
119*> \verbatim
120*> BSAV is COMPLEX*16 array, dimension (NMAX*NRHS)
121*> \endverbatim
122*>
123*> \param[out] X
124*> \verbatim
125*> X is COMPLEX*16 array, dimension (NMAX*NRHS)
126*> \endverbatim
127*>
128*> \param[out] XACT
129*> \verbatim
130*> XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
131*> \endverbatim
132*>
133*> \param[out] S
134*> \verbatim
135*> S is DOUBLE PRECISION array, dimension (2*NMAX)
136*> \endverbatim
137*>
138*> \param[out] WORK
139*> \verbatim
140*> WORK is COMPLEX*16 array, dimension
141*> (NMAX*max(3,NRHS,NMAX))
142*> \endverbatim
143*>
144*> \param[out] RWORK
145*> \verbatim
146*> RWORK is DOUBLE PRECISION array, dimension
147*> (max(2*NMAX,NMAX+2*NRHS))
148*> \endverbatim
149*>
150*> \param[out] IWORK
151*> \verbatim
152*> IWORK is INTEGER array, dimension (NMAX)
153*> \endverbatim
154*>
155*> \param[in] NOUT
156*> \verbatim
157*> NOUT is INTEGER
158*> The unit number for output.
159*> \endverbatim
160*
161* Authors:
162* ========
163*
164*> \author Univ. of Tennessee
165*> \author Univ. of California Berkeley
166*> \author Univ. of Colorado Denver
167*> \author NAG Ltd.
168*
169*> \ingroup complex16_lin
170*
171* =====================================================================
172 SUBROUTINE zdrvgb( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, LA,
173 $ AFB, LAFB, ASAV, B, BSAV, X, XACT, S, WORK,
174 $ RWORK, IWORK, NOUT )
175*
176* -- LAPACK test routine --
177* -- LAPACK is a software package provided by Univ. of Tennessee, --
178* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
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 RWORK( * ), S( * )
189 COMPLEX*16 A( * ), AFB( * ), ASAV( * ), B( * ), BSAV( * ),
190 $ work( * ), x( * ), 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 $ n_err_bnds
214 DOUBLE PRECISION AINVNM, AMAX, ANORM, ANORMI, ANORMO, ANRMPV,
215 $ CNDNUM, COLCND, RCOND, RCONDC, RCONDI, RCONDO,
216 $ roldc, roldi, roldo, rowcnd, rpvgrw,
217 $ rpvgrw_svxx
218* ..
219* .. Local Arrays ..
220 CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN )
221 INTEGER ISEED( 4 ), ISEEDY( 4 )
222 DOUBLE PRECISION RDUM( 1 ), RESULT( NTESTS ), BERR( NRHS ),
223 $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
224* ..
225* .. External Functions ..
226 LOGICAL LSAME
227 DOUBLE PRECISION DGET06, DLAMCH, ZLANGB, ZLANGE, ZLANTB,
229 EXTERNAL lsame, dget06, dlamch, zlangb, zlange, zlantb,
231* ..
232* .. External Subroutines ..
233 EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx, zgbequ,
237* ..
238* .. Intrinsic Functions ..
239 INTRINSIC abs, dcmplx, max, min
240* ..
241* .. Scalars in Common ..
242 LOGICAL LERR, OK
243 CHARACTER*32 SRNAMT
244 INTEGER INFOT, NUNIT
245* ..
246* .. Common blocks ..
247 COMMON / infoc / infot, nunit, ok, lerr
248 COMMON / srnamc / srnamt
249* ..
250* .. Data statements ..
251 DATA iseedy / 1988, 1989, 1990, 1991 /
252 DATA transs / 'N', 'T', 'C' /
253 DATA facts / 'F', 'N', 'E' /
254 DATA equeds / 'N', 'R', 'C', 'B' /
255* ..
256* .. Executable Statements ..
257*
258* Initialize constants and the random number seed.
259*
260 path( 1: 1 ) = 'Zomplex precision'
261 path( 2: 3 ) = 'GB'
262 nrun = 0
263 nfail = 0
264 nerrs = 0
265 DO 10 i = 1, 4
266 iseed( i ) = iseedy( i )
267 10 CONTINUE
268*
269* Test the error exits
270*
271 IF( tsterr )
272 $ CALL zerrvx( path, nout )
273 infot = 0
274*
275* Set the block size and minimum block size for testing.
276*
277 nb = 1
278 nbmin = 2
279 CALL xlaenv( 1, nb )
280 CALL xlaenv( 2, nbmin )
281*
282* Do for each value of N in NVAL
283*
284 DO 150 in = 1, nn
285 n = nval( in )
286 ldb = max( n, 1 )
287 xtype = 'N'
288*
289* Set limits on the number of loop iterations.
290*
291 nkl = max( 1, min( n, 4 ) )
292 IF( n.EQ.0 )
293 $ nkl = 1
294 nku = nkl
295 nimat = ntypes
296 IF( n.LE.0 )
297 $ nimat = 1
298*
299 DO 140 ikl = 1, nkl
300*
301* Do for KL = 0, N-1, (3N-1)/4, and (N+1)/4. This order makes
302* it easier to skip redundant values for small values of N.
303*
304 IF( ikl.EQ.1 ) THEN
305 kl = 0
306 ELSE IF( ikl.EQ.2 ) THEN
307 kl = max( n-1, 0 )
308 ELSE IF( ikl.EQ.3 ) THEN
309 kl = ( 3*n-1 ) / 4
310 ELSE IF( ikl.EQ.4 ) THEN
311 kl = ( n+1 ) / 4
312 END IF
313 DO 130 iku = 1, nku
314*
315* Do for KU = 0, N-1, (3N-1)/4, and (N+1)/4. This order
316* makes it easier to skip redundant values for small
317* values of N.
318*
319 IF( iku.EQ.1 ) THEN
320 ku = 0
321 ELSE IF( iku.EQ.2 ) THEN
322 ku = max( n-1, 0 )
323 ELSE IF( iku.EQ.3 ) THEN
324 ku = ( 3*n-1 ) / 4
325 ELSE IF( iku.EQ.4 ) THEN
326 ku = ( n+1 ) / 4
327 END IF
328*
329* Check that A and AFB are big enough to generate this
330* matrix.
331*
332 lda = kl + ku + 1
333 ldafb = 2*kl + ku + 1
334 IF( lda*n.GT.la .OR. ldafb*n.GT.lafb ) THEN
335 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
336 $ CALL aladhd( nout, path )
337 IF( lda*n.GT.la ) THEN
338 WRITE( nout, fmt = 9999 )la, n, kl, ku,
339 $ n*( kl+ku+1 )
340 nerrs = nerrs + 1
341 END IF
342 IF( ldafb*n.GT.lafb ) THEN
343 WRITE( nout, fmt = 9998 )lafb, n, kl, ku,
344 $ n*( 2*kl+ku+1 )
345 nerrs = nerrs + 1
346 END IF
347 GO TO 130
348 END IF
349*
350 DO 120 imat = 1, nimat
351*
352* Do the tests only if DOTYPE( IMAT ) is true.
353*
354 IF( .NOT.dotype( imat ) )
355 $ GO TO 120
356*
357* Skip types 2, 3, or 4 if the matrix is too small.
358*
359 zerot = imat.GE.2 .AND. imat.LE.4
360 IF( zerot .AND. n.LT.imat-1 )
361 $ GO TO 120
362*
363* Set up parameters with ZLATB4 and generate a
364* test matrix with ZLATMS.
365*
366 CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm,
367 $ mode, cndnum, dist )
368 rcondc = one / cndnum
369*
370 srnamt = 'ZLATMS'
371 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
372 $ cndnum, anorm, kl, ku, 'Z', a, lda, work,
373 $ info )
374*
375* Check the error code from ZLATMS.
376*
377 IF( info.NE.0 ) THEN
378 CALL alaerh( path, 'ZLATMS', info, 0, ' ', n, n,
379 $ kl, ku, -1, imat, nfail, nerrs, nout )
380 GO TO 120
381 END IF
382*
383* For types 2, 3, and 4, zero one or more columns of
384* the matrix to test that INFO is returned correctly.
385*
386 izero = 0
387 IF( zerot ) THEN
388 IF( imat.EQ.2 ) THEN
389 izero = 1
390 ELSE IF( imat.EQ.3 ) THEN
391 izero = n
392 ELSE
393 izero = n / 2 + 1
394 END IF
395 ioff = ( izero-1 )*lda
396 IF( imat.LT.4 ) THEN
397 i1 = max( 1, ku+2-izero )
398 i2 = min( kl+ku+1, ku+1+( n-izero ) )
399 DO 20 i = i1, i2
400 a( ioff+i ) = zero
401 20 CONTINUE
402 ELSE
403 DO 40 j = izero, n
404 DO 30 i = max( 1, ku+2-j ),
405 $ min( kl+ku+1, ku+1+( n-j ) )
406 a( ioff+i ) = zero
407 30 CONTINUE
408 ioff = ioff + lda
409 40 CONTINUE
410 END IF
411 END IF
412*
413* Save a copy of the matrix A in ASAV.
414*
415 CALL zlacpy( 'Full', kl+ku+1, n, a, lda, asav, lda )
416*
417 DO 110 iequed = 1, 4
418 equed = equeds( iequed )
419 IF( iequed.EQ.1 ) THEN
420 nfact = 3
421 ELSE
422 nfact = 1
423 END IF
424*
425 DO 100 ifact = 1, nfact
426 fact = facts( ifact )
427 prefac = lsame( fact, 'F' )
428 nofact = lsame( fact, 'N' )
429 equil = lsame( fact, 'E' )
430*
431 IF( zerot ) THEN
432 IF( prefac )
433 $ GO TO 100
434 rcondo = zero
435 rcondi = zero
436*
437 ELSE IF( .NOT.nofact ) THEN
438*
439* Compute the condition number for comparison
440* with the value returned by DGESVX (FACT =
441* 'N' reuses the condition number from the
442* previous iteration with FACT = 'F').
443*
444 CALL zlacpy( 'Full', kl+ku+1, n, asav, lda,
445 $ afb( kl+1 ), ldafb )
446 IF( equil .OR. iequed.GT.1 ) THEN
447*
448* Compute row and column scale factors to
449* equilibrate the matrix A.
450*
451 CALL zgbequ( n, n, kl, ku, afb( kl+1 ),
452 $ ldafb, s, s( n+1 ), rowcnd,
453 $ colcnd, amax, info )
454 IF( info.EQ.0 .AND. n.GT.0 ) THEN
455 IF( lsame( equed, 'R' ) ) THEN
456 rowcnd = zero
457 colcnd = one
458 ELSE IF( lsame( equed, 'C' ) ) THEN
459 rowcnd = one
460 colcnd = zero
461 ELSE IF( lsame( equed, 'B' ) ) THEN
462 rowcnd = zero
463 colcnd = zero
464 END IF
465*
466* Equilibrate the matrix.
467*
468 CALL zlaqgb( n, n, kl, ku, afb( kl+1 ),
469 $ ldafb, s, s( n+1 ),
470 $ rowcnd, colcnd, amax,
471 $ equed )
472 END IF
473 END IF
474*
475* Save the condition number of the
476* non-equilibrated system for use in ZGET04.
477*
478 IF( equil ) THEN
479 roldo = rcondo
480 roldi = rcondi
481 END IF
482*
483* Compute the 1-norm and infinity-norm of A.
484*
485 anormo = zlangb( '1', n, kl, ku, afb( kl+1 ),
486 $ ldafb, rwork )
487 anormi = zlangb( 'I', n, kl, ku, afb( kl+1 ),
488 $ ldafb, rwork )
489*
490* Factor the matrix A.
491*
492 CALL zgbtrf( n, n, kl, ku, afb, ldafb, iwork,
493 $ info )
494*
495* Form the inverse of A.
496*
497 CALL zlaset( 'Full', n, n, dcmplx( zero ),
498 $ dcmplx( one ), work, ldb )
499 srnamt = 'ZGBTRS'
500 CALL zgbtrs( 'No transpose', n, kl, ku, n,
501 $ afb, ldafb, iwork, work, ldb,
502 $ info )
503*
504* Compute the 1-norm condition number of A.
505*
506 ainvnm = zlange( '1', n, n, work, ldb,
507 $ rwork )
508 IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
509 rcondo = one
510 ELSE
511 rcondo = ( one / anormo ) / ainvnm
512 END IF
513*
514* Compute the infinity-norm condition number
515* of A.
516*
517 ainvnm = zlange( 'I', n, n, work, ldb,
518 $ rwork )
519 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
520 rcondi = one
521 ELSE
522 rcondi = ( one / anormi ) / ainvnm
523 END IF
524 END IF
525*
526 DO 90 itran = 1, ntran
527*
528* Do for each value of TRANS.
529*
530 trans = transs( itran )
531 IF( itran.EQ.1 ) THEN
532 rcondc = rcondo
533 ELSE
534 rcondc = rcondi
535 END IF
536*
537* Restore the matrix A.
538*
539 CALL zlacpy( 'Full', kl+ku+1, n, asav, lda,
540 $ a, lda )
541*
542* Form an exact solution and set the right hand
543* side.
544*
545 srnamt = 'ZLARHS'
546 CALL zlarhs( path, xtype, 'Full', trans, n,
547 $ n, kl, ku, nrhs, a, lda, xact,
548 $ ldb, b, ldb, iseed, info )
549 xtype = 'C'
550 CALL zlacpy( 'Full', n, nrhs, b, ldb, bsav,
551 $ ldb )
552*
553 IF( nofact .AND. itran.EQ.1 ) THEN
554*
555* --- Test ZGBSV ---
556*
557* Compute the LU factorization of the matrix
558* and solve the system.
559*
560 CALL zlacpy( 'Full', kl+ku+1, n, a, lda,
561 $ afb( kl+1 ), ldafb )
562 CALL zlacpy( 'Full', n, nrhs, b, ldb, x,
563 $ ldb )
564*
565 srnamt = 'ZGBSV '
566 CALL zgbsv( n, kl, ku, nrhs, afb, ldafb,
567 $ iwork, x, ldb, info )
568*
569* Check error code from ZGBSV .
570*
571 IF( info.NE.izero )
572 $ CALL alaerh( path, 'ZGBSV ', info,
573 $ izero, ' ', n, n, kl, ku,
574 $ nrhs, imat, nfail, nerrs,
575 $ nout )
576*
577* Reconstruct matrix from factors and
578* compute residual.
579*
580 CALL zgbt01( n, n, kl, ku, a, lda, afb,
581 $ ldafb, iwork, work,
582 $ result( 1 ) )
583 nt = 1
584 IF( izero.EQ.0 ) THEN
585*
586* Compute residual of the computed
587* solution.
588*
589 CALL zlacpy( 'Full', n, nrhs, b, ldb,
590 $ work, ldb )
591 CALL zgbt02( 'No transpose', n, n, kl,
592 $ ku, nrhs, a, lda, x, ldb,
593 $ work, ldb, rwork,
594 $ result( 2 ) )
595*
596* Check solution from generated exact
597* solution.
598*
599 CALL zget04( n, nrhs, x, ldb, xact,
600 $ ldb, rcondc, result( 3 ) )
601 nt = 3
602 END IF
603*
604* Print information about the tests that did
605* not pass the threshold.
606*
607 DO 50 k = 1, nt
608 IF( result( k ).GE.thresh ) THEN
609 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
610 $ CALL aladhd( nout, path )
611 WRITE( nout, fmt = 9997 )'ZGBSV ',
612 $ n, kl, ku, imat, k, result( k )
613 nfail = nfail + 1
614 END IF
615 50 CONTINUE
616 nrun = nrun + nt
617 END IF
618*
619* --- Test ZGBSVX ---
620*
621 IF( .NOT.prefac )
622 $ CALL zlaset( 'Full', 2*kl+ku+1, n,
623 $ dcmplx( zero ),
624 $ dcmplx( zero ), afb, ldafb )
625 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
626 $ dcmplx( zero ), x, ldb )
627 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
628*
629* Equilibrate the matrix if FACT = 'F' and
630* EQUED = 'R', 'C', or 'B'.
631*
632 CALL zlaqgb( n, n, kl, ku, a, lda, s,
633 $ s( n+1 ), rowcnd, colcnd,
634 $ amax, equed )
635 END IF
636*
637* Solve the system and compute the condition
638* number and error bounds using ZGBSVX.
639*
640 srnamt = 'ZGBSVX'
641 CALL zgbsvx( fact, trans, n, kl, ku, nrhs, a,
642 $ lda, afb, ldafb, iwork, equed,
643 $ s, s( ldb+1 ), b, ldb, x, ldb,
644 $ rcond, rwork, rwork( nrhs+1 ),
645 $ work, rwork( 2*nrhs+1 ), info )
646*
647* Check the error code from ZGBSVX.
648*
649 IF( info.NE.izero )
650 $ CALL alaerh( path, 'ZGBSVX', info, izero,
651 $ fact // trans, n, n, kl, ku,
652 $ nrhs, imat, nfail, nerrs,
653 $ nout )
654*
655* Compare RWORK(2*NRHS+1) from ZGBSVX with the
656* computed reciprocal pivot growth RPVGRW
657*
658 IF( info.NE.0 ) THEN
659 anrmpv = zero
660 DO 70 j = 1, info
661 DO 60 i = max( ku+2-j, 1 ),
662 $ min( n+ku+1-j, kl+ku+1 )
663 anrmpv = max( anrmpv,
664 $ abs( a( i+( j-1 )*lda ) ) )
665 60 CONTINUE
666 70 CONTINUE
667 rpvgrw = zlantb( 'M', 'U', 'N', info,
668 $ min( info-1, kl+ku ),
669 $ afb( max( 1, kl+ku+2-info ) ),
670 $ ldafb, rdum )
671 IF( rpvgrw.EQ.zero ) THEN
672 rpvgrw = one
673 ELSE
674 rpvgrw = anrmpv / rpvgrw
675 END IF
676 ELSE
677 rpvgrw = zlantb( 'M', 'U', 'N', n, kl+ku,
678 $ afb, ldafb, rdum )
679 IF( rpvgrw.EQ.zero ) THEN
680 rpvgrw = one
681 ELSE
682 rpvgrw = zlangb( 'M', n, kl, ku, a,
683 $ lda, rdum ) / rpvgrw
684 END IF
685 END IF
686 result( 7 ) = abs( rpvgrw-rwork( 2*nrhs+1 ) )
687 $ / max( rwork( 2*nrhs+1 ),
688 $ rpvgrw ) / dlamch( 'E' )
689*
690 IF( .NOT.prefac ) THEN
691*
692* Reconstruct matrix from factors and
693* compute residual.
694*
695 CALL zgbt01( n, n, kl, ku, a, lda, afb,
696 $ ldafb, iwork, work,
697 $ result( 1 ) )
698 k1 = 1
699 ELSE
700 k1 = 2
701 END IF
702*
703 IF( info.EQ.0 ) THEN
704 trfcon = .false.
705*
706* Compute residual of the computed solution.
707*
708 CALL zlacpy( 'Full', n, nrhs, bsav, ldb,
709 $ work, ldb )
710 CALL zgbt02( trans, n, n, kl, ku, nrhs,
711 $ asav, lda, x, ldb, work, ldb,
712 $ rwork( 2*nrhs+1 ),
713 $ result( 2 ) )
714*
715* Check solution from generated exact
716* solution.
717*
718 IF( nofact .OR. ( prefac .AND.
719 $ lsame( equed, 'N' ) ) ) THEN
720 CALL zget04( n, nrhs, x, ldb, xact,
721 $ ldb, rcondc, result( 3 ) )
722 ELSE
723 IF( itran.EQ.1 ) THEN
724 roldc = roldo
725 ELSE
726 roldc = roldi
727 END IF
728 CALL zget04( n, nrhs, x, ldb, xact,
729 $ ldb, roldc, result( 3 ) )
730 END IF
731*
732* Check the error bounds from iterative
733* refinement.
734*
735 CALL zgbt05( trans, n, kl, ku, nrhs, asav,
736 $ lda, bsav, ldb, x, ldb, xact,
737 $ ldb, rwork, rwork( nrhs+1 ),
738 $ result( 4 ) )
739 ELSE
740 trfcon = .true.
741 END IF
742*
743* Compare RCOND from ZGBSVX with the computed
744* value in RCONDC.
745*
746 result( 6 ) = dget06( rcond, rcondc )
747*
748* Print information about the tests that did
749* not pass the threshold.
750*
751 IF( .NOT.trfcon ) THEN
752 DO 80 k = k1, ntests
753 IF( result( k ).GE.thresh ) THEN
754 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
755 $ CALL aladhd( nout, path )
756 IF( prefac ) THEN
757 WRITE( nout, fmt = 9995 )
758 $ 'ZGBSVX', fact, trans, n, kl,
759 $ ku, equed, imat, k,
760 $ result( k )
761 ELSE
762 WRITE( nout, fmt = 9996 )
763 $ 'ZGBSVX', fact, trans, n, kl,
764 $ ku, imat, k, result( k )
765 END IF
766 nfail = nfail + 1
767 END IF
768 80 CONTINUE
769 nrun = nrun + 7 - k1
770 ELSE
771 IF( result( 1 ).GE.thresh .AND. .NOT.
772 $ prefac ) THEN
773 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
774 $ CALL aladhd( nout, path )
775 IF( prefac ) THEN
776 WRITE( nout, fmt = 9995 )'ZGBSVX',
777 $ fact, trans, n, kl, ku, equed,
778 $ imat, 1, result( 1 )
779 ELSE
780 WRITE( nout, fmt = 9996 )'ZGBSVX',
781 $ fact, trans, n, kl, ku, imat, 1,
782 $ result( 1 )
783 END IF
784 nfail = nfail + 1
785 nrun = nrun + 1
786 END IF
787 IF( result( 6 ).GE.thresh ) THEN
788 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
789 $ CALL aladhd( nout, path )
790 IF( prefac ) THEN
791 WRITE( nout, fmt = 9995 )'ZGBSVX',
792 $ fact, trans, n, kl, ku, equed,
793 $ imat, 6, result( 6 )
794 ELSE
795 WRITE( nout, fmt = 9996 )'ZGBSVX',
796 $ fact, trans, n, kl, ku, imat, 6,
797 $ result( 6 )
798 END IF
799 nfail = nfail + 1
800 nrun = nrun + 1
801 END IF
802 IF( result( 7 ).GE.thresh ) THEN
803 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
804 $ CALL aladhd( nout, path )
805 IF( prefac ) THEN
806 WRITE( nout, fmt = 9995 )'ZGBSVX',
807 $ fact, trans, n, kl, ku, equed,
808 $ imat, 7, result( 7 )
809 ELSE
810 WRITE( nout, fmt = 9996 )'ZGBSVX',
811 $ fact, trans, n, kl, ku, imat, 7,
812 $ result( 7 )
813 END IF
814 nfail = nfail + 1
815 nrun = nrun + 1
816 END IF
817 END IF
818
819* --- Test ZGBSVXX ---
820
821* Restore the matrices A and B.
822
823c write(*,*) 'begin zgbsvxx testing'
824
825 CALL zlacpy( 'Full', kl+ku+1, n, asav, lda, a,
826 $ lda )
827 CALL zlacpy( 'Full', n, nrhs, bsav, ldb, b, ldb )
828
829 IF( .NOT.prefac )
830 $ CALL zlaset( 'Full', 2*kl+ku+1, n,
831 $ dcmplx( zero ), dcmplx( zero ),
832 $ afb, ldafb )
833 CALL zlaset( 'Full', n, nrhs,
834 $ dcmplx( zero ), dcmplx( zero ),
835 $ x, ldb )
836 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
837*
838* Equilibrate the matrix if FACT = 'F' and
839* EQUED = 'R', 'C', or 'B'.
840*
841 CALL zlaqgb( n, n, kl, ku, a, lda, s,
842 $ s( n+1 ), rowcnd, colcnd, amax, equed )
843 END IF
844*
845* Solve the system and compute the condition number
846* and error bounds using ZGBSVXX.
847*
848 srnamt = 'ZGBSVXX'
849 n_err_bnds = 3
850 CALL zgbsvxx( fact, trans, n, kl, ku, nrhs, a, lda,
851 $ afb, ldafb, iwork, equed, s, s( n+1 ), b, ldb,
852 $ x, ldb, rcond, rpvgrw_svxx, berr, n_err_bnds,
853 $ errbnds_n, errbnds_c, 0, zero, work,
854 $ rwork, info )
855*
856* Check the error code from ZGBSVXX.
857*
858 IF( info.EQ.n+1 ) GOTO 90
859 IF( info.NE.izero ) THEN
860 CALL alaerh( path, 'ZGBSVXX', info, izero,
861 $ fact // trans, n, n, -1, -1, nrhs,
862 $ imat, nfail, nerrs, nout )
863 GOTO 90
864 END IF
865*
866* Compare rpvgrw_svxx from ZGESVXX with the computed
867* reciprocal pivot growth factor RPVGRW
868*
869
870 IF ( info .GT. 0 .AND. info .LT. n+1 ) THEN
871 rpvgrw = zla_gbrpvgrw(n, kl, ku, info, a, lda,
872 $ afb, ldafb)
873 ELSE
874 rpvgrw = zla_gbrpvgrw(n, kl, ku, n, a, lda,
875 $ afb, ldafb)
876 ENDIF
877
878 result( 7 ) = abs( rpvgrw-rpvgrw_svxx ) /
879 $ max( rpvgrw_svxx, rpvgrw ) /
880 $ dlamch( 'E' )
881*
882 IF( .NOT.prefac ) THEN
883*
884* Reconstruct matrix from factors and compute
885* residual.
886*
887 CALL zgbt01( n, n, kl, ku, a, lda, afb, ldafb,
888 $ iwork, work( 2*nrhs+1 ), result( 1 ) )
889 k1 = 1
890 ELSE
891 k1 = 2
892 END IF
893*
894 IF( info.EQ.0 ) THEN
895 trfcon = .false.
896*
897* Compute residual of the computed solution.
898*
899 CALL zlacpy( 'Full', n, nrhs, bsav, ldb, work,
900 $ ldb )
901 CALL zgbt02( trans, n, n, kl, ku, nrhs, asav,
902 $ lda, x, ldb, work, ldb, rwork,
903 $ result( 2 ) )
904*
905* Check solution from generated exact solution.
906*
907 IF( nofact .OR. ( prefac .AND. lsame( equed,
908 $ 'N' ) ) ) THEN
909 CALL zget04( n, nrhs, x, ldb, xact, ldb,
910 $ rcondc, result( 3 ) )
911 ELSE
912 IF( itran.EQ.1 ) THEN
913 roldc = roldo
914 ELSE
915 roldc = roldi
916 END IF
917 CALL zget04( n, nrhs, x, ldb, xact, ldb,
918 $ roldc, result( 3 ) )
919 END IF
920 ELSE
921 trfcon = .true.
922 END IF
923*
924* Compare RCOND from ZGBSVXX with the computed value
925* in RCONDC.
926*
927 result( 6 ) = dget06( rcond, rcondc )
928*
929* Print information about the tests that did not pass
930* the threshold.
931*
932 IF( .NOT.trfcon ) THEN
933 DO 45 k = k1, ntests
934 IF( result( k ).GE.thresh ) THEN
935 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
936 $ CALL aladhd( nout, path )
937 IF( prefac ) THEN
938 WRITE( nout, fmt = 9995 )'ZGBSVXX',
939 $ fact, trans, n, kl, ku, equed,
940 $ imat, k, result( k )
941 ELSE
942 WRITE( nout, fmt = 9996 )'ZGBSVXX',
943 $ fact, trans, n, kl, ku, imat, k,
944 $ result( k )
945 END IF
946 nfail = nfail + 1
947 END IF
948 45 CONTINUE
949 nrun = nrun + 7 - k1
950 ELSE
951 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
952 $ THEN
953 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
954 $ CALL aladhd( nout, path )
955 IF( prefac ) THEN
956 WRITE( nout, fmt = 9995 )'ZGBSVXX', fact,
957 $ trans, n, kl, ku, equed, imat, 1,
958 $ result( 1 )
959 ELSE
960 WRITE( nout, fmt = 9996 )'ZGBSVXX', fact,
961 $ trans, n, kl, ku, imat, 1,
962 $ result( 1 )
963 END IF
964 nfail = nfail + 1
965 nrun = nrun + 1
966 END IF
967 IF( result( 6 ).GE.thresh ) THEN
968 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
969 $ CALL aladhd( nout, path )
970 IF( prefac ) THEN
971 WRITE( nout, fmt = 9995 )'ZGBSVXX', fact,
972 $ trans, n, kl, ku, equed, imat, 6,
973 $ result( 6 )
974 ELSE
975 WRITE( nout, fmt = 9996 )'ZGBSVXX', fact,
976 $ trans, n, kl, ku, imat, 6,
977 $ result( 6 )
978 END IF
979 nfail = nfail + 1
980 nrun = nrun + 1
981 END IF
982 IF( result( 7 ).GE.thresh ) THEN
983 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
984 $ CALL aladhd( nout, path )
985 IF( prefac ) THEN
986 WRITE( nout, fmt = 9995 )'ZGBSVXX', fact,
987 $ trans, n, kl, ku, equed, imat, 7,
988 $ result( 7 )
989 ELSE
990 WRITE( nout, fmt = 9996 )'ZGBSVXX', fact,
991 $ trans, n, kl, ku, imat, 7,
992 $ result( 7 )
993 END IF
994 nfail = nfail + 1
995 nrun = nrun + 1
996 END IF
997*
998 END IF
999*
1000 90 CONTINUE
1001 100 CONTINUE
1002 110 CONTINUE
1003 120 CONTINUE
1004 130 CONTINUE
1005 140 CONTINUE
1006 150 CONTINUE
1007*
1008* Print a summary of the results.
1009*
1010 CALL alasvm( path, nout, nfail, nrun, nerrs )
1011*
1012
1013* Test Error Bounds from ZGBSVXX
1014
1015 CALL zebchvxx(thresh, path)
1016
1017 9999 FORMAT( ' *** In ZDRVGB, LA=', i5, ' is too small for N=', i5,
1018 $ ', KU=', i5, ', KL=', i5, / ' ==> Increase LA to at least ',
1019 $ i5 )
1020 9998 FORMAT( ' *** In ZDRVGB, LAFB=', i5, ' is too small for N=', i5,
1021 $ ', KU=', i5, ', KL=', i5, /
1022 $ ' ==> Increase LAFB to at least ', i5 )
1023 9997 FORMAT( 1x, a, ', N=', i5, ', KL=', i5, ', KU=', i5, ', type ',
1024 $ i1, ', test(', i1, ')=', g12.5 )
1025 9996 FORMAT( 1x, a, '( ''', a1, ''',''', a1, ''',', i5, ',', i5, ',',
1026 $ i5, ',...), type ', i1, ', test(', i1, ')=', g12.5 )
1027 9995 FORMAT( 1x, a, '( ''', a1, ''',''', a1, ''',', i5, ',', i5, ',',
1028 $ i5, ',...), EQUED=''', a1, ''', type ', i1, ', test(', i1,
1029 $ ')=', g12.5 )
1030*
1031 RETURN
1032*
1033* End of ZDRVGBX
1034*
1035 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.f:208
subroutine aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine zgbequ(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, info)
ZGBEQU
Definition zgbequ.f:154
subroutine zgbsv(n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
ZGBSV computes the solution to system of linear equations A * X = B for GB matrices (simple driver)
Definition zgbsv.f:162
subroutine zgbsvx(fact, trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, equed, r, c, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZGBSVX computes the solution to system of linear equations A * X = B for GB matrices
Definition zgbsvx.f:370
subroutine zgbsvxx(fact, trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, equed, r, c, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
ZGBSVXX computes the solution to system of linear equations A * X = B for GB matrices
Definition zgbsvxx.f:560
subroutine zgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
ZGBTRF
Definition zgbtrf.f:144
subroutine zgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
ZGBTRS
Definition zgbtrs.f:138
double precision function zla_gbrpvgrw(n, kl, ku, ncols, ab, ldab, afb, ldafb)
ZLA_GBRPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a general banded matrix.
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlaqgb(m, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd, amax, equed)
ZLAQGB scales a general band matrix, using row and column scaling factors computed by sgbequ.
Definition zlaqgb.f:160
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zdrvgb(dotype, nn, nval, nrhs, thresh, tsterr, a, la, afb, lafb, asav, b, bsav, x, xact, s, work, rwork, iwork, nout)
ZDRVGB
Definition zdrvgb.f:172
subroutine zebchvxx(thresh, path)
ZEBCHVXX
Definition zebchvxx.f:96
subroutine zerrvx(path, nunit)
ZERRVX
Definition zerrvx.f:55
subroutine zgbt01(m, n, kl, ku, a, lda, afac, ldafac, ipiv, work, resid)
ZGBT01
Definition zgbt01.f:126
subroutine zgbt02(trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZGBT02
Definition zgbt02.f:148
subroutine zgbt05(trans, n, kl, ku, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZGBT05
Definition zgbt05.f:176
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332