LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cdrvge ( 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, dimension( * )  IWORK,
integer  NOUT 
)

CDRVGE

CDRVGEX

Purpose:
 CDRVGE tests the driver routines CGESV 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 column 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 (2*NMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is REAL array, dimension (2*NRHS+NMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[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 2015
Purpose:
 CDRVGE tests the driver routines CGESV, -SVX, and -SVXX.

 Note that this file is used only when the XBLAS are available,
 otherwise cdrvge.f defines this subroutine.
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 column 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 (2*NMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is REAL array, dimension (2*NRHS+NMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[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
April 2012

Definition at line 166 of file cdrvge.f.

166 *
167 * -- LAPACK test routine (version 3.6.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 2015
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 rwork( * ), s( * )
181  COMPLEX a( * ), afac( * ), asav( * ), b( * ),
182  $ bsav( * ), work( * ), 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
191  parameter ( ntypes = 11 )
192  INTEGER ntests
193  parameter ( ntests = 7 )
194  INTEGER ntran
195  parameter ( ntran = 3 )
196 * ..
197 * .. Local Scalars ..
198  LOGICAL equil, nofact, prefac, trfcon, zerot
199  CHARACTER dist, equed, fact, trans, TYPE, xtype
200  CHARACTER*3 path
201  INTEGER i, iequed, ifact, imat, in, info, ioff, itran,
202  $ izero, k, k1, kl, ku, lda, lwork, mode, n, nb,
203  $ nbmin, nerrs, nfact, nfail, nimat, nrun, nt
204  REAL ainvnm, amax, anorm, anormi, anormo, cndnum,
205  $ colcnd, rcond, rcondc, rcondi, rcondo, roldc,
206  $ roldi, roldo, rowcnd, rpvgrw
207 * ..
208 * .. Local Arrays ..
209  CHARACTER equeds( 4 ), facts( 3 ), transs( ntran )
210  INTEGER iseed( 4 ), iseedy( 4 )
211  REAL rdum( 1 ), result( ntests )
212 * ..
213 * .. External Functions ..
214  LOGICAL lsame
215  REAL clange, clantr, sget06, slamch
216  EXTERNAL lsame, clange, clantr, sget06, slamch
217 * ..
218 * .. External Subroutines ..
219  EXTERNAL aladhd, alaerh, alasvm, cerrvx, cgeequ, cgesv,
222  $ clatms, xlaenv
223 * ..
224 * .. Intrinsic Functions ..
225  INTRINSIC abs, cmplx, max
226 * ..
227 * .. Scalars in Common ..
228  LOGICAL lerr, ok
229  CHARACTER*32 srnamt
230  INTEGER infot, nunit
231 * ..
232 * .. Common blocks ..
233  COMMON / infoc / infot, nunit, ok, lerr
234  COMMON / srnamc / srnamt
235 * ..
236 * .. Data statements ..
237  DATA iseedy / 1988, 1989, 1990, 1991 /
238  DATA transs / 'N', 'T', 'C' /
239  DATA facts / 'F', 'N', 'E' /
240  DATA equeds / 'N', 'R', 'C', 'B' /
241 * ..
242 * .. Executable Statements ..
243 *
244 * Initialize constants and the random number seed.
245 *
246  path( 1: 1 ) = 'Complex precision'
247  path( 2: 3 ) = 'GE'
248  nrun = 0
249  nfail = 0
250  nerrs = 0
251  DO 10 i = 1, 4
252  iseed( i ) = iseedy( i )
253  10 CONTINUE
254 *
255 * Test the error exits
256 *
257  IF( tsterr )
258  $ CALL cerrvx( path, nout )
259  infot = 0
260 *
261 * Set the block size and minimum block size for testing.
262 *
263  nb = 1
264  nbmin = 2
265  CALL xlaenv( 1, nb )
266  CALL xlaenv( 2, nbmin )
267 *
268 * Do for each value of N in NVAL
269 *
270  DO 90 in = 1, nn
271  n = nval( in )
272  lda = max( n, 1 )
273  xtype = 'N'
274  nimat = ntypes
275  IF( n.LE.0 )
276  $ nimat = 1
277 *
278  DO 80 imat = 1, nimat
279 *
280 * Do the tests only if DOTYPE( IMAT ) is true.
281 *
282  IF( .NOT.dotype( imat ) )
283  $ GO TO 80
284 *
285 * Skip types 5, 6, or 7 if the matrix size is too small.
286 *
287  zerot = imat.GE.5 .AND. imat.LE.7
288  IF( zerot .AND. n.LT.imat-4 )
289  $ GO TO 80
290 *
291 * Set up parameters with CLATB4 and generate a test matrix
292 * with CLATMS.
293 *
294  CALL clatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
295  $ cndnum, dist )
296  rcondc = one / cndnum
297 *
298  srnamt = 'CLATMS'
299  CALL clatms( n, n, dist, iseed, TYPE, rwork, mode, cndnum,
300  $ anorm, kl, ku, 'No packing', a, lda, work,
301  $ info )
302 *
303 * Check error code from CLATMS.
304 *
305  IF( info.NE.0 ) THEN
306  CALL alaerh( path, 'CLATMS', info, 0, ' ', n, n, -1, -1,
307  $ -1, imat, nfail, nerrs, nout )
308  GO TO 80
309  END IF
310 *
311 * For types 5-7, zero one or more columns of the matrix to
312 * test that INFO is returned correctly.
313 *
314  IF( zerot ) THEN
315  IF( imat.EQ.5 ) THEN
316  izero = 1
317  ELSE IF( imat.EQ.6 ) THEN
318  izero = n
319  ELSE
320  izero = n / 2 + 1
321  END IF
322  ioff = ( izero-1 )*lda
323  IF( imat.LT.7 ) THEN
324  DO 20 i = 1, n
325  a( ioff+i ) = zero
326  20 CONTINUE
327  ELSE
328  CALL claset( 'Full', n, n-izero+1, cmplx( zero ),
329  $ cmplx( zero ), a( ioff+1 ), lda )
330  END IF
331  ELSE
332  izero = 0
333  END IF
334 *
335 * Save a copy of the matrix A in ASAV.
336 *
337  CALL clacpy( 'Full', n, n, a, lda, asav, lda )
338 *
339  DO 70 iequed = 1, 4
340  equed = equeds( iequed )
341  IF( iequed.EQ.1 ) THEN
342  nfact = 3
343  ELSE
344  nfact = 1
345  END IF
346 *
347  DO 60 ifact = 1, nfact
348  fact = facts( ifact )
349  prefac = lsame( fact, 'F' )
350  nofact = lsame( fact, 'N' )
351  equil = lsame( fact, 'E' )
352 *
353  IF( zerot ) THEN
354  IF( prefac )
355  $ GO TO 60
356  rcondo = zero
357  rcondi = zero
358 *
359  ELSE IF( .NOT.nofact ) THEN
360 *
361 * Compute the condition number for comparison with
362 * the value returned by CGESVX (FACT = 'N' reuses
363 * the condition number from the previous iteration
364 * with FACT = 'F').
365 *
366  CALL clacpy( 'Full', n, n, asav, lda, afac, lda )
367  IF( equil .OR. iequed.GT.1 ) THEN
368 *
369 * Compute row and column scale factors to
370 * equilibrate the matrix A.
371 *
372  CALL cgeequ( n, n, afac, lda, s, s( n+1 ),
373  $ rowcnd, colcnd, amax, info )
374  IF( info.EQ.0 .AND. n.GT.0 ) THEN
375  IF( lsame( equed, 'R' ) ) THEN
376  rowcnd = zero
377  colcnd = one
378  ELSE IF( lsame( equed, 'C' ) ) THEN
379  rowcnd = one
380  colcnd = zero
381  ELSE IF( lsame( equed, 'B' ) ) THEN
382  rowcnd = zero
383  colcnd = zero
384  END IF
385 *
386 * Equilibrate the matrix.
387 *
388  CALL claqge( n, n, afac, lda, s, s( n+1 ),
389  $ rowcnd, colcnd, amax, equed )
390  END IF
391  END IF
392 *
393 * Save the condition number of the non-equilibrated
394 * system for use in CGET04.
395 *
396  IF( equil ) THEN
397  roldo = rcondo
398  roldi = rcondi
399  END IF
400 *
401 * Compute the 1-norm and infinity-norm of A.
402 *
403  anormo = clange( '1', n, n, afac, lda, rwork )
404  anormi = clange( 'I', n, n, afac, lda, rwork )
405 *
406 * Factor the matrix A.
407 *
408  srnamt = 'CGETRF'
409  CALL cgetrf( n, n, afac, lda, iwork, info )
410 *
411 * Form the inverse of A.
412 *
413  CALL clacpy( 'Full', n, n, afac, lda, a, lda )
414  lwork = nmax*max( 3, nrhs )
415  srnamt = 'CGETRI'
416  CALL cgetri( n, a, lda, iwork, work, lwork, info )
417 *
418 * Compute the 1-norm condition number of A.
419 *
420  ainvnm = clange( '1', n, n, a, lda, rwork )
421  IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
422  rcondo = one
423  ELSE
424  rcondo = ( one / anormo ) / ainvnm
425  END IF
426 *
427 * Compute the infinity-norm condition number of A.
428 *
429  ainvnm = clange( 'I', n, n, a, lda, rwork )
430  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
431  rcondi = one
432  ELSE
433  rcondi = ( one / anormi ) / ainvnm
434  END IF
435  END IF
436 *
437  DO 50 itran = 1, ntran
438 *
439 * Do for each value of TRANS.
440 *
441  trans = transs( itran )
442  IF( itran.EQ.1 ) THEN
443  rcondc = rcondo
444  ELSE
445  rcondc = rcondi
446  END IF
447 *
448 * Restore the matrix A.
449 *
450  CALL clacpy( 'Full', n, n, asav, lda, a, lda )
451 *
452 * Form an exact solution and set the right hand side.
453 *
454  srnamt = 'CLARHS'
455  CALL clarhs( path, xtype, 'Full', trans, n, n, kl,
456  $ ku, nrhs, a, lda, xact, lda, b, lda,
457  $ iseed, info )
458  xtype = 'C'
459  CALL clacpy( 'Full', n, nrhs, b, lda, bsav, lda )
460 *
461  IF( nofact .AND. itran.EQ.1 ) THEN
462 *
463 * --- Test CGESV ---
464 *
465 * Compute the LU factorization of the matrix and
466 * solve the system.
467 *
468  CALL clacpy( 'Full', n, n, a, lda, afac, lda )
469  CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
470 *
471  srnamt = 'CGESV '
472  CALL cgesv( n, nrhs, afac, lda, iwork, x, lda,
473  $ info )
474 *
475 * Check error code from CGESV .
476 *
477  IF( info.NE.izero )
478  $ CALL alaerh( path, 'CGESV ', info, izero,
479  $ ' ', n, n, -1, -1, nrhs, imat,
480  $ nfail, nerrs, nout )
481 *
482 * Reconstruct matrix from factors and compute
483 * residual.
484 *
485  CALL cget01( n, n, a, lda, afac, lda, iwork,
486  $ rwork, result( 1 ) )
487  nt = 1
488  IF( izero.EQ.0 ) THEN
489 *
490 * Compute residual of the computed solution.
491 *
492  CALL clacpy( 'Full', n, nrhs, b, lda, work,
493  $ lda )
494  CALL cget02( 'No transpose', n, n, nrhs, a,
495  $ lda, x, lda, work, lda, rwork,
496  $ result( 2 ) )
497 *
498 * Check solution from generated exact solution.
499 *
500  CALL cget04( n, nrhs, x, lda, xact, lda,
501  $ rcondc, result( 3 ) )
502  nt = 3
503  END IF
504 *
505 * Print information about the tests that did not
506 * pass the threshold.
507 *
508  DO 30 k = 1, nt
509  IF( result( k ).GE.thresh ) THEN
510  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
511  $ CALL aladhd( nout, path )
512  WRITE( nout, fmt = 9999 )'CGESV ', n,
513  $ imat, k, result( k )
514  nfail = nfail + 1
515  END IF
516  30 CONTINUE
517  nrun = nrun + nt
518  END IF
519 *
520 * --- Test CGESVX ---
521 *
522  IF( .NOT.prefac )
523  $ CALL claset( 'Full', n, n, cmplx( zero ),
524  $ cmplx( zero ), afac, lda )
525  CALL claset( 'Full', n, nrhs, cmplx( zero ),
526  $ cmplx( zero ), x, lda )
527  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
528 *
529 * Equilibrate the matrix if FACT = 'F' and
530 * EQUED = 'R', 'C', or 'B'.
531 *
532  CALL claqge( n, n, a, lda, s, s( n+1 ), rowcnd,
533  $ colcnd, amax, equed )
534  END IF
535 *
536 * Solve the system and compute the condition number
537 * and error bounds using CGESVX.
538 *
539  srnamt = 'CGESVX'
540  CALL cgesvx( fact, trans, n, nrhs, a, lda, afac,
541  $ lda, iwork, equed, s, s( n+1 ), b,
542  $ lda, x, lda, rcond, rwork,
543  $ rwork( nrhs+1 ), work,
544  $ rwork( 2*nrhs+1 ), info )
545 *
546 * Check the error code from CGESVX.
547 *
548  IF( info.NE.izero )
549  $ CALL alaerh( path, 'CGESVX', info, izero,
550  $ fact // trans, n, n, -1, -1, nrhs,
551  $ imat, nfail, nerrs, nout )
552 *
553 * Compare RWORK(2*NRHS+1) from CGESVX with the
554 * computed reciprocal pivot growth factor RPVGRW
555 *
556  IF( info.NE.0 .AND. info.LE.n) THEN
557  rpvgrw = clantr( 'M', 'U', 'N', info, info,
558  $ afac, lda, rdum )
559  IF( rpvgrw.EQ.zero ) THEN
560  rpvgrw = one
561  ELSE
562  rpvgrw = clange( 'M', n, info, a, lda,
563  $ rdum ) / rpvgrw
564  END IF
565  ELSE
566  rpvgrw = clantr( 'M', 'U', 'N', n, n, afac, lda,
567  $ rdum )
568  IF( rpvgrw.EQ.zero ) THEN
569  rpvgrw = one
570  ELSE
571  rpvgrw = clange( 'M', n, n, a, lda, rdum ) /
572  $ rpvgrw
573  END IF
574  END IF
575  result( 7 ) = abs( rpvgrw-rwork( 2*nrhs+1 ) ) /
576  $ max( rwork( 2*nrhs+1 ), rpvgrw ) /
577  $ slamch( 'E' )
578 *
579  IF( .NOT.prefac ) THEN
580 *
581 * Reconstruct matrix from factors and compute
582 * residual.
583 *
584  CALL cget01( n, n, a, lda, afac, lda, iwork,
585  $ rwork( 2*nrhs+1 ), result( 1 ) )
586  k1 = 1
587  ELSE
588  k1 = 2
589  END IF
590 *
591  IF( info.EQ.0 ) THEN
592  trfcon = .false.
593 *
594 * Compute residual of the computed solution.
595 *
596  CALL clacpy( 'Full', n, nrhs, bsav, lda, work,
597  $ lda )
598  CALL cget02( trans, n, n, nrhs, asav, lda, x,
599  $ lda, work, lda, rwork( 2*nrhs+1 ),
600  $ result( 2 ) )
601 *
602 * Check solution from generated exact solution.
603 *
604  IF( nofact .OR. ( prefac .AND. lsame( equed,
605  $ 'N' ) ) ) THEN
606  CALL cget04( n, nrhs, x, lda, xact, lda,
607  $ rcondc, result( 3 ) )
608  ELSE
609  IF( itran.EQ.1 ) THEN
610  roldc = roldo
611  ELSE
612  roldc = roldi
613  END IF
614  CALL cget04( n, nrhs, x, lda, xact, lda,
615  $ roldc, result( 3 ) )
616  END IF
617 *
618 * Check the error bounds from iterative
619 * refinement.
620 *
621  CALL cget07( trans, n, nrhs, asav, lda, b, lda,
622  $ x, lda, xact, lda, rwork, .true.,
623  $ rwork( nrhs+1 ), result( 4 ) )
624  ELSE
625  trfcon = .true.
626  END IF
627 *
628 * Compare RCOND from CGESVX with the computed value
629 * in RCONDC.
630 *
631  result( 6 ) = sget06( rcond, rcondc )
632 *
633 * Print information about the tests that did not pass
634 * the threshold.
635 *
636  IF( .NOT.trfcon ) THEN
637  DO 40 k = k1, ntests
638  IF( result( k ).GE.thresh ) THEN
639  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
640  $ CALL aladhd( nout, path )
641  IF( prefac ) THEN
642  WRITE( nout, fmt = 9997 )'CGESVX',
643  $ fact, trans, n, equed, imat, k,
644  $ result( k )
645  ELSE
646  WRITE( nout, fmt = 9998 )'CGESVX',
647  $ fact, trans, n, imat, k, result( k )
648  END IF
649  nfail = nfail + 1
650  END IF
651  40 CONTINUE
652  nrun = nrun + ntests - k1 + 1
653  ELSE
654  IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
655  $ THEN
656  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
657  $ CALL aladhd( nout, path )
658  IF( prefac ) THEN
659  WRITE( nout, fmt = 9997 )'CGESVX', fact,
660  $ trans, n, equed, imat, 1, result( 1 )
661  ELSE
662  WRITE( nout, fmt = 9998 )'CGESVX', fact,
663  $ trans, n, imat, 1, result( 1 )
664  END IF
665  nfail = nfail + 1
666  nrun = nrun + 1
667  END IF
668  IF( result( 6 ).GE.thresh ) THEN
669  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
670  $ CALL aladhd( nout, path )
671  IF( prefac ) THEN
672  WRITE( nout, fmt = 9997 )'CGESVX', fact,
673  $ trans, n, equed, imat, 6, result( 6 )
674  ELSE
675  WRITE( nout, fmt = 9998 )'CGESVX', fact,
676  $ trans, n, imat, 6, result( 6 )
677  END IF
678  nfail = nfail + 1
679  nrun = nrun + 1
680  END IF
681  IF( result( 7 ).GE.thresh ) THEN
682  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
683  $ CALL aladhd( nout, path )
684  IF( prefac ) THEN
685  WRITE( nout, fmt = 9997 )'CGESVX', fact,
686  $ trans, n, equed, imat, 7, result( 7 )
687  ELSE
688  WRITE( nout, fmt = 9998 )'CGESVX', fact,
689  $ trans, n, imat, 7, result( 7 )
690  END IF
691  nfail = nfail + 1
692  nrun = nrun + 1
693  END IF
694 *
695  END IF
696 *
697  50 CONTINUE
698  60 CONTINUE
699  70 CONTINUE
700  80 CONTINUE
701  90 CONTINUE
702 *
703 * Print a summary of the results.
704 *
705  CALL alasvm( path, nout, nfail, nrun, nerrs )
706 *
707  9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test(', i2, ') =',
708  $ g12.5 )
709  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
710  $ ', type ', i2, ', test(', i1, ')=', g12.5 )
711  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
712  $ ', EQUED=''', a1, ''', type ', i2, ', test(', i1, ')=',
713  $ g12.5 )
714  RETURN
715 *
716 * End of CDRVGE
717 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine cget01(M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK, RESID)
CGET01
Definition: cget01.f:110
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 cget07(TRANS, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, CHKFERR, BERR, RESLTS)
CGET07
Definition: cget07.f:168
subroutine cgetri(N, A, LDA, IPIV, WORK, LWORK, INFO)
CGETRI
Definition: cgetri.f:116
subroutine cerrvx(PATH, NUNIT)
CERRVX
Definition: cerrvx.f:57
real function clantr(NORM, UPLO, DIAG, M, N, A, LDA, WORK)
CLANTR returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a trapezoidal or triangular matrix.
Definition: clantr.f:144
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 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:145
subroutine cget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
CGET02
Definition: cget02.f:135
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 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:352
subroutine cgeequ(M, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFO)
CGEEQU
Definition: cgeequ.f:142
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 cgetrf(M, N, A, LDA, IPIV, INFO)
CGETRF
Definition: cgetrf.f:110
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
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine clatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
CLATB4
Definition: clatb4.f:123
subroutine cgesv(N, NRHS, A, LDA, IPIV, B, LDB, INFO)
CGESV computes the solution to system of linear equations A * X = B for GE matrices (simple driver) ...
Definition: cgesv.f:124
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

Here is the call graph for this function:

Here is the caller graph for this function: