LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cdrvge.f
Go to the documentation of this file.
1*> \brief \b CDRVGE
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 CDRVGE( 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 RWORK( * ), S( * )
24* COMPLEX A( * ), AFAC( * ), ASAV( * ), B( * ),
25* $ BSAV( * ), WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> CDRVGE tests the driver routines CGESV 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 column 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 COMPLEX array, dimension (NMAX*NMAX)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*> AFAC is COMPLEX array, dimension (NMAX*NMAX)
96*> \endverbatim
97*>
98*> \param[out] ASAV
99*> \verbatim
100*> ASAV is COMPLEX array, dimension (NMAX*NMAX)
101*> \endverbatim
102*>
103*> \param[out] B
104*> \verbatim
105*> B is COMPLEX array, dimension (NMAX*NRHS)
106*> \endverbatim
107*>
108*> \param[out] BSAV
109*> \verbatim
110*> BSAV is COMPLEX array, dimension (NMAX*NRHS)
111*> \endverbatim
112*>
113*> \param[out] X
114*> \verbatim
115*> X is COMPLEX array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] XACT
119*> \verbatim
120*> XACT is COMPLEX array, dimension (NMAX*NRHS)
121*> \endverbatim
122*>
123*> \param[out] S
124*> \verbatim
125*> S is REAL array, dimension (2*NMAX)
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> WORK is COMPLEX array, dimension
131*> (NMAX*max(3,NRHS))
132*> \endverbatim
133*>
134*> \param[out] RWORK
135*> \verbatim
136*> RWORK is REAL array, dimension (2*NRHS+NMAX)
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*> \ingroup complex_lin
159*
160* =====================================================================
161 SUBROUTINE cdrvge( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
162 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
163 $ RWORK, IWORK, NOUT )
164*
165* -- LAPACK test routine --
166* -- LAPACK is a software package provided by Univ. of Tennessee, --
167* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
168*
169* .. Scalar Arguments ..
170 LOGICAL TSTERR
171 INTEGER NMAX, NN, NOUT, NRHS
172 REAL THRESH
173* ..
174* .. Array Arguments ..
175 LOGICAL DOTYPE( * )
176 INTEGER IWORK( * ), NVAL( * )
177 REAL RWORK( * ), S( * )
178 COMPLEX A( * ), AFAC( * ), ASAV( * ), B( * ),
179 $ bsav( * ), work( * ), x( * ), xact( * )
180* ..
181*
182* =====================================================================
183*
184* .. Parameters ..
185 REAL ONE, ZERO
186 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
187 INTEGER NTYPES
188 parameter( ntypes = 11 )
189 INTEGER NTESTS
190 parameter( ntests = 7 )
191 INTEGER NTRAN
192 parameter( ntran = 3 )
193* ..
194* .. Local Scalars ..
195 LOGICAL EQUIL, NOFACT, PREFAC, TRFCON, ZEROT
196 CHARACTER DIST, EQUED, FACT, TRANS, TYPE, XTYPE
197 CHARACTER*3 PATH
198 INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, ITRAN,
199 $ izero, k, k1, kl, ku, lda, lwork, mode, n, nb,
200 $ nbmin, nerrs, nfact, nfail, nimat, nrun, nt
201 REAL AINVNM, AMAX, ANORM, ANORMI, ANORMO, CNDNUM,
202 $ COLCND, RCOND, RCONDC, RCONDI, RCONDO, ROLDC,
203 $ roldi, roldo, rowcnd, rpvgrw
204* ..
205* .. Local Arrays ..
206 CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN )
207 INTEGER ISEED( 4 ), ISEEDY( 4 )
208 REAL RDUM( 1 ), RESULT( NTESTS )
209* ..
210* .. External Functions ..
211 LOGICAL LSAME
212 REAL CLANGE, CLANTR, SGET06, SLAMCH
213 EXTERNAL lsame, clange, clantr, sget06, slamch
214* ..
215* .. External Subroutines ..
216 EXTERNAL aladhd, alaerh, alasvm, cerrvx, cgeequ, cgesv,
219 $ clatms, xlaenv
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC abs, cmplx, max
223* ..
224* .. Scalars in Common ..
225 LOGICAL LERR, OK
226 CHARACTER*32 SRNAMT
227 INTEGER INFOT, NUNIT
228* ..
229* .. Common blocks ..
230 COMMON / infoc / infot, nunit, ok, lerr
231 COMMON / srnamc / srnamt
232* ..
233* .. Data statements ..
234 DATA iseedy / 1988, 1989, 1990, 1991 /
235 DATA transs / 'N', 'T', 'C' /
236 DATA facts / 'F', 'N', 'E' /
237 DATA equeds / 'N', 'R', 'C', 'B' /
238* ..
239* .. Executable Statements ..
240*
241* Initialize constants and the random number seed.
242*
243 path( 1: 1 ) = 'Complex precision'
244 path( 2: 3 ) = 'GE'
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 cerrvx( path, nout )
256 infot = 0
257*
258* Set the block size and minimum block size for testing.
259*
260 nb = 1
261 nbmin = 2
262 CALL xlaenv( 1, nb )
263 CALL xlaenv( 2, nbmin )
264*
265* Do for each value of N in NVAL
266*
267 DO 90 in = 1, nn
268 n = nval( in )
269 lda = max( n, 1 )
270 xtype = 'N'
271 nimat = ntypes
272 IF( n.LE.0 )
273 $ nimat = 1
274*
275 DO 80 imat = 1, nimat
276*
277* Do the tests only if DOTYPE( IMAT ) is true.
278*
279 IF( .NOT.dotype( imat ) )
280 $ GO TO 80
281*
282* Skip types 5, 6, or 7 if the matrix size is too small.
283*
284 zerot = imat.GE.5 .AND. imat.LE.7
285 IF( zerot .AND. n.LT.imat-4 )
286 $ GO TO 80
287*
288* Set up parameters with CLATB4 and generate a test matrix
289* with CLATMS.
290*
291 CALL clatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
292 $ cndnum, dist )
293 rcondc = one / cndnum
294*
295 srnamt = 'CLATMS'
296 CALL clatms( n, n, dist, iseed, TYPE, rwork, mode, cndnum,
297 $ anorm, kl, ku, 'No packing', a, lda, work,
298 $ info )
299*
300* Check error code from CLATMS.
301*
302 IF( info.NE.0 ) THEN
303 CALL alaerh( path, 'CLATMS', info, 0, ' ', n, n, -1, -1,
304 $ -1, imat, nfail, nerrs, nout )
305 GO TO 80
306 END IF
307*
308* For types 5-7, zero one or more columns of the matrix to
309* test that INFO is returned correctly.
310*
311 IF( zerot ) THEN
312 IF( imat.EQ.5 ) THEN
313 izero = 1
314 ELSE IF( imat.EQ.6 ) THEN
315 izero = n
316 ELSE
317 izero = n / 2 + 1
318 END IF
319 ioff = ( izero-1 )*lda
320 IF( imat.LT.7 ) THEN
321 DO 20 i = 1, n
322 a( ioff+i ) = zero
323 20 CONTINUE
324 ELSE
325 CALL claset( 'Full', n, n-izero+1, cmplx( zero ),
326 $ cmplx( zero ), a( ioff+1 ), lda )
327 END IF
328 ELSE
329 izero = 0
330 END IF
331*
332* Save a copy of the matrix A in ASAV.
333*
334 CALL clacpy( 'Full', n, n, a, lda, asav, lda )
335*
336 DO 70 iequed = 1, 4
337 equed = equeds( iequed )
338 IF( iequed.EQ.1 ) THEN
339 nfact = 3
340 ELSE
341 nfact = 1
342 END IF
343*
344 DO 60 ifact = 1, nfact
345 fact = facts( ifact )
346 prefac = lsame( fact, 'F' )
347 nofact = lsame( fact, 'N' )
348 equil = lsame( fact, 'E' )
349*
350 IF( zerot ) THEN
351 IF( prefac )
352 $ GO TO 60
353 rcondo = zero
354 rcondi = zero
355*
356 ELSE IF( .NOT.nofact ) THEN
357*
358* Compute the condition number for comparison with
359* the value returned by CGESVX (FACT = 'N' reuses
360* the condition number from the previous iteration
361* with FACT = 'F').
362*
363 CALL clacpy( 'Full', n, n, asav, lda, afac, lda )
364 IF( equil .OR. iequed.GT.1 ) THEN
365*
366* Compute row and column scale factors to
367* equilibrate the matrix A.
368*
369 CALL cgeequ( n, n, afac, lda, s, s( n+1 ),
370 $ rowcnd, colcnd, amax, info )
371 IF( info.EQ.0 .AND. n.GT.0 ) THEN
372 IF( lsame( equed, 'R' ) ) THEN
373 rowcnd = zero
374 colcnd = one
375 ELSE IF( lsame( equed, 'C' ) ) THEN
376 rowcnd = one
377 colcnd = zero
378 ELSE IF( lsame( equed, 'B' ) ) THEN
379 rowcnd = zero
380 colcnd = zero
381 END IF
382*
383* Equilibrate the matrix.
384*
385 CALL claqge( n, n, afac, lda, s, s( n+1 ),
386 $ rowcnd, colcnd, amax, equed )
387 END IF
388 END IF
389*
390* Save the condition number of the non-equilibrated
391* system for use in CGET04.
392*
393 IF( equil ) THEN
394 roldo = rcondo
395 roldi = rcondi
396 END IF
397*
398* Compute the 1-norm and infinity-norm of A.
399*
400 anormo = clange( '1', n, n, afac, lda, rwork )
401 anormi = clange( 'I', n, n, afac, lda, rwork )
402*
403* Factor the matrix A.
404*
405 srnamt = 'CGETRF'
406 CALL cgetrf( n, n, afac, lda, iwork, info )
407*
408* Form the inverse of A.
409*
410 CALL clacpy( 'Full', n, n, afac, lda, a, lda )
411 lwork = nmax*max( 3, nrhs )
412 srnamt = 'CGETRI'
413 CALL cgetri( n, a, lda, iwork, work, lwork, info )
414*
415* Compute the 1-norm condition number of A.
416*
417 ainvnm = clange( '1', n, n, a, lda, rwork )
418 IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
419 rcondo = one
420 ELSE
421 rcondo = ( one / anormo ) / ainvnm
422 END IF
423*
424* Compute the infinity-norm condition number of A.
425*
426 ainvnm = clange( 'I', n, n, a, lda, rwork )
427 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
428 rcondi = one
429 ELSE
430 rcondi = ( one / anormi ) / ainvnm
431 END IF
432 END IF
433*
434 DO 50 itran = 1, ntran
435*
436* Do for each value of TRANS.
437*
438 trans = transs( itran )
439 IF( itran.EQ.1 ) THEN
440 rcondc = rcondo
441 ELSE
442 rcondc = rcondi
443 END IF
444*
445* Restore the matrix A.
446*
447 CALL clacpy( 'Full', n, n, asav, lda, a, lda )
448*
449* Form an exact solution and set the right hand side.
450*
451 srnamt = 'CLARHS'
452 CALL clarhs( path, xtype, 'Full', trans, n, n, kl,
453 $ ku, nrhs, a, lda, xact, lda, b, lda,
454 $ iseed, info )
455 xtype = 'C'
456 CALL clacpy( 'Full', n, nrhs, b, lda, bsav, lda )
457*
458 IF( nofact .AND. itran.EQ.1 ) THEN
459*
460* --- Test CGESV ---
461*
462* Compute the LU factorization of the matrix and
463* solve the system.
464*
465 CALL clacpy( 'Full', n, n, a, lda, afac, lda )
466 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
467*
468 srnamt = 'CGESV '
469 CALL cgesv( n, nrhs, afac, lda, iwork, x, lda,
470 $ info )
471*
472* Check error code from CGESV .
473*
474 IF( info.NE.izero )
475 $ CALL alaerh( path, 'CGESV ', info, izero,
476 $ ' ', n, n, -1, -1, nrhs, imat,
477 $ nfail, nerrs, nout )
478*
479* Reconstruct matrix from factors and compute
480* residual.
481*
482 CALL cget01( n, n, a, lda, afac, lda, iwork,
483 $ rwork, result( 1 ) )
484 nt = 1
485 IF( izero.EQ.0 ) THEN
486*
487* Compute residual of the computed solution.
488*
489 CALL clacpy( 'Full', n, nrhs, b, lda, work,
490 $ lda )
491 CALL cget02( 'No transpose', n, n, nrhs, a,
492 $ lda, x, lda, work, lda, rwork,
493 $ result( 2 ) )
494*
495* Check solution from generated exact solution.
496*
497 CALL cget04( n, nrhs, x, lda, xact, lda,
498 $ rcondc, result( 3 ) )
499 nt = 3
500 END IF
501*
502* Print information about the tests that did not
503* pass the threshold.
504*
505 DO 30 k = 1, nt
506 IF( result( k ).GE.thresh ) THEN
507 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
508 $ CALL aladhd( nout, path )
509 WRITE( nout, fmt = 9999 )'CGESV ', n,
510 $ imat, k, result( k )
511 nfail = nfail + 1
512 END IF
513 30 CONTINUE
514 nrun = nrun + nt
515 END IF
516*
517* --- Test CGESVX ---
518*
519 IF( .NOT.prefac )
520 $ CALL claset( 'Full', n, n, cmplx( zero ),
521 $ cmplx( zero ), afac, lda )
522 CALL claset( 'Full', n, nrhs, cmplx( zero ),
523 $ cmplx( zero ), x, lda )
524 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
525*
526* Equilibrate the matrix if FACT = 'F' and
527* EQUED = 'R', 'C', or 'B'.
528*
529 CALL claqge( n, n, a, lda, s, s( n+1 ), rowcnd,
530 $ colcnd, amax, equed )
531 END IF
532*
533* Solve the system and compute the condition number
534* and error bounds using CGESVX.
535*
536 srnamt = 'CGESVX'
537 CALL cgesvx( fact, trans, n, nrhs, a, lda, afac,
538 $ lda, iwork, equed, s, s( n+1 ), b,
539 $ lda, x, lda, rcond, rwork,
540 $ rwork( nrhs+1 ), work,
541 $ rwork( 2*nrhs+1 ), info )
542*
543* Check the error code from CGESVX.
544*
545 IF( info.NE.izero )
546 $ CALL alaerh( path, 'CGESVX', info, izero,
547 $ fact // trans, n, n, -1, -1, nrhs,
548 $ imat, nfail, nerrs, nout )
549*
550* Compare RWORK(2*NRHS+1) from CGESVX with the
551* computed reciprocal pivot growth factor RPVGRW
552*
553 IF( info.NE.0 .AND. info.LE.n) THEN
554 rpvgrw = clantr( 'M', 'U', 'N', info, info,
555 $ afac, lda, rdum )
556 IF( rpvgrw.EQ.zero ) THEN
557 rpvgrw = one
558 ELSE
559 rpvgrw = clange( 'M', n, info, a, lda,
560 $ rdum ) / rpvgrw
561 END IF
562 ELSE
563 rpvgrw = clantr( 'M', 'U', 'N', n, n, afac, lda,
564 $ rdum )
565 IF( rpvgrw.EQ.zero ) THEN
566 rpvgrw = one
567 ELSE
568 rpvgrw = clange( 'M', n, n, a, lda, rdum ) /
569 $ rpvgrw
570 END IF
571 END IF
572 result( 7 ) = abs( rpvgrw-rwork( 2*nrhs+1 ) ) /
573 $ max( rwork( 2*nrhs+1 ), rpvgrw ) /
574 $ slamch( 'E' )
575*
576 IF( .NOT.prefac ) THEN
577*
578* Reconstruct matrix from factors and compute
579* residual.
580*
581 CALL cget01( n, n, a, lda, afac, lda, iwork,
582 $ rwork( 2*nrhs+1 ), result( 1 ) )
583 k1 = 1
584 ELSE
585 k1 = 2
586 END IF
587*
588 IF( info.EQ.0 ) THEN
589 trfcon = .false.
590*
591* Compute residual of the computed solution.
592*
593 CALL clacpy( 'Full', n, nrhs, bsav, lda, work,
594 $ lda )
595 CALL cget02( trans, n, n, nrhs, asav, lda, x,
596 $ lda, work, lda, rwork( 2*nrhs+1 ),
597 $ result( 2 ) )
598*
599* Check solution from generated exact solution.
600*
601 IF( nofact .OR. ( prefac .AND. lsame( equed,
602 $ 'N' ) ) ) THEN
603 CALL cget04( n, nrhs, x, lda, xact, lda,
604 $ rcondc, result( 3 ) )
605 ELSE
606 IF( itran.EQ.1 ) THEN
607 roldc = roldo
608 ELSE
609 roldc = roldi
610 END IF
611 CALL cget04( n, nrhs, x, lda, xact, lda,
612 $ roldc, result( 3 ) )
613 END IF
614*
615* Check the error bounds from iterative
616* refinement.
617*
618 CALL cget07( trans, n, nrhs, asav, lda, b, lda,
619 $ x, lda, xact, lda, rwork, .true.,
620 $ rwork( nrhs+1 ), result( 4 ) )
621 ELSE
622 trfcon = .true.
623 END IF
624*
625* Compare RCOND from CGESVX with the computed value
626* in RCONDC.
627*
628 result( 6 ) = sget06( rcond, rcondc )
629*
630* Print information about the tests that did not pass
631* the threshold.
632*
633 IF( .NOT.trfcon ) THEN
634 DO 40 k = k1, ntests
635 IF( result( k ).GE.thresh ) THEN
636 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
637 $ CALL aladhd( nout, path )
638 IF( prefac ) THEN
639 WRITE( nout, fmt = 9997 )'CGESVX',
640 $ fact, trans, n, equed, imat, k,
641 $ result( k )
642 ELSE
643 WRITE( nout, fmt = 9998 )'CGESVX',
644 $ fact, trans, n, imat, k, result( k )
645 END IF
646 nfail = nfail + 1
647 END IF
648 40 CONTINUE
649 nrun = nrun + ntests - k1 + 1
650 ELSE
651 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
652 $ THEN
653 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
654 $ CALL aladhd( nout, path )
655 IF( prefac ) THEN
656 WRITE( nout, fmt = 9997 )'CGESVX', fact,
657 $ trans, n, equed, imat, 1, result( 1 )
658 ELSE
659 WRITE( nout, fmt = 9998 )'CGESVX', fact,
660 $ trans, n, imat, 1, result( 1 )
661 END IF
662 nfail = nfail + 1
663 nrun = nrun + 1
664 END IF
665 IF( result( 6 ).GE.thresh ) THEN
666 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
667 $ CALL aladhd( nout, path )
668 IF( prefac ) THEN
669 WRITE( nout, fmt = 9997 )'CGESVX', fact,
670 $ trans, n, equed, imat, 6, result( 6 )
671 ELSE
672 WRITE( nout, fmt = 9998 )'CGESVX', fact,
673 $ trans, n, imat, 6, result( 6 )
674 END IF
675 nfail = nfail + 1
676 nrun = nrun + 1
677 END IF
678 IF( result( 7 ).GE.thresh ) THEN
679 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
680 $ CALL aladhd( nout, path )
681 IF( prefac ) THEN
682 WRITE( nout, fmt = 9997 )'CGESVX', fact,
683 $ trans, n, equed, imat, 7, result( 7 )
684 ELSE
685 WRITE( nout, fmt = 9998 )'CGESVX', fact,
686 $ trans, n, imat, 7, result( 7 )
687 END IF
688 nfail = nfail + 1
689 nrun = nrun + 1
690 END IF
691*
692 END IF
693*
694 50 CONTINUE
695 60 CONTINUE
696 70 CONTINUE
697 80 CONTINUE
698 90 CONTINUE
699*
700* Print a summary of the results.
701*
702 CALL alasvm( path, nout, nfail, nrun, nerrs )
703*
704 9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test(', i2, ') =',
705 $ g12.5 )
706 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
707 $ ', type ', i2, ', test(', i1, ')=', g12.5 )
708 9997 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
709 $ ', EQUED=''', a1, ''', type ', i2, ', test(', i1, ')=',
710 $ g12.5 )
711 RETURN
712*
713* End of CDRVGE
714*
715 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine cget02(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
CGET02
Definition cget02.f:134
subroutine clarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
CLARHS
Definition clarhs.f:208
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
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 cdrvge(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, asav, b, bsav, x, xact, s, work, rwork, iwork, nout)
CDRVGE
Definition cdrvge.f:164
subroutine cerrvx(path, nunit)
CERRVX
Definition cerrvx.f:55
subroutine cget01(m, n, a, lda, afac, ldafac, ipiv, rwork, resid)
CGET01
Definition cget01.f:108
subroutine cget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
CGET04
Definition cget04.f:102
subroutine cget07(trans, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, chkferr, berr, reslts)
CGET07
Definition cget07.f:166
subroutine clatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
CLATB4
Definition clatb4.f:121
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine cgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
CGEEQU
Definition cgeequ.f:140
subroutine cgesv(n, nrhs, a, lda, ipiv, b, ldb, info)
Download CGESV + dependencies <a href="http://www.netlib.org/cgi-bin/netlibfiles....
Definition cgesv.f:124
subroutine cgesvx(fact, trans, n, nrhs, a, lda, af, ldaf, ipiv, equed, r, c, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
CGESVX computes the solution to system of linear equations A * X = B for GE matrices
Definition cgesvx.f:350
subroutine cgetrf(m, n, a, lda, ipiv, info)
CGETRF
Definition cgetrf.f:108
subroutine cgetri(n, a, lda, ipiv, work, lwork, info)
CGETRI
Definition cgetri.f:114
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
subroutine claqge(m, n, a, lda, r, c, rowcnd, colcnd, amax, equed)
CLAQGE scales a general rectangular matrix, using row and column scaling factors computed by sgeequ.
Definition claqge.f:143
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:106