161 SUBROUTINE ddrvge( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
162 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
163 $ RWORK, IWORK, NOUT )
171 INTEGER NMAX, NN, NOUT, NRHS
172 DOUBLE PRECISION THRESH
176 INTEGER IWORK( * ), NVAL( * )
177 DOUBLE PRECISION A( * ), AFAC( * ), ASAV( * ), B( * ),
178 $ bsav( * ), rwork( * ), s( * ), work( * ),
185 DOUBLE PRECISION ONE, ZERO
186 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
188 parameter( ntypes = 11 )
190 parameter( ntests = 7 )
192 parameter( ntran = 3 )
195 LOGICAL EQUIL, NOFACT, PREFAC, TRFCON, ZEROT
196 CHARACTER DIST, EQUED, FACT, TRANS,
TYPE, XTYPE
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 DOUBLE PRECISION AINVNM, AMAX, ANORM, ANORMI, ANORMO, CNDNUM,
202 $ COLCND, RCOND, RCONDC, RCONDI, RCONDO, ROLDC,
203 $ roldi, roldo, rowcnd, rpvgrw
206 CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN )
207 INTEGER ISEED( 4 ), ISEEDY( 4 )
208 DOUBLE PRECISION RESULT( NTESTS )
212 DOUBLE PRECISION DGET06, DLAMCH, DLANGE, DLANTR
213 EXTERNAL lsame, dget06, dlamch, dlange, dlantr
230 COMMON / infoc / infot, nunit, ok, lerr
231 COMMON / srnamc / srnamt
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' /
243 path( 1: 1 ) =
'Double precision'
249 iseed( i ) = iseedy( i )
255 $
CALL derrvx( path, nout )
275 DO 80 imat = 1, nimat
279 IF( .NOT.dotype( imat ) )
284 zerot = imat.GE.5 .AND. imat.LE.7
285 IF( zerot .AND. n.LT.imat-4 )
291 CALL dlatb4( path, imat, n, n,
TYPE, kl, ku, anorm, mode,
293 rcondc = one / cndnum
296 CALL dlatms( n, n, dist, iseed,
TYPE, rwork, mode, cndnum,
297 $ anorm, kl, ku,
'No packing', a, lda, work,
303 CALL alaerh( path,
'DLATMS', info, 0,
' ', n, n, -1, -1,
304 $ -1, imat, nfail, nerrs, nout )
314 ELSE IF( imat.EQ.6 )
THEN
319 ioff = ( izero-1 )*lda
325 CALL dlaset(
'Full', n, n-izero+1, zero, zero,
334 CALL dlacpy(
'Full', n, n, a, lda, asav, lda )
337 equed = equeds( iequed )
338 IF( iequed.EQ.1 )
THEN
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' )
356 ELSE IF( .NOT.nofact )
THEN
363 CALL dlacpy(
'Full', n, n, asav, lda, afac, lda )
364 IF( equil .OR. iequed.GT.1 )
THEN
369 CALL dgeequ( 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
375 ELSE IF( lsame( equed,
'C' ) )
THEN
378 ELSE IF( lsame( equed,
'B' ) )
THEN
385 CALL dlaqge( n, n, afac, lda, s, s( n+1 ),
386 $ rowcnd, colcnd, amax, equed )
400 anormo = dlange(
'1', n, n, afac, lda, rwork )
401 anormi = dlange(
'I', n, n, afac, lda, rwork )
406 CALL dgetrf( n, n, afac, lda, iwork, info )
410 CALL dlacpy(
'Full', n, n, afac, lda, a, lda )
411 lwork = nmax*max( 3, nrhs )
413 CALL dgetri( n, a, lda, iwork, work, lwork, info )
417 ainvnm = dlange(
'1', n, n, a, lda, rwork )
418 IF( anormo.LE.zero .OR. ainvnm.LE.zero )
THEN
421 rcondo = ( one / anormo ) / ainvnm
426 ainvnm = dlange(
'I', n, n, a, lda, rwork )
427 IF( anormi.LE.zero .OR. ainvnm.LE.zero )
THEN
430 rcondi = ( one / anormi ) / ainvnm
434 DO 50 itran = 1, ntran
438 trans = transs( itran )
439 IF( itran.EQ.1 )
THEN
447 CALL dlacpy(
'Full', n, n, asav, lda, a, lda )
452 CALL dlarhs( path, xtype,
'Full', trans, n, n, kl,
453 $ ku, nrhs, a, lda, xact, lda, b, lda,
456 CALL dlacpy(
'Full', n, nrhs, b, lda, bsav, lda )
458 IF( nofact .AND. itran.EQ.1 )
THEN
465 CALL dlacpy(
'Full', n, n, a, lda, afac, lda )
466 CALL dlacpy(
'Full', n, nrhs, b, lda, x, lda )
469 CALL dgesv( n, nrhs, afac, lda, iwork, x, lda,
475 $
CALL alaerh( path,
'DGESV ', info, izero,
476 $
' ', n, n, -1, -1, nrhs, imat,
477 $ nfail, nerrs, nout )
482 CALL dget01( n, n, a, lda, afac, lda, iwork,
483 $ rwork, result( 1 ) )
485 IF( izero.EQ.0 )
THEN
489 CALL dlacpy(
'Full', n, nrhs, b, lda, work,
491 CALL dget02(
'No transpose', n, n, nrhs, a,
492 $ lda, x, lda, work, lda, rwork,
497 CALL dget04( n, nrhs, x, lda, xact, lda,
498 $ rcondc, result( 3 ) )
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 )
'DGESV ', n,
510 $ imat, k, result( k )
520 $
CALL dlaset(
'Full', n, n, zero, zero, afac,
522 CALL dlaset(
'Full', n, nrhs, zero, zero, x, lda )
523 IF( iequed.GT.1 .AND. n.GT.0 )
THEN
528 CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
529 $ colcnd, amax, equed )
536 CALL dgesvx( fact, trans, n, nrhs, a, lda, afac,
537 $ lda, iwork, equed, s, s( n+1 ), b,
538 $ lda, x, lda, rcond, rwork,
539 $ rwork( nrhs+1 ), work, iwork( n+1 ),
545 $
CALL alaerh( path,
'DGESVX', info, izero,
546 $ fact // trans, n, n, -1, -1, nrhs,
547 $ imat, nfail, nerrs, nout )
552 IF( info.NE.0 .AND. info.LE.n)
THEN
553 rpvgrw = dlantr(
'M',
'U',
'N', info, info,
555 IF( rpvgrw.EQ.zero )
THEN
558 rpvgrw = dlange(
'M', n, info, a, lda,
562 rpvgrw = dlantr(
'M',
'U',
'N', n, n, afac, lda,
564 IF( rpvgrw.EQ.zero )
THEN
567 rpvgrw = dlange(
'M', n, n, a, lda, work ) /
571 result( 7 ) = abs( rpvgrw-work( 1 ) ) /
572 $ max( work( 1 ), rpvgrw ) /
575 IF( .NOT.prefac )
THEN
580 CALL dget01( n, n, a, lda, afac, lda, iwork,
581 $ rwork( 2*nrhs+1 ), result( 1 ) )
592 CALL dlacpy(
'Full', n, nrhs, bsav, lda, work,
594 CALL dget02( trans, n, n, nrhs, asav, lda, x,
595 $ lda, work, lda, rwork( 2*nrhs+1 ),
600 IF( nofact .OR. ( prefac .AND. lsame( equed,
602 CALL dget04( n, nrhs, x, lda, xact, lda,
603 $ rcondc, result( 3 ) )
605 IF( itran.EQ.1 )
THEN
610 CALL dget04( n, nrhs, x, lda, xact, lda,
611 $ roldc, result( 3 ) )
617 CALL dget07( trans, n, nrhs, asav, lda, b, lda,
618 $ x, lda, xact, lda, rwork, .true.,
619 $ rwork( nrhs+1 ), result( 4 ) )
627 result( 6 ) = dget06( rcond, rcondc )
632 IF( .NOT.trfcon )
THEN
634 IF( result( k ).GE.thresh )
THEN
635 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
636 $
CALL aladhd( nout, path )
638 WRITE( nout, fmt = 9997 )
'DGESVX',
639 $ fact, trans, n, equed, imat, k,
642 WRITE( nout, fmt = 9998 )
'DGESVX',
643 $ fact, trans, n, imat, k, result( k )
648 nrun = nrun + ntests - k1 + 1
650 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
652 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
653 $
CALL aladhd( nout, path )
655 WRITE( nout, fmt = 9997 )
'DGESVX', fact,
656 $ trans, n, equed, imat, 1, result( 1 )
658 WRITE( nout, fmt = 9998 )
'DGESVX', fact,
659 $ trans, n, imat, 1, result( 1 )
664 IF( result( 6 ).GE.thresh )
THEN
665 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
666 $
CALL aladhd( nout, path )
668 WRITE( nout, fmt = 9997 )
'DGESVX', fact,
669 $ trans, n, equed, imat, 6, result( 6 )
671 WRITE( nout, fmt = 9998 )
'DGESVX', fact,
672 $ trans, n, imat, 6, result( 6 )
677 IF( result( 7 ).GE.thresh )
THEN
678 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
679 $
CALL aladhd( nout, path )
681 WRITE( nout, fmt = 9997 )
'DGESVX', fact,
682 $ trans, n, equed, imat, 7, result( 7 )
684 WRITE( nout, fmt = 9998 )
'DGESVX', fact,
685 $ trans, n, imat, 7, result( 7 )
701 CALL alasvm( path, nout, nfail, nrun, nerrs )
703 9999
FORMAT( 1x, a,
', N =', i5,
', type ', i2,
', test(', i2,
') =',
705 9998
FORMAT( 1x, a,
', FACT=''', a1,
''', TRANS=''', a1,
''', N=', i5,
706 $
', type ', i2,
', test(', i1,
')=', g12.5 )
707 9997
FORMAT( 1x, a,
', FACT=''', a1,
''', TRANS=''', a1,
''', N=', i5,
708 $
', EQUED=''', a1,
''', type ', i2,
', test(', i1,
')=',
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
subroutine dget02(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DGET02
subroutine dlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
DLARHS
subroutine xlaenv(ispec, nvalue)
XLAENV
subroutine aladhd(iounit, path)
ALADHD
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
subroutine ddrvge(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, asav, b, bsav, x, xact, s, work, rwork, iwork, nout)
DDRVGE
subroutine derrvx(path, nunit)
DERRVX
subroutine dget01(m, n, a, lda, afac, ldafac, ipiv, rwork, resid)
DGET01
subroutine dget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
DGET04
subroutine dget07(trans, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, chkferr, berr, reslts)
DGET07
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
subroutine dgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
DGEEQU
subroutine dgesv(n, nrhs, a, lda, ipiv, b, ldb, info)
Download DGESV + dependencies <a href="http://www.netlib.org/cgi-bin/netlibfiles....
subroutine dgesvx(fact, trans, n, nrhs, a, lda, af, ldaf, ipiv, equed, r, c, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
DGESVX computes the solution to system of linear equations A * X = B for GE matrices
subroutine dgetrf(m, n, a, lda, ipiv, info)
DGETRF
subroutine dgetri(n, a, lda, ipiv, work, lwork, info)
DGETRI
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
subroutine dlaqge(m, n, a, lda, r, c, rowcnd, colcnd, amax, equed)
DLAQGE scales a general rectangular matrix, using row and column scaling factors computed by sgeequ.
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.