OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
ddrvpp.f
Go to the documentation of this file.
1*> \brief \b DDRVPP
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 DDRVPP( 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* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NVAL( * )
23* DOUBLE PRECISION A( * ), AFAC( * ), ASAV( * ), B( * ),
24* $ BSAV( * ), RWORK( * ), S( * ), WORK( * ),
25* $ X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> DDRVPP tests the driver routines DPPSV 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 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 DOUBLE PRECISION
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 DOUBLE PRECISION array, dimension
91*> (NMAX*(NMAX+1)/2)
92*> \endverbatim
93*>
94*> \param[out] AFAC
95*> \verbatim
96*> AFAC is DOUBLE PRECISION array, dimension
97*> (NMAX*(NMAX+1)/2)
98*> \endverbatim
99*>
100*> \param[out] ASAV
101*> \verbatim
102*> ASAV is DOUBLE PRECISION array, dimension
103*> (NMAX*(NMAX+1)/2)
104*> \endverbatim
105*>
106*> \param[out] B
107*> \verbatim
108*> B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
109*> \endverbatim
110*>
111*> \param[out] BSAV
112*> \verbatim
113*> BSAV is DOUBLE PRECISION array, dimension (NMAX*NRHS)
114*> \endverbatim
115*>
116*> \param[out] X
117*> \verbatim
118*> X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
119*> \endverbatim
120*>
121*> \param[out] XACT
122*> \verbatim
123*> XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
124*> \endverbatim
125*>
126*> \param[out] S
127*> \verbatim
128*> S is DOUBLE PRECISION array, dimension (NMAX)
129*> \endverbatim
130*>
131*> \param[out] WORK
132*> \verbatim
133*> WORK is DOUBLE PRECISION array, dimension
134*> (NMAX*max(3,NRHS))
135*> \endverbatim
136*>
137*> \param[out] RWORK
138*> \verbatim
139*> RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
140*> \endverbatim
141*>
142*> \param[out] IWORK
143*> \verbatim
144*> IWORK is INTEGER array, dimension (NMAX)
145*> \endverbatim
146*>
147*> \param[in] NOUT
148*> \verbatim
149*> NOUT is INTEGER
150*> The unit number for output.
151*> \endverbatim
152*
153* Authors:
154* ========
155*
156*> \author Univ. of Tennessee
157*> \author Univ. of California Berkeley
158*> \author Univ. of Colorado Denver
159*> \author NAG Ltd.
160*
161*> \ingroup double_lin
162*
163* =====================================================================
164 SUBROUTINE ddrvpp( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
165 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
166 $ RWORK, IWORK, NOUT )
167*
168* -- LAPACK test routine --
169* -- LAPACK is a software package provided by Univ. of Tennessee, --
170* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
171*
172* .. Scalar Arguments ..
173 LOGICAL TSTERR
174 INTEGER NMAX, NN, NOUT, NRHS
175 DOUBLE PRECISION THRESH
176* ..
177* .. Array Arguments ..
178 LOGICAL DOTYPE( * )
179 INTEGER IWORK( * ), NVAL( * )
180 DOUBLE PRECISION A( * ), AFAC( * ), ASAV( * ), B( * ),
181 $ bsav( * ), rwork( * ), s( * ), work( * ),
182 $ x( * ), xact( * )
183* ..
184*
185* =====================================================================
186*
187* .. Parameters ..
188 DOUBLE PRECISION ONE, ZERO
189 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
190 INTEGER NTYPES
191 parameter( ntypes = 9 )
192 INTEGER NTESTS
193 parameter( ntests = 6 )
194* ..
195* .. Local Scalars ..
196 LOGICAL EQUIL, NOFACT, PREFAC, ZEROT
197 CHARACTER DIST, EQUED, FACT, PACKIT, TYPE, UPLO, XTYPE
198 CHARACTER*3 PATH
199 INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
200 $ izero, k, k1, kl, ku, lda, mode, n, nerrs,
201 $ nfact, nfail, nimat, npp, nrun, nt
202 DOUBLE PRECISION AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
203 $ ROLDC, SCOND
204* ..
205* .. Local Arrays ..
206 CHARACTER EQUEDS( 2 ), FACTS( 3 ), PACKS( 2 ), UPLOS( 2 )
207 INTEGER ISEED( 4 ), ISEEDY( 4 )
208 DOUBLE PRECISION RESULT( NTESTS )
209* ..
210* .. External Functions ..
211 LOGICAL LSAME
212 DOUBLE PRECISION DGET06, DLANSP
213 EXTERNAL lsame, dget06, dlansp
214* ..
215* .. External Subroutines ..
216 EXTERNAL aladhd, alaerh, alasvm, dcopy, derrvx, dget04,
219 $ dpptrf, dpptri
220* ..
221* .. Scalars in Common ..
222 LOGICAL LERR, OK
223 CHARACTER*32 SRNAMT
224 INTEGER INFOT, NUNIT
225* ..
226* .. Common blocks ..
227 COMMON / infoc / infot, nunit, ok, lerr
228 COMMON / srnamc / srnamt
229* ..
230* .. Intrinsic Functions ..
231 INTRINSIC max
232* ..
233* .. Data statements ..
234 DATA iseedy / 1988, 1989, 1990, 1991 /
235 DATA uplos / 'U', 'l' / , FACTS / 'f', 'n', 'e' / ,
236 $ PACKS / 'c', 'r' / , EQUEDS / 'n', 'y' /
237* ..
238* .. Executable Statements ..
239*
240* Initialize constants and the random number seed.
241*
242 PATH( 1: 1 ) = 'double precision'
243 PATH( 2: 3 ) = 'pp'
244 NRUN = 0
245 NFAIL = 0
246 NERRS = 0
247 DO 10 I = 1, 4
248 ISEED( I ) = ISEEDY( I )
249 10 CONTINUE
250*
251* Test the error exits
252*
253 IF( TSTERR )
254 $ CALL DERRVX( PATH, NOUT )
255 INFOT = 0
256*
257* Do for each value of N in NVAL
258*
259 DO 140 IN = 1, NN
260 N = NVAL( IN )
261 LDA = MAX( N, 1 )
262 NPP = N*( N+1 ) / 2
263 XTYPE = 'n'
264 NIMAT = NTYPES
265.LE. IF( N0 )
266 $ NIMAT = 1
267*
268 DO 130 IMAT = 1, NIMAT
269*
270* Do the tests only if DOTYPE( IMAT ) is true.
271*
272.NOT. IF( DOTYPE( IMAT ) )
273 $ GO TO 130
274*
275* Skip types 3, 4, or 5 if the matrix size is too small.
276*
277.GE..AND..LE. ZEROT = IMAT3 IMAT5
278.AND..LT. IF( ZEROT NIMAT-2 )
279 $ GO TO 130
280*
281* Do first for UPLO = 'U', then for UPLO = 'L'
282*
283 DO 120 IUPLO = 1, 2
284 UPLO = UPLOS( IUPLO )
285 PACKIT = PACKS( IUPLO )
286*
287* Set up parameters with DLATB4 and generate a test matrix
288* with DLATMS.
289*
290 CALL DLATB4( PATH, IMAT, N, N, TYPE, KL, KU, ANORM, MODE,
291 $ CNDNUM, DIST )
292 RCONDC = ONE / CNDNUM
293*
294 SRNAMT = 'dlatms'
295 CALL DLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
296 $ CNDNUM, ANORM, KL, KU, PACKIT, A, LDA, WORK,
297 $ INFO )
298*
299* Check error code from DLATMS.
300*
301.NE. IF( INFO0 ) THEN
302 CALL ALAERH( PATH, 'dlatms', INFO, 0, UPLO, N, N, -1,
303 $ -1, -1, IMAT, NFAIL, NERRS, NOUT )
304 GO TO 120
305 END IF
306*
307* For types 3-5, zero one row and column of the matrix to
308* test that INFO is returned correctly.
309*
310 IF( ZEROT ) THEN
311.EQ. IF( IMAT3 ) THEN
312 IZERO = 1
313.EQ. ELSE IF( IMAT4 ) THEN
314 IZERO = N
315 ELSE
316 IZERO = N / 2 + 1
317 END IF
318*
319* Set row and column IZERO of A to 0.
320*
321.EQ. IF( IUPLO1 ) THEN
322 IOFF = ( IZERO-1 )*IZERO / 2
323 DO 20 I = 1, IZERO - 1
324 A( IOFF+I ) = ZERO
325 20 CONTINUE
326 IOFF = IOFF + IZERO
327 DO 30 I = IZERO, N
328 A( IOFF ) = ZERO
329 IOFF = IOFF + I
330 30 CONTINUE
331 ELSE
332 IOFF = IZERO
333 DO 40 I = 1, IZERO - 1
334 A( IOFF ) = ZERO
335 IOFF = IOFF + N - I
336 40 CONTINUE
337 IOFF = IOFF - IZERO
338 DO 50 I = IZERO, N
339 A( IOFF+I ) = ZERO
340 50 CONTINUE
341 END IF
342 ELSE
343 IZERO = 0
344 END IF
345*
346* Save a copy of the matrix A in ASAV.
347*
348 CALL DCOPY( NPP, A, 1, ASAV, 1 )
349*
350 DO 110 IEQUED = 1, 2
351 EQUED = EQUEDS( IEQUED )
352.EQ. IF( IEQUED1 ) THEN
353 NFACT = 3
354 ELSE
355 NFACT = 1
356 END IF
357*
358 DO 100 IFACT = 1, NFACT
359 FACT = FACTS( IFACT )
360 PREFAC = LSAME( FACT, 'f' )
361 NOFACT = LSAME( FACT, 'n' )
362 EQUIL = LSAME( FACT, 'e' )
363*
364 IF( ZEROT ) THEN
365 IF( PREFAC )
366 $ GO TO 100
367 RCONDC = ZERO
368*
369.NOT. ELSE IF( LSAME( FACT, 'n' ) ) THEN
370*
371* Compute the condition number for comparison with
372* the value returned by DPPSVX (FACT = 'N' reuses
373* the condition number from the previous iteration
374* with FACT = 'F').
375*
376 CALL DCOPY( NPP, ASAV, 1, AFAC, 1 )
377.OR..GT. IF( EQUIL IEQUED1 ) THEN
378*
379* Compute row and column scale factors to
380* equilibrate the matrix A.
381*
382 CALL DPPEQU( UPLO, N, AFAC, S, SCOND, AMAX,
383 $ INFO )
384.EQ..AND..GT. IF( INFO0 N0 ) THEN
385.GT. IF( IEQUED1 )
386 $ SCOND = ZERO
387*
388* Equilibrate the matrix.
389*
390 CALL DLAQSP( UPLO, N, AFAC, S, SCOND,
391 $ AMAX, EQUED )
392 END IF
393 END IF
394*
395* Save the condition number of the
396* non-equilibrated system for use in DGET04.
397*
398 IF( EQUIL )
399 $ ROLDC = RCONDC
400*
401* Compute the 1-norm of A.
402*
403 ANORM = DLANSP( '1', UPLO, N, AFAC, RWORK )
404*
405* Factor the matrix A.
406*
407 CALL DPPTRF( UPLO, N, AFAC, INFO )
408*
409* Form the inverse of A.
410*
411 CALL DCOPY( NPP, AFAC, 1, A, 1 )
412 CALL DPPTRI( UPLO, N, A, INFO )
413*
414* Compute the 1-norm condition number of A.
415*
416 AINVNM = DLANSP( '1', UPLO, N, A, RWORK )
417.LE..OR..LE. IF( ANORMZERO AINVNMZERO ) THEN
418 RCONDC = ONE
419 ELSE
420 RCONDC = ( ONE / ANORM ) / AINVNM
421 END IF
422 END IF
423*
424* Restore the matrix A.
425*
426 CALL DCOPY( NPP, ASAV, 1, A, 1 )
427*
428* Form an exact solution and set the right hand side.
429*
430 SRNAMT = 'dlarhs'
431 CALL DLARHS( PATH, XTYPE, UPLO, ' ', N, N, KL, KU,
432 $ NRHS, A, LDA, XACT, LDA, B, LDA,
433 $ ISEED, INFO )
434 XTYPE = 'c'
435 CALL DLACPY( 'full', N, NRHS, B, LDA, BSAV, LDA )
436*
437 IF( NOFACT ) THEN
438*
439* --- Test DPPSV ---
440*
441* Compute the L*L' or U'*U factorization of the
442* matrix and solve the system.
443*
444 CALL DCOPY( NPP, A, 1, AFAC, 1 )
445 CALL DLACPY( 'full', N, NRHS, B, LDA, X, LDA )
446*
447 SRNAMT = 'dppsv '
448 CALL DPPSV( UPLO, N, NRHS, AFAC, X, LDA, INFO )
449*
450* Check error code from DPPSV .
451*
452.NE. IF( INFOIZERO ) THEN
453 CALL ALAERH( PATH, 'dppsv ', INFO, IZERO,
454 $ UPLO, N, N, -1, -1, NRHS, IMAT,
455 $ NFAIL, NERRS, NOUT )
456 GO TO 70
457.NE. ELSE IF( INFO0 ) THEN
458 GO TO 70
459 END IF
460*
461* Reconstruct matrix from factors and compute
462* residual.
463*
464 CALL DPPT01( UPLO, N, A, AFAC, RWORK,
465 $ RESULT( 1 ) )
466*
467* Compute residual of the computed solution.
468*
469 CALL DLACPY( 'full', N, NRHS, B, LDA, WORK,
470 $ LDA )
471 CALL DPPT02( UPLO, N, NRHS, A, X, LDA, WORK,
472 $ LDA, RWORK, RESULT( 2 ) )
473*
474* Check solution from generated exact solution.
475*
476 CALL DGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
477 $ RESULT( 3 ) )
478 NT = 3
479*
480* Print information about the tests that did not
481* pass the threshold.
482*
483 DO 60 K = 1, NT
484.GE. IF( RESULT( K )THRESH ) THEN
485.EQ..AND..EQ. IF( NFAIL0 NERRS0 )
486 $ CALL ALADHD( NOUT, PATH )
487 WRITE( NOUT, FMT = 9999 )'dppsv ', UPLO,
488 $ N, IMAT, K, RESULT( K )
489 NFAIL = NFAIL + 1
490 END IF
491 60 CONTINUE
492 NRUN = NRUN + NT
493 70 CONTINUE
494 END IF
495*
496* --- Test DPPSVX ---
497*
498.NOT..AND..GT. IF( PREFAC NPP0 )
499 $ CALL DLASET( 'full', NPP, 1, ZERO, ZERO, AFAC,
500 $ NPP )
501 CALL DLASET( 'full', N, NRHS, ZERO, ZERO, X, LDA )
502.GT..AND..GT. IF( IEQUED1 N0 ) THEN
503*
504* Equilibrate the matrix if FACT='F' and
505* EQUED='Y'.
506*
507 CALL DLAQSP( UPLO, N, A, S, SCOND, AMAX, EQUED )
508 END IF
509*
510* Solve the system and compute the condition number
511* and error bounds using DPPSVX.
512*
513 SRNAMT = 'dppsvx'
514 CALL DPPSVX( FACT, UPLO, N, NRHS, A, AFAC, EQUED,
515 $ S, B, LDA, X, LDA, RCOND, RWORK,
516 $ RWORK( NRHS+1 ), WORK, IWORK, INFO )
517*
518* Check the error code from DPPSVX.
519*
520.NE. IF( INFOIZERO ) THEN
521 CALL ALAERH( PATH, 'dppsvx', INFO, IZERO,
522 $ FACT // UPLO, N, N, -1, -1, NRHS,
523 $ IMAT, NFAIL, NERRS, NOUT )
524 GO TO 90
525 END IF
526*
527.EQ. IF( INFO0 ) THEN
528.NOT. IF( PREFAC ) THEN
529*
530* Reconstruct matrix from factors and compute
531* residual.
532*
533 CALL DPPT01( UPLO, N, A, AFAC,
534 $ RWORK( 2*NRHS+1 ), RESULT( 1 ) )
535 K1 = 1
536 ELSE
537 K1 = 2
538 END IF
539*
540* Compute residual of the computed solution.
541*
542 CALL DLACPY( 'full', N, NRHS, BSAV, LDA, WORK,
543 $ LDA )
544 CALL DPPT02( UPLO, N, NRHS, ASAV, X, LDA, WORK,
545 $ LDA, RWORK( 2*NRHS+1 ),
546 $ RESULT( 2 ) )
547*
548* Check solution from generated exact solution.
549*
550.OR..AND. IF( NOFACT ( PREFAC LSAME( EQUED,
551 $ 'n' ) ) ) THEN
552 CALL DGET04( N, NRHS, X, LDA, XACT, LDA,
553 $ RCONDC, RESULT( 3 ) )
554 ELSE
555 CALL DGET04( N, NRHS, X, LDA, XACT, LDA,
556 $ ROLDC, RESULT( 3 ) )
557 END IF
558*
559* Check the error bounds from iterative
560* refinement.
561*
562 CALL DPPT05( UPLO, N, NRHS, ASAV, B, LDA, X,
563 $ LDA, XACT, LDA, RWORK,
564 $ RWORK( NRHS+1 ), RESULT( 4 ) )
565 ELSE
566 K1 = 6
567 END IF
568*
569* Compare RCOND from DPPSVX with the computed value
570* in RCONDC.
571*
572 RESULT( 6 ) = DGET06( RCOND, RCONDC )
573*
574* Print information about the tests that did not pass
575* the threshold.
576*
577 DO 80 K = K1, 6
578.GE. IF( RESULT( K )THRESH ) THEN
579.EQ..AND..EQ. IF( NFAIL0 NERRS0 )
580 $ CALL ALADHD( NOUT, PATH )
581 IF( PREFAC ) THEN
582 WRITE( NOUT, FMT = 9997 )'dppsvx', FACT,
583 $ UPLO, N, EQUED, IMAT, K, RESULT( K )
584 ELSE
585 WRITE( NOUT, FMT = 9998 )'dppsvx', FACT,
586 $ UPLO, N, IMAT, K, RESULT( K )
587 END IF
588 NFAIL = NFAIL + 1
589 END IF
590 80 CONTINUE
591 NRUN = NRUN + 7 - K1
592 90 CONTINUE
593 100 CONTINUE
594 110 CONTINUE
595 120 CONTINUE
596 130 CONTINUE
597 140 CONTINUE
598*
599* Print a summary of the results.
600*
601 CALL ALASVM( PATH, NOUT, NFAIL, NRUN, NERRS )
602*
603 9999 FORMAT( 1X, A, ', uplo=''', A1, ''', n =', I5, ', type ', I1,
604 $ ', test(', I1, ')=', G12.5 )
605 9998 FORMAT( 1X, A, ', fact=''', A1, ''', uplo=''', A1, ''', n=', I5,
606 $ ', type ', I1, ', test(', I1, ')=', G12.5 )
607 9997 FORMAT( 1X, A, ', fact=''', A1, ''', uplo=''', A1, ''', n=', I5,
608 $ ', equed=''', A1, ''', type ', I1, ', test(', I1, ')=',
609 $ G12.5 )
610 RETURN
611*
612* End of DDRVPP
613*
614 END
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
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
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
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 dlaqsp(uplo, n, ap, s, scond, amax, equed)
DLAQSP scales a symmetric/Hermitian matrix in packed storage, using scaling factors computed by sppeq...
Definition dlaqsp.f:125
subroutine dpptri(uplo, n, ap, info)
DPPTRI
Definition dpptri.f:93
subroutine dpptrf(uplo, n, ap, info)
DPPTRF
Definition dpptrf.f:119
subroutine dppequ(uplo, n, ap, s, scond, amax, info)
DPPEQU
Definition dppequ.f:116
subroutine dppsvx(fact, uplo, n, nrhs, ap, afp, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
DPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices
Definition dppsvx.f:311
subroutine dppsv(uplo, n, nrhs, ap, b, ldb, info)
DPPSV computes the solution to system of linear equations A * X = B for OTHER matrices
Definition dppsv.f:144
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
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 dppt01(uplo, n, a, afac, rwork, resid)
DPPT01
Definition dppt01.f:93
subroutine dppt02(uplo, n, nrhs, a, x, ldx, b, ldb, rwork, resid)
DPPT02
Definition dppt02.f:122
subroutine ddrvpp(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, asav, b, bsav, x, xact, s, work, rwork, iwork, nout)
DDRVPP
Definition ddrvpp.f:167
subroutine derrvx(path, nunit)
DERRVX
Definition derrvx.f:55
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
Definition dlatb4.f:120
subroutine dppt05(uplo, n, nrhs, ap, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
DPPT05
Definition dppt05.f:156
subroutine dget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
DGET04
Definition dget04.f:102
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
#define max(a, b)
Definition macros.h:21