164 SUBROUTINE ddrvge( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
165 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
166 $ RWORK, IWORK, NOUT )
174 INTEGER NMAX, NN, NOUT, NRHS
175 DOUBLE PRECISION THRESH
179 INTEGER IWORK( * ), NVAL( * )
180 DOUBLE PRECISION A( * ), AFAC( * ), ASAV( * ), B( * ),
181 $ bsav( * ), rwork( * ), s( * ), work( * ),
188 DOUBLE PRECISION ONE, ZERO
189 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
191 parameter( ntypes = 11 )
193 parameter( ntests = 7 )
195 parameter( ntran = 3 )
198 LOGICAL EQUIL, NOFACT, PREFAC, TRFCON, ZEROT
199 CHARACTER DIST, EQUED, FACT, TRANS,
TYPE, XTYPE
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,
205 DOUBLE PRECISION AINVNM, AMAX, ANORM, ANORMI, ANORMO, CNDNUM,
206 $ COLCND, RCOND, RCONDC, RCONDI, RCONDO, ROLDC,
207 $ roldi, roldo, rowcnd, rpvgrw, rpvgrw_svxx
210 CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN )
211 INTEGER ISEED( 4 ), ISEEDY( 4 )
212 DOUBLE PRECISION RESULT( NTESTS ), BERR( NRHS ),
213 $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
217 DOUBLE PRECISION DGET06, DLAMCH, DLANGE, DLANTR, DLA_GERPVGRW
218 EXTERNAL lsame, dget06, dlamch, dlange, dlantr,
236 COMMON / infoc / infot, nunit, ok, lerr
237 COMMON / srnamc / srnamt
240 DATA iseedy / 1988, 1989, 1990, 1991 /
241 DATA transs /
'N',
'T',
'C' /
242 DATA facts /
'F',
'N',
'E' /
243 DATA equeds /
'N',
'R',
'C',
'B' /
249 path( 1: 1 ) =
'Double precision'
255 iseed( i ) = iseedy( i )
261 $
CALL derrvx( path, nout )
281 DO 80 imat = 1, nimat
285 IF( .NOT.dotype( imat ) )
290 zerot = imat.GE.5 .AND. imat.LE.7
291 IF( zerot .AND. n.LT.imat-4 )
297 CALL dlatb4( path, imat, n, n,
TYPE, kl, ku, anorm, mode,
299 rcondc = one / cndnum
302 CALL dlatms( n, n, dist, iseed,
TYPE, rwork, mode, cndnum,
303 $ anorm, kl, ku,
'No packing', a, lda, work,
309 CALL alaerh( path,
'DLATMS', info, 0,
' ', n, n, -1, -1,
310 $ -1, imat, nfail, nerrs, nout )
320 ELSE IF( imat.EQ.6 )
THEN
325 ioff = ( izero-1 )*lda
331 CALL dlaset(
'Full', n, n-izero+1, zero, zero,
340 CALL dlacpy(
'Full', n, n, a, lda, asav, lda )
343 equed = equeds( iequed )
344 IF( iequed.EQ.1 )
THEN
350 DO 60 ifact = 1, nfact
351 fact = facts( ifact )
352 prefac = lsame( fact,
'F' )
353 nofact = lsame( fact,
'N' )
354 equil = lsame( fact,
'E' )
362 ELSE IF( .NOT.nofact )
THEN
369 CALL dlacpy(
'Full', n, n, asav, lda, afac, lda )
370 IF( equil .OR. iequed.GT.1 )
THEN
375 CALL dgeequ( n, n, afac, lda, s, s( n+1 ),
376 $ rowcnd, colcnd, amax, info )
377 IF( info.EQ.0 .AND. n.GT.0 )
THEN
378 IF( lsame( equed,
'R' ) )
THEN
381 ELSE IF( lsame( equed,
'C' ) )
THEN
384 ELSE IF( lsame( equed,
'B' ) )
THEN
391 CALL dlaqge( n, n, afac, lda, s, s( n+1 ),
392 $ rowcnd, colcnd, amax, equed )
406 anormo = dlange(
'1', n, n, afac, lda, rwork )
407 anormi = dlange(
'I', n, n, afac, lda, rwork )
411 CALL dgetrf( n, n, afac, lda, iwork, info )
415 CALL dlacpy(
'Full', n, n, afac, lda, a, lda )
416 lwork = nmax*max( 3, nrhs )
417 CALL dgetri( n, a, lda, iwork, work, lwork, info )
421 ainvnm = dlange(
'1', n, n, a, lda, rwork )
422 IF( anormo.LE.zero .OR. ainvnm.LE.zero )
THEN
425 rcondo = ( one / anormo ) / ainvnm
430 ainvnm = dlange(
'I', n, n, a, lda, rwork )
431 IF( anormi.LE.zero .OR. ainvnm.LE.zero )
THEN
434 rcondi = ( one / anormi ) / ainvnm
438 DO 50 itran = 1, ntran
442 trans = transs( itran )
443 IF( itran.EQ.1 )
THEN
451 CALL dlacpy(
'Full', n, n, asav, lda, a, lda )
456 CALL dlarhs( path, xtype,
'Full', trans, n, n, kl,
457 $ ku, nrhs, a, lda, xact, lda, b, lda,
460 CALL dlacpy(
'Full', n, nrhs, b, lda, bsav, lda )
462 IF( nofact .AND. itran.EQ.1 )
THEN
469 CALL dlacpy(
'Full', n, n, a, lda, afac, lda )
470 CALL dlacpy(
'Full', n, nrhs, b, lda, x, lda )
473 CALL dgesv( n, nrhs, afac, lda, iwork, x, lda,
479 $
CALL alaerh( path,
'DGESV ', info, izero,
480 $
' ', n, n, -1, -1, nrhs, imat,
481 $ nfail, nerrs, nout )
486 CALL dget01( n, n, a, lda, afac, lda, iwork,
487 $ rwork, result( 1 ) )
489 IF( izero.EQ.0 )
THEN
493 CALL dlacpy(
'Full', n, nrhs, b, lda, work,
495 CALL dget02(
'No transpose', n, n, nrhs, a,
496 $ lda, x, lda, work, lda, rwork,
501 CALL dget04( n, nrhs, x, lda, xact, lda,
502 $ rcondc, result( 3 ) )
510 IF( result( k ).GE.thresh )
THEN
511 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
512 $
CALL aladhd( nout, path )
513 WRITE( nout, fmt = 9999 )
'DGESV ', n,
514 $ imat, k, result( k )
524 $
CALL dlaset(
'Full', n, n, zero, zero, afac,
526 CALL dlaset(
'Full', n, nrhs, zero, zero, x, lda )
527 IF( iequed.GT.1 .AND. n.GT.0 )
THEN
532 CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
533 $ colcnd, amax, equed )
540 CALL dgesvx( 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, iwork( n+1 ),
549 $
CALL alaerh( path,
'DGESVX', info, izero,
550 $ fact // trans, n, n, -1, -1, nrhs,
551 $ imat, nfail, nerrs, nout )
557 rpvgrw = dlantr(
'M',
'U',
'N', info, info,
559 IF( rpvgrw.EQ.zero )
THEN
562 rpvgrw = dlange(
'M', n, info, a, lda,
566 rpvgrw = dlantr(
'M',
'U',
'N', n, n, afac, lda,
568 IF( rpvgrw.EQ.zero )
THEN
571 rpvgrw = dlange(
'M', n, n, a, lda, work ) /
575 result( 7 ) = abs( rpvgrw-work( 1 ) ) /
576 $ max( work( 1 ), rpvgrw ) /
579 IF( .NOT.prefac )
THEN
584 CALL dget01( n, n, a, lda, afac, lda, iwork,
585 $ rwork( 2*nrhs+1 ), result( 1 ) )
596 CALL dlacpy(
'Full', n, nrhs, bsav, lda, work,
598 CALL dget02( trans, n, n, nrhs, asav, lda, x,
599 $ lda, work, lda, rwork( 2*nrhs+1 ),
604 IF( nofact .OR. ( prefac .AND. lsame( equed,
606 CALL dget04( n, nrhs, x, lda, xact, lda,
607 $ rcondc, result( 3 ) )
609 IF( itran.EQ.1 )
THEN
614 CALL dget04( n, nrhs, x, lda, xact, lda,
615 $ roldc, result( 3 ) )
621 CALL dget07( trans, n, nrhs, asav, lda, b, lda,
622 $ x, lda, xact, lda, rwork, .true.,
623 $ rwork( nrhs+1 ), result( 4 ) )
631 result( 6 ) = dget06( rcond, rcondc )
636 IF( .NOT.trfcon )
THEN
638 IF( result( k ).GE.thresh )
THEN
639 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
640 $
CALL aladhd( nout, path )
642 WRITE( nout, fmt = 9997 )
'DGESVX',
643 $ fact, trans, n, equed, imat, k,
646 WRITE( nout, fmt = 9998 )
'DGESVX',
647 $ fact, trans, n, imat, k, result( k )
654 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
656 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
657 $
CALL aladhd( nout, path )
659 WRITE( nout, fmt = 9997 )
'DGESVX', fact,
660 $ trans, n, equed, imat, 1, result( 1 )
662 WRITE( nout, fmt = 9998 )
'DGESVX', fact,
663 $ trans, n, imat, 1, result( 1 )
668 IF( result( 6 ).GE.thresh )
THEN
669 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
670 $
CALL aladhd( nout, path )
672 WRITE( nout, fmt = 9997 )
'DGESVX', fact,
673 $ trans, n, equed, imat, 6, result( 6 )
675 WRITE( nout, fmt = 9998 )
'DGESVX', fact,
676 $ trans, n, imat, 6, result( 6 )
681 IF( result( 7 ).GE.thresh )
THEN
682 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
683 $
CALL aladhd( nout, path )
685 WRITE( nout, fmt = 9997 )
'DGESVX', fact,
686 $ trans, n, equed, imat, 7, result( 7 )
688 WRITE( nout, fmt = 9998 )
'DGESVX', fact,
689 $ trans, n, imat, 7, result( 7 )
701 CALL dlacpy(
'Full', n, n, asav, lda, a, lda )
702 CALL dlacpy(
'Full', n, nrhs, bsav, lda, b, lda )
705 $
CALL dlaset(
'Full', n, n, zero, zero, afac,
707 CALL dlaset(
'Full', n, nrhs, zero, zero, x, lda )
708 IF( iequed.GT.1 .AND. n.GT.0 )
THEN
713 CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
714 $ colcnd, amax, equed )
722 CALL dgesvxx( fact, trans, n, nrhs, a, lda, afac,
723 $ lda, iwork, equed, s, s( n+1 ), b, lda, x,
724 $ lda, rcond, rpvgrw_svxx, berr, n_err_bnds,
725 $ errbnds_n, errbnds_c, 0, zero, work,
726 $ iwork( n+1 ), info )
730 IF( info.EQ.n+1 )
GOTO 50
731 IF( info.NE.izero )
THEN
732 CALL alaerh( path,
'DGESVXX', info, izero,
733 $ fact // trans, n, n, -1, -1, nrhs,
734 $ imat, nfail, nerrs, nout )
742 IF ( info .GT. 0 .AND. info .LT. n+1 )
THEN
743 rpvgrw = dla_gerpvgrw
744 $ (n, info, a, lda, afac, lda)
746 rpvgrw = dla_gerpvgrw
747 $ (n, n, a, lda, afac, lda)
750 result( 7 ) = abs( rpvgrw-rpvgrw_svxx ) /
751 $ max( rpvgrw_svxx, rpvgrw ) /
754 IF( .NOT.prefac )
THEN
759 CALL dget01( n, n, a, lda, afac, lda, iwork,
760 $ rwork( 2*nrhs+1 ), result( 1 ) )
771 CALL dlacpy(
'Full', n, nrhs, bsav, lda, work,
773 CALL dget02( trans, n, n, nrhs, asav, lda, x,
774 $ lda, work, lda, rwork( 2*nrhs+1 ),
779 IF( nofact .OR. ( prefac .AND. lsame( equed,
781 CALL dget04( n, nrhs, x, lda, xact, lda,
782 $ rcondc, result( 3 ) )
784 IF( itran.EQ.1 )
THEN
789 CALL dget04( n, nrhs, x, lda, xact, lda,
790 $ roldc, result( 3 ) )
799 result( 6 ) = dget06( rcond, rcondc )
804 IF( .NOT.trfcon )
THEN
806 IF( result( k ).GE.thresh )
THEN
807 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
808 $
CALL aladhd( nout, path )
810 WRITE( nout, fmt = 9997 )
'DGESVXX',
811 $ fact, trans, n, equed, imat, k,
814 WRITE( nout, fmt = 9998 )
'DGESVXX',
815 $ fact, trans, n, imat, k, result( k )
822 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
824 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
825 $
CALL aladhd( nout, path )
827 WRITE( nout, fmt = 9997 )
'DGESVXX', fact,
828 $ trans, n, equed, imat, 1, result( 1 )
830 WRITE( nout, fmt = 9998 )
'DGESVXX', fact,
831 $ trans, n, imat, 1, result( 1 )
836 IF( result( 6 ).GE.thresh )
THEN
837 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
838 $
CALL aladhd( nout, path )
840 WRITE( nout, fmt = 9997 )
'DGESVXX', fact,
841 $ trans, n, equed, imat, 6, result( 6 )
843 WRITE( nout, fmt = 9998 )
'DGESVXX', fact,
844 $ trans, n, imat, 6, result( 6 )
849 IF( result( 7 ).GE.thresh )
THEN
850 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
851 $
CALL aladhd( nout, path )
853 WRITE( nout, fmt = 9997 )
'DGESVXX', fact,
854 $ trans, n, equed, imat, 7, result( 7 )
856 WRITE( nout, fmt = 9998 )
'DGESVXX', fact,
857 $ trans, n, imat, 7, result( 7 )
873 CALL alasvm( path, nout, nfail, nrun, nerrs )
880 9999
FORMAT( 1x, a,
', N =', i5,
', type ', i2,
', test(', i2,
') =',
882 9998
FORMAT( 1x, a,
', FACT=''', a1,
''', TRANS=''', a1,
''', N=', i5,
883 $
', type ', i2,
', test(', i1,
')=', g12.5 )
884 9997
FORMAT( 1x, a,
', FACT=''', a1,
''', TRANS=''', a1,
''', N=', i5,
885 $
', 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 debchvxx(thresh, path)
DEBCHVXX
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 dgesvxx(fact, trans, n, nrhs, a, lda, af, ldaf, ipiv, equed, r, c, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
DGESVXX 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.