LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ ddrvge()

subroutine ddrvge ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nrhs,
double precision  thresh,
logical  tsterr,
integer  nmax,
double precision, dimension( * )  a,
double precision, dimension( * )  afac,
double precision, dimension( * )  asav,
double precision, dimension( * )  b,
double precision, dimension( * )  bsav,
double precision, dimension( * )  x,
double precision, dimension( * )  xact,
double precision, dimension( * )  s,
double precision, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

DDRVGE

Purpose:
 DDRVGE tests the driver routines DGESV 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 DOUBLE PRECISION
          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 DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]ASAV
          ASAV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]S
          S is DOUBLE PRECISION array, dimension (2*NMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (2*NRHS+NMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (2*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.

Definition at line 161 of file ddrvge.f.

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 DOUBLE PRECISION THRESH
173* ..
174* .. Array Arguments ..
175 LOGICAL DOTYPE( * )
176 INTEGER IWORK( * ), NVAL( * )
177 DOUBLE PRECISION A( * ), AFAC( * ), ASAV( * ), B( * ),
178 $ BSAV( * ), RWORK( * ), S( * ), WORK( * ),
179 $ X( * ), XACT( * )
180* ..
181*
182* =====================================================================
183*
184* .. Parameters ..
185 DOUBLE PRECISION ONE, ZERO
186 parameter( one = 1.0d+0, zero = 0.0d+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 DOUBLE PRECISION 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 DOUBLE PRECISION RESULT( NTESTS )
209* ..
210* .. External Functions ..
211 LOGICAL LSAME
212 DOUBLE PRECISION DGET06, DLAMCH, DLANGE, DLANTR
213 EXTERNAL lsame, dget06, dlamch, dlange, dlantr
214* ..
215* .. External Subroutines ..
216 EXTERNAL aladhd, alaerh, alasvm, derrvx, dgeequ, dgesv,
219 $ dlatms, xlaenv
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC abs, 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 ) = 'Double 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 derrvx( 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 DLATB4 and generate a test matrix
289* with DLATMS.
290*
291 CALL dlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
292 $ CNDNUM, DIST )
293 rcondc = one / cndnum
294*
295 srnamt = 'DLATMS'
296 CALL dlatms( 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 DLATMS.
301*
302 IF( info.NE.0 ) THEN
303 CALL alaerh( path, 'DLATMS', 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 dlaset( 'Full', n, n-izero+1, zero, zero,
326 $ 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 dlacpy( '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 DGESVX (FACT = 'N' reuses
360* the condition number from the previous iteration
361* with FACT = 'F').
362*
363 CALL dlacpy( '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 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
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 dlaqge( 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 DGET04.
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 = dlange( '1', n, n, afac, lda, rwork )
401 anormi = dlange( 'I', n, n, afac, lda, rwork )
402*
403* Factor the matrix A.
404*
405 srnamt = 'DGETRF'
406 CALL dgetrf( n, n, afac, lda, iwork, info )
407*
408* Form the inverse of A.
409*
410 CALL dlacpy( 'Full', n, n, afac, lda, a, lda )
411 lwork = nmax*max( 3, nrhs )
412 srnamt = 'DGETRI'
413 CALL dgetri( n, a, lda, iwork, work, lwork, info )
414*
415* Compute the 1-norm condition number of A.
416*
417 ainvnm = dlange( '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 = dlange( '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 dlacpy( 'Full', n, n, asav, lda, a, lda )
448*
449* Form an exact solution and set the right hand side.
450*
451 srnamt = 'DLARHS'
452 CALL dlarhs( path, xtype, 'Full', trans, n, n, kl,
453 $ ku, nrhs, a, lda, xact, lda, b, lda,
454 $ iseed, info )
455 xtype = 'C'
456 CALL dlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
457*
458 IF( nofact .AND. itran.EQ.1 ) THEN
459*
460* --- Test DGESV ---
461*
462* Compute the LU factorization of the matrix and
463* solve the system.
464*
465 CALL dlacpy( 'Full', n, n, a, lda, afac, lda )
466 CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
467*
468 srnamt = 'DGESV '
469 CALL dgesv( n, nrhs, afac, lda, iwork, x, lda,
470 $ info )
471*
472* Check error code from DGESV .
473*
474 IF( info.NE.izero )
475 $ CALL alaerh( path, 'DGESV ', 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 dget01( 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 dlacpy( 'Full', n, nrhs, b, lda, work,
490 $ lda )
491 CALL dget02( '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 dget04( 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 )'DGESV ', 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 DGESVX ---
518*
519 IF( .NOT.prefac )
520 $ CALL dlaset( 'Full', n, n, zero, zero, afac,
521 $ lda )
522 CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
523 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
524*
525* Equilibrate the matrix if FACT = 'F' and
526* EQUED = 'R', 'C', or 'B'.
527*
528 CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
529 $ colcnd, amax, equed )
530 END IF
531*
532* Solve the system and compute the condition number
533* and error bounds using DGESVX.
534*
535 srnamt = 'DGESVX'
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 ),
540 $ info )
541*
542* Check the error code from DGESVX.
543*
544 IF( info.NE.izero )
545 $ CALL alaerh( path, 'DGESVX', info, izero,
546 $ fact // trans, n, n, -1, -1, nrhs,
547 $ imat, nfail, nerrs, nout )
548*
549* Compare WORK(1) from DGESVX with the computed
550* reciprocal pivot growth factor RPVGRW
551*
552 IF( info.NE.0 .AND. info.LE.n) THEN
553 rpvgrw = dlantr( 'M', 'U', 'N', info, info,
554 $ afac, lda, work )
555 IF( rpvgrw.EQ.zero ) THEN
556 rpvgrw = one
557 ELSE
558 rpvgrw = dlange( 'M', n, info, a, lda,
559 $ work ) / rpvgrw
560 END IF
561 ELSE
562 rpvgrw = dlantr( 'M', 'U', 'N', n, n, afac, lda,
563 $ work )
564 IF( rpvgrw.EQ.zero ) THEN
565 rpvgrw = one
566 ELSE
567 rpvgrw = dlange( 'M', n, n, a, lda, work ) /
568 $ rpvgrw
569 END IF
570 END IF
571 result( 7 ) = abs( rpvgrw-work( 1 ) ) /
572 $ max( work( 1 ), rpvgrw ) /
573 $ dlamch( 'E' )
574*
575 IF( .NOT.prefac ) THEN
576*
577* Reconstruct matrix from factors and compute
578* residual.
579*
580 CALL dget01( n, n, a, lda, afac, lda, iwork,
581 $ rwork( 2*nrhs+1 ), result( 1 ) )
582 k1 = 1
583 ELSE
584 k1 = 2
585 END IF
586*
587 IF( info.EQ.0 ) THEN
588 trfcon = .false.
589*
590* Compute residual of the computed solution.
591*
592 CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
593 $ lda )
594 CALL dget02( trans, n, n, nrhs, asav, lda, x,
595 $ lda, work, lda, rwork( 2*nrhs+1 ),
596 $ result( 2 ) )
597*
598* Check solution from generated exact solution.
599*
600 IF( nofact .OR. ( prefac .AND. lsame( equed,
601 $ 'N' ) ) ) THEN
602 CALL dget04( n, nrhs, x, lda, xact, lda,
603 $ rcondc, result( 3 ) )
604 ELSE
605 IF( itran.EQ.1 ) THEN
606 roldc = roldo
607 ELSE
608 roldc = roldi
609 END IF
610 CALL dget04( n, nrhs, x, lda, xact, lda,
611 $ roldc, result( 3 ) )
612 END IF
613*
614* Check the error bounds from iterative
615* refinement.
616*
617 CALL dget07( trans, n, nrhs, asav, lda, b, lda,
618 $ x, lda, xact, lda, rwork, .true.,
619 $ rwork( nrhs+1 ), result( 4 ) )
620 ELSE
621 trfcon = .true.
622 END IF
623*
624* Compare RCOND from DGESVX with the computed value
625* in RCONDC.
626*
627 result( 6 ) = dget06( rcond, rcondc )
628*
629* Print information about the tests that did not pass
630* the threshold.
631*
632 IF( .NOT.trfcon ) THEN
633 DO 40 k = k1, ntests
634 IF( result( k ).GE.thresh ) THEN
635 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
636 $ CALL aladhd( nout, path )
637 IF( prefac ) THEN
638 WRITE( nout, fmt = 9997 )'DGESVX',
639 $ fact, trans, n, equed, imat, k,
640 $ result( k )
641 ELSE
642 WRITE( nout, fmt = 9998 )'DGESVX',
643 $ fact, trans, n, imat, k, result( k )
644 END IF
645 nfail = nfail + 1
646 END IF
647 40 CONTINUE
648 nrun = nrun + ntests - k1 + 1
649 ELSE
650 IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
651 $ THEN
652 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
653 $ CALL aladhd( nout, path )
654 IF( prefac ) THEN
655 WRITE( nout, fmt = 9997 )'DGESVX', fact,
656 $ trans, n, equed, imat, 1, result( 1 )
657 ELSE
658 WRITE( nout, fmt = 9998 )'DGESVX', fact,
659 $ trans, n, imat, 1, result( 1 )
660 END IF
661 nfail = nfail + 1
662 nrun = nrun + 1
663 END IF
664 IF( result( 6 ).GE.thresh ) THEN
665 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
666 $ CALL aladhd( nout, path )
667 IF( prefac ) THEN
668 WRITE( nout, fmt = 9997 )'DGESVX', fact,
669 $ trans, n, equed, imat, 6, result( 6 )
670 ELSE
671 WRITE( nout, fmt = 9998 )'DGESVX', fact,
672 $ trans, n, imat, 6, result( 6 )
673 END IF
674 nfail = nfail + 1
675 nrun = nrun + 1
676 END IF
677 IF( result( 7 ).GE.thresh ) THEN
678 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
679 $ CALL aladhd( nout, path )
680 IF( prefac ) THEN
681 WRITE( nout, fmt = 9997 )'DGESVX', fact,
682 $ trans, n, equed, imat, 7, result( 7 )
683 ELSE
684 WRITE( nout, fmt = 9998 )'DGESVX', fact,
685 $ trans, n, imat, 7, result( 7 )
686 END IF
687 nfail = nfail + 1
688 nrun = nrun + 1
689 END IF
690*
691 END IF
692*
693 50 CONTINUE
694 60 CONTINUE
695 70 CONTINUE
696 80 CONTINUE
697 90 CONTINUE
698*
699* Print a summary of the results.
700*
701 CALL alasvm( path, nout, nfail, nrun, nerrs )
702*
703 9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test(', i2, ') =',
704 $ g12.5 )
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, ')=',
709 $ g12.5 )
710 RETURN
711*
712* End of DDRVGE
713*
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine dget02(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DGET02
Definition dget02.f:135
subroutine dlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
DLARHS
Definition dlarhs.f:205
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 derrvx(path, nunit)
DERRVX
Definition derrvx.f:55
subroutine dget01(m, n, a, lda, afac, ldafac, ipiv, rwork, resid)
DGET01
Definition dget01.f:107
subroutine dget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
DGET04
Definition dget04.f:102
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
subroutine dget07(trans, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, chkferr, berr, reslts)
DGET07
Definition dget07.f:165
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
Definition dlatb4.f:120
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
DGEEQU
Definition dgeequ.f:139
subroutine dgesv(n, nrhs, a, lda, ipiv, b, ldb, info)
Download DGESV + dependencies <a href="http://www.netlib.org/cgi-bin/netlibfiles....
Definition dgesv.f:124
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
Definition dgesvx.f:349
subroutine dgetrf(m, n, a, lda, ipiv, info)
DGETRF
Definition dgetrf.f:108
subroutine dgetri(n, a, lda, ipiv, work, lwork, info)
DGETRI
Definition dgetri.f:114
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
double precision function dlantr(norm, uplo, diag, m, n, a, lda, work)
DLANTR returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlantr.f:141
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.
Definition dlaqge.f:142
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.
Definition dlaset.f:110
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: