LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvpox.f
Go to the documentation of this file.
1*> \brief \b ZDRVPOX
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 ZDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12* A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
13* RWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NOUT, NRHS
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER NVAL( * )
23* DOUBLE PRECISION RWORK( * ), S( * )
24* COMPLEX*16 A( * ), AFAC( * ), ASAV( * ), B( * ),
25* $ BSAV( * ), WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZDRVPO tests the driver routines ZPOSV, -SVX, and -SVXX.
35*>
36*> Note that this file is used only when the XBLAS are available,
37*> otherwise zdrvpo.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 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[in] NMAX
85*> \verbatim
86*> NMAX is INTEGER
87*> The maximum value permitted for N, used in dimensioning the
88*> work arrays.
89*> \endverbatim
90*>
91*> \param[out] A
92*> \verbatim
93*> A is COMPLEX*16 array, dimension (NMAX*NMAX)
94*> \endverbatim
95*>
96*> \param[out] AFAC
97*> \verbatim
98*> AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
99*> \endverbatim
100*>
101*> \param[out] ASAV
102*> \verbatim
103*> ASAV is COMPLEX*16 array, dimension (NMAX*NMAX)
104*> \endverbatim
105*>
106*> \param[out] B
107*> \verbatim
108*> B is COMPLEX*16 array, dimension (NMAX*NRHS)
109*> \endverbatim
110*>
111*> \param[out] BSAV
112*> \verbatim
113*> BSAV is COMPLEX*16 array, dimension (NMAX*NRHS)
114*> \endverbatim
115*>
116*> \param[out] X
117*> \verbatim
118*> X is COMPLEX*16 array, dimension (NMAX*NRHS)
119*> \endverbatim
120*>
121*> \param[out] XACT
122*> \verbatim
123*> XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
124*> \endverbatim
125*>
126*> \param[out] S
127*> \verbatim
128*> S is DOUBLE PRECISION array, dimension (NMAX)
129*> \endverbatim
130*>
131*> \param[out] WORK
132*> \verbatim
133*> WORK is COMPLEX*16 array, dimension
134*> (NMAX*max(3,NRHS))
135*> \endverbatim
136*>
137*> \param[out] RWORK
138*> \verbatim
139*> RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
140*> \endverbatim
141*>
142*> \param[in] NOUT
143*> \verbatim
144*> NOUT is INTEGER
145*> The unit number for output.
146*> \endverbatim
147*
148* Authors:
149* ========
150*
151*> \author Univ. of Tennessee
152*> \author Univ. of California Berkeley
153*> \author Univ. of Colorado Denver
154*> \author NAG Ltd.
155*
156*> \ingroup complex16_lin
157*
158* =====================================================================
159 SUBROUTINE zdrvpo( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
160 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
161 $ RWORK, NOUT )
162*
163* -- LAPACK test routine --
164* -- LAPACK is a software package provided by Univ. of Tennessee, --
165* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
166*
167* .. Scalar Arguments ..
168 LOGICAL TSTERR
169 INTEGER NMAX, NN, NOUT, NRHS
170 DOUBLE PRECISION THRESH
171* ..
172* .. Array Arguments ..
173 LOGICAL DOTYPE( * )
174 INTEGER NVAL( * )
175 DOUBLE PRECISION RWORK( * ), S( * )
176 COMPLEX*16 A( * ), AFAC( * ), ASAV( * ), B( * ),
177 $ bsav( * ), work( * ), x( * ), xact( * )
178* ..
179*
180* =====================================================================
181*
182* .. Parameters ..
183 DOUBLE PRECISION ONE, ZERO
184 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
185 INTEGER NTYPES
186 parameter( ntypes = 9 )
187 INTEGER NTESTS
188 parameter( ntests = 6 )
189* ..
190* .. Local Scalars ..
191 LOGICAL EQUIL, NOFACT, PREFAC, ZEROT
192 CHARACTER DIST, EQUED, FACT, TYPE, UPLO, XTYPE
193 CHARACTER*3 PATH
194 INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
195 $ izero, k, k1, kl, ku, lda, mode, n, nb, nbmin,
196 $ nerrs, nfact, nfail, nimat, nrun, nt,
197 $ n_err_bnds
198 DOUBLE PRECISION AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
199 $ ROLDC, SCOND, RPVGRW_SVXX
200* ..
201* .. Local Arrays ..
202 CHARACTER EQUEDS( 2 ), FACTS( 3 ), UPLOS( 2 )
203 INTEGER ISEED( 4 ), ISEEDY( 4 )
204 DOUBLE PRECISION RESULT( NTESTS ), BERR( NRHS ),
205 $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
206* ..
207* .. External Functions ..
208 LOGICAL LSAME
209 DOUBLE PRECISION DGET06, ZLANHE
210 EXTERNAL lsame, dget06, zlanhe
211* ..
212* .. External Subroutines ..
213 EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx, zget04,
217* ..
218* .. Scalars in Common ..
219 LOGICAL LERR, OK
220 CHARACTER*32 SRNAMT
221 INTEGER INFOT, NUNIT
222* ..
223* .. Common blocks ..
224 COMMON / infoc / infot, nunit, ok, lerr
225 COMMON / srnamc / srnamt
226* ..
227* .. Intrinsic Functions ..
228 INTRINSIC dcmplx, max
229* ..
230* .. Data statements ..
231 DATA iseedy / 1988, 1989, 1990, 1991 /
232 DATA uplos / 'U', 'L' /
233 DATA facts / 'F', 'N', 'E' /
234 DATA equeds / 'N', 'Y' /
235* ..
236* .. Executable Statements ..
237*
238* Initialize constants and the random number seed.
239*
240 path( 1: 1 ) = 'Zomplex precision'
241 path( 2: 3 ) = 'PO'
242 nrun = 0
243 nfail = 0
244 nerrs = 0
245 DO 10 i = 1, 4
246 iseed( i ) = iseedy( i )
247 10 CONTINUE
248*
249* Test the error exits
250*
251 IF( tsterr )
252 $ CALL zerrvx( path, nout )
253 infot = 0
254*
255* Set the block size and minimum block size for testing.
256*
257 nb = 1
258 nbmin = 2
259 CALL xlaenv( 1, nb )
260 CALL xlaenv( 2, nbmin )
261*
262* Do for each value of N in NVAL
263*
264 DO 130 in = 1, nn
265 n = nval( in )
266 lda = max( n, 1 )
267 xtype = 'N'
268 nimat = ntypes
269 IF( n.LE.0 )
270 $ nimat = 1
271*
272 DO 120 imat = 1, nimat
273*
274* Do the tests only if DOTYPE( IMAT ) is true.
275*
276 IF( .NOT.dotype( imat ) )
277 $ GO TO 120
278*
279* Skip types 3, 4, or 5 if the matrix size is too small.
280*
281 zerot = imat.GE.3 .AND. imat.LE.5
282 IF( zerot .AND. n.LT.imat-2 )
283 $ GO TO 120
284*
285* Do first for UPLO = 'U', then for UPLO = 'L'
286*
287 DO 110 iuplo = 1, 2
288 uplo = uplos( iuplo )
289*
290* Set up parameters with ZLATB4 and generate a test matrix
291* with ZLATMS.
292*
293 CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
294 $ cndnum, dist )
295*
296 srnamt = 'ZLATMS'
297 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
298 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
299 $ info )
300*
301* Check error code from ZLATMS.
302*
303 IF( info.NE.0 ) THEN
304 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
305 $ -1, -1, imat, nfail, nerrs, nout )
306 GO TO 110
307 END IF
308*
309* For types 3-5, zero one row and column of the matrix to
310* test that INFO is returned correctly.
311*
312 IF( zerot ) THEN
313 IF( imat.EQ.3 ) THEN
314 izero = 1
315 ELSE IF( imat.EQ.4 ) THEN
316 izero = n
317 ELSE
318 izero = n / 2 + 1
319 END IF
320 ioff = ( izero-1 )*lda
321*
322* Set row and column IZERO of A to 0.
323*
324 IF( iuplo.EQ.1 ) THEN
325 DO 20 i = 1, izero - 1
326 a( ioff+i ) = zero
327 20 CONTINUE
328 ioff = ioff + izero
329 DO 30 i = izero, n
330 a( ioff ) = zero
331 ioff = ioff + lda
332 30 CONTINUE
333 ELSE
334 ioff = izero
335 DO 40 i = 1, izero - 1
336 a( ioff ) = zero
337 ioff = ioff + lda
338 40 CONTINUE
339 ioff = ioff - izero
340 DO 50 i = izero, n
341 a( ioff+i ) = zero
342 50 CONTINUE
343 END IF
344 ELSE
345 izero = 0
346 END IF
347*
348* Set the imaginary part of the diagonals.
349*
350 CALL zlaipd( n, a, lda+1, 0 )
351*
352* Save a copy of the matrix A in ASAV.
353*
354 CALL zlacpy( uplo, n, n, a, lda, asav, lda )
355*
356 DO 100 iequed = 1, 2
357 equed = equeds( iequed )
358 IF( iequed.EQ.1 ) THEN
359 nfact = 3
360 ELSE
361 nfact = 1
362 END IF
363*
364 DO 90 ifact = 1, nfact
365 fact = facts( ifact )
366 prefac = lsame( fact, 'F' )
367 nofact = lsame( fact, 'N' )
368 equil = lsame( fact, 'E' )
369*
370 IF( zerot ) THEN
371 IF( prefac )
372 $ GO TO 90
373 rcondc = zero
374*
375 ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
376*
377* Compute the condition number for comparison with
378* the value returned by ZPOSVX (FACT = 'N' reuses
379* the condition number from the previous iteration
380* with FACT = 'F').
381*
382 CALL zlacpy( uplo, n, n, asav, lda, afac, lda )
383 IF( equil .OR. iequed.GT.1 ) THEN
384*
385* Compute row and column scale factors to
386* equilibrate the matrix A.
387*
388 CALL zpoequ( n, afac, lda, s, scond, amax,
389 $ info )
390 IF( info.EQ.0 .AND. n.GT.0 ) THEN
391 IF( iequed.GT.1 )
392 $ scond = zero
393*
394* Equilibrate the matrix.
395*
396 CALL zlaqhe( uplo, n, afac, lda, s, scond,
397 $ amax, equed )
398 END IF
399 END IF
400*
401* Save the condition number of the
402* non-equilibrated system for use in ZGET04.
403*
404 IF( equil )
405 $ roldc = rcondc
406*
407* Compute the 1-norm of A.
408*
409 anorm = zlanhe( '1', uplo, n, afac, lda, rwork )
410*
411* Factor the matrix A.
412*
413 CALL zpotrf( uplo, n, afac, lda, info )
414*
415* Form the inverse of A.
416*
417 CALL zlacpy( uplo, n, n, afac, lda, a, lda )
418 CALL zpotri( uplo, n, a, lda, info )
419*
420* Compute the 1-norm condition number of A.
421*
422 ainvnm = zlanhe( '1', uplo, n, a, lda, rwork )
423 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
424 rcondc = one
425 ELSE
426 rcondc = ( one / anorm ) / ainvnm
427 END IF
428 END IF
429*
430* Restore the matrix A.
431*
432 CALL zlacpy( uplo, n, n, asav, lda, a, lda )
433*
434* Form an exact solution and set the right hand side.
435*
436 srnamt = 'ZLARHS'
437 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
438 $ nrhs, a, lda, xact, lda, b, lda,
439 $ iseed, info )
440 xtype = 'C'
441 CALL zlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
442*
443 IF( nofact ) THEN
444*
445* --- Test ZPOSV ---
446*
447* Compute the L*L' or U'*U factorization of the
448* matrix and solve the system.
449*
450 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
451 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
452*
453 srnamt = 'ZPOSV '
454 CALL zposv( uplo, n, nrhs, afac, lda, x, lda,
455 $ info )
456*
457* Check error code from ZPOSV .
458*
459 IF( info.NE.izero ) THEN
460 CALL alaerh( path, 'ZPOSV ', info, izero,
461 $ uplo, n, n, -1, -1, nrhs, imat,
462 $ nfail, nerrs, nout )
463 GO TO 70
464 ELSE IF( info.NE.0 ) THEN
465 GO TO 70
466 END IF
467*
468* Reconstruct matrix from factors and compute
469* residual.
470*
471 CALL zpot01( uplo, n, a, lda, afac, lda, rwork,
472 $ result( 1 ) )
473*
474* Compute residual of the computed solution.
475*
476 CALL zlacpy( 'Full', n, nrhs, b, lda, work,
477 $ lda )
478 CALL zpot02( uplo, n, nrhs, a, lda, x, lda,
479 $ work, lda, rwork, result( 2 ) )
480*
481* Check solution from generated exact solution.
482*
483 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
484 $ result( 3 ) )
485 nt = 3
486*
487* Print information about the tests that did not
488* pass the threshold.
489*
490 DO 60 k = 1, nt
491 IF( result( k ).GE.thresh ) THEN
492 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
493 $ CALL aladhd( nout, path )
494 WRITE( nout, fmt = 9999 )'ZPOSV ', uplo,
495 $ n, imat, k, result( k )
496 nfail = nfail + 1
497 END IF
498 60 CONTINUE
499 nrun = nrun + nt
500 70 CONTINUE
501 END IF
502*
503* --- Test ZPOSVX ---
504*
505 IF( .NOT.prefac )
506 $ CALL zlaset( uplo, n, n, dcmplx( zero ),
507 $ dcmplx( zero ), afac, lda )
508 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
509 $ dcmplx( zero ), x, lda )
510 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
511*
512* Equilibrate the matrix if FACT='F' and
513* EQUED='Y'.
514*
515 CALL zlaqhe( uplo, n, a, lda, s, scond, amax,
516 $ equed )
517 END IF
518*
519* Solve the system and compute the condition number
520* and error bounds using ZPOSVX.
521*
522 srnamt = 'ZPOSVX'
523 CALL zposvx( fact, uplo, n, nrhs, a, lda, afac,
524 $ lda, equed, s, b, lda, x, lda, rcond,
525 $ rwork, rwork( nrhs+1 ), work,
526 $ rwork( 2*nrhs+1 ), info )
527*
528* Check the error code from ZPOSVX.
529*
530 IF( info.NE.izero ) THEN
531 CALL alaerh( path, 'ZPOSVX', info, izero,
532 $ fact // uplo, n, n, -1, -1, nrhs,
533 $ imat, nfail, nerrs, nout )
534 GO TO 90
535 END IF
536*
537 IF( info.EQ.0 ) THEN
538 IF( .NOT.prefac ) THEN
539*
540* Reconstruct matrix from factors and compute
541* residual.
542*
543 CALL zpot01( uplo, n, a, lda, afac, lda,
544 $ rwork( 2*nrhs+1 ), result( 1 ) )
545 k1 = 1
546 ELSE
547 k1 = 2
548 END IF
549*
550* Compute residual of the computed solution.
551*
552 CALL zlacpy( 'Full', n, nrhs, bsav, lda, work,
553 $ lda )
554 CALL zpot02( uplo, n, nrhs, asav, lda, x, lda,
555 $ work, lda, rwork( 2*nrhs+1 ),
556 $ result( 2 ) )
557*
558* Check solution from generated exact solution.
559*
560 IF( nofact .OR. ( prefac .AND. lsame( equed,
561 $ 'N' ) ) ) THEN
562 CALL zget04( n, nrhs, x, lda, xact, lda,
563 $ rcondc, result( 3 ) )
564 ELSE
565 CALL zget04( n, nrhs, x, lda, xact, lda,
566 $ roldc, result( 3 ) )
567 END IF
568*
569* Check the error bounds from iterative
570* refinement.
571*
572 CALL zpot05( uplo, n, nrhs, asav, lda, b, lda,
573 $ x, lda, xact, lda, rwork,
574 $ rwork( nrhs+1 ), result( 4 ) )
575 ELSE
576 k1 = 6
577 END IF
578*
579* Compare RCOND from ZPOSVX with the computed value
580* in RCONDC.
581*
582 result( 6 ) = dget06( rcond, rcondc )
583*
584* Print information about the tests that did not pass
585* the threshold.
586*
587 DO 80 k = k1, 6
588 IF( result( k ).GE.thresh ) THEN
589 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
590 $ CALL aladhd( nout, path )
591 IF( prefac ) THEN
592 WRITE( nout, fmt = 9997 )'ZPOSVX', fact,
593 $ uplo, n, equed, imat, k, result( k )
594 ELSE
595 WRITE( nout, fmt = 9998 )'ZPOSVX', fact,
596 $ uplo, n, imat, k, result( k )
597 END IF
598 nfail = nfail + 1
599 END IF
600 80 CONTINUE
601 nrun = nrun + 7 - k1
602*
603* --- Test ZPOSVXX ---
604*
605* Restore the matrices A and B.
606*
607 CALL zlacpy( 'Full', n, n, asav, lda, a, lda )
608 CALL zlacpy( 'Full', n, nrhs, bsav, lda, b, lda )
609
610 IF( .NOT.prefac )
611 $ CALL zlaset( uplo, n, n, dcmplx( zero ),
612 $ dcmplx( zero ), afac, lda )
613 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
614 $ dcmplx( zero ), x, lda )
615 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
616*
617* Equilibrate the matrix if FACT='F' and
618* EQUED='Y'.
619*
620 CALL zlaqhe( uplo, n, a, lda, s, scond, amax,
621 $ equed )
622 END IF
623*
624* Solve the system and compute the condition number
625* and error bounds using ZPOSVXX.
626*
627 srnamt = 'ZPOSVXX'
628 n_err_bnds = 3
629 CALL zposvxx( fact, uplo, n, nrhs, a, lda, afac,
630 $ lda, equed, s, b, lda, x,
631 $ lda, rcond, rpvgrw_svxx, berr, n_err_bnds,
632 $ errbnds_n, errbnds_c, 0, zero, work,
633 $ rwork( 2*nrhs+1 ), info )
634*
635* Check the error code from ZPOSVXX.
636*
637 IF( info.EQ.n+1 ) GOTO 90
638 IF( info.NE.izero ) THEN
639 CALL alaerh( path, 'ZPOSVXX', info, izero,
640 $ fact // uplo, n, n, -1, -1, nrhs,
641 $ imat, nfail, nerrs, nout )
642 GO TO 90
643 END IF
644*
645 IF( info.EQ.0 ) THEN
646 IF( .NOT.prefac ) THEN
647*
648* Reconstruct matrix from factors and compute
649* residual.
650*
651 CALL zpot01( uplo, n, a, lda, afac, lda,
652 $ rwork( 2*nrhs+1 ), result( 1 ) )
653 k1 = 1
654 ELSE
655 k1 = 2
656 END IF
657*
658* Compute residual of the computed solution.
659*
660 CALL zlacpy( 'Full', n, nrhs, bsav, lda, work,
661 $ lda )
662 CALL zpot02( uplo, n, nrhs, asav, lda, x, lda,
663 $ work, lda, rwork( 2*nrhs+1 ),
664 $ result( 2 ) )
665*
666* Check solution from generated exact solution.
667*
668 IF( nofact .OR. ( prefac .AND. lsame( equed,
669 $ 'N' ) ) ) THEN
670 CALL zget04( n, nrhs, x, lda, xact, lda,
671 $ rcondc, result( 3 ) )
672 ELSE
673 CALL zget04( n, nrhs, x, lda, xact, lda,
674 $ roldc, result( 3 ) )
675 END IF
676*
677* Check the error bounds from iterative
678* refinement.
679*
680 CALL zpot05( uplo, n, nrhs, asav, lda, b, lda,
681 $ x, lda, xact, lda, rwork,
682 $ rwork( nrhs+1 ), result( 4 ) )
683 ELSE
684 k1 = 6
685 END IF
686*
687* Compare RCOND from ZPOSVXX with the computed value
688* in RCONDC.
689*
690 result( 6 ) = dget06( rcond, rcondc )
691*
692* Print information about the tests that did not pass
693* the threshold.
694*
695 DO 85 k = k1, 6
696 IF( result( k ).GE.thresh ) THEN
697 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
698 $ CALL aladhd( nout, path )
699 IF( prefac ) THEN
700 WRITE( nout, fmt = 9997 )'ZPOSVXX', fact,
701 $ uplo, n, equed, imat, k, result( k )
702 ELSE
703 WRITE( nout, fmt = 9998 )'ZPOSVXX', fact,
704 $ uplo, n, imat, k, result( k )
705 END IF
706 nfail = nfail + 1
707 END IF
708 85 CONTINUE
709 nrun = nrun + 7 - k1
710 90 CONTINUE
711 100 CONTINUE
712 110 CONTINUE
713 120 CONTINUE
714 130 CONTINUE
715*
716* Print a summary of the results.
717*
718 CALL alasvm( path, nout, nfail, nrun, nerrs )
719*
720
721* Test Error Bounds for ZGESVXX
722
723 CALL zebchvxx(thresh, path)
724
725 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
726 $ ', test(', i1, ')=', g12.5 )
727 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
728 $ ', type ', i1, ', test(', i1, ')=', g12.5 )
729 9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
730 $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ') =',
731 $ g12.5 )
732 RETURN
733*
734* End of ZDRVPOX
735*
736 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 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 zlaqhe(uplo, n, a, lda, s, scond, amax, equed)
ZLAQHE scales a Hermitian matrix.
Definition zlaqhe.f:134
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 zpoequ(n, a, lda, s, scond, amax, info)
ZPOEQU
Definition zpoequ.f:113
subroutine zposv(uplo, n, nrhs, a, lda, b, ldb, info)
ZPOSV computes the solution to system of linear equations A * X = B for PO matrices
Definition zposv.f:130
subroutine zposvx(fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZPOSVX computes the solution to system of linear equations A * X = B for PO matrices
Definition zposvx.f:306
subroutine zposvxx(fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
ZPOSVXX computes the solution to system of linear equations A * X = B for PO matrices
Definition zposvxx.f:493
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:107
subroutine zpotri(uplo, n, a, lda, info)
ZPOTRI
Definition zpotri.f:95
subroutine zdrvpo(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, asav, b, bsav, x, xact, s, work, rwork, nout)
ZDRVPO
Definition zdrvpo.f:159
subroutine zebchvxx(thresh, path)
ZEBCHVXX
Definition zebchvxx.f:96
subroutine zerrvx(path, nunit)
ZERRVX
Definition zerrvx.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlaipd(n, a, inda, vinda)
ZLAIPD
Definition zlaipd.f:83
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
subroutine zpot01(uplo, n, a, lda, afac, ldafac, rwork, resid)
ZPOT01
Definition zpot01.f:106
subroutine zpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZPOT02
Definition zpot02.f:127
subroutine zpot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPOT05
Definition zpot05.f:165