LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dchkge.f
Go to the documentation of this file.
1*> \brief \b DCHKGE
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 DCHKGE( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
12* NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B,
13* X, XACT, WORK, RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NM, NMAX, NN, NNB, NNS, NOUT
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
23* $ NVAL( * )
24* DOUBLE PRECISION A( * ), AFAC( * ), AINV( * ), B( * ),
25* $ RWORK( * ), WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> DCHKGE tests DGETRF, -TRI, -TRS, -RFS, and -CON.
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] NM
49*> \verbatim
50*> NM is INTEGER
51*> The number of values of M contained in the vector MVAL.
52*> \endverbatim
53*>
54*> \param[in] MVAL
55*> \verbatim
56*> MVAL is INTEGER array, dimension (NM)
57*> The values of the matrix row dimension M.
58*> \endverbatim
59*>
60*> \param[in] NN
61*> \verbatim
62*> NN is INTEGER
63*> The number of values of N contained in the vector NVAL.
64*> \endverbatim
65*>
66*> \param[in] NVAL
67*> \verbatim
68*> NVAL is INTEGER array, dimension (NN)
69*> The values of the matrix column dimension N.
70*> \endverbatim
71*>
72*> \param[in] NNB
73*> \verbatim
74*> NNB is INTEGER
75*> The number of values of NB contained in the vector NBVAL.
76*> \endverbatim
77*>
78*> \param[in] NBVAL
79*> \verbatim
80*> NBVAL is INTEGER array, dimension (NNB)
81*> The values of the blocksize NB.
82*> \endverbatim
83*>
84*> \param[in] NNS
85*> \verbatim
86*> NNS is INTEGER
87*> The number of values of NRHS contained in the vector NSVAL.
88*> \endverbatim
89*>
90*> \param[in] NSVAL
91*> \verbatim
92*> NSVAL is INTEGER array, dimension (NNS)
93*> The values of the number of right hand sides NRHS.
94*> \endverbatim
95*>
96*> \param[in] THRESH
97*> \verbatim
98*> THRESH is DOUBLE PRECISION
99*> The threshold value for the test ratios. A result is
100*> included in the output file if RESULT >= THRESH. To have
101*> every test ratio printed, use THRESH = 0.
102*> \endverbatim
103*>
104*> \param[in] TSTERR
105*> \verbatim
106*> TSTERR is LOGICAL
107*> Flag that indicates whether error exits are to be tested.
108*> \endverbatim
109*>
110*> \param[in] NMAX
111*> \verbatim
112*> NMAX is INTEGER
113*> The maximum value permitted for M or N, used in dimensioning
114*> the work arrays.
115*> \endverbatim
116*>
117*> \param[out] A
118*> \verbatim
119*> A is DOUBLE PRECISION array, dimension (NMAX*NMAX)
120*> \endverbatim
121*>
122*> \param[out] AFAC
123*> \verbatim
124*> AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
125*> \endverbatim
126*>
127*> \param[out] AINV
128*> \verbatim
129*> AINV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
130*> \endverbatim
131*>
132*> \param[out] B
133*> \verbatim
134*> B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
135*> where NSMAX is the largest entry in NSVAL.
136*> \endverbatim
137*>
138*> \param[out] X
139*> \verbatim
140*> X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
141*> \endverbatim
142*>
143*> \param[out] XACT
144*> \verbatim
145*> XACT is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
146*> \endverbatim
147*>
148*> \param[out] WORK
149*> \verbatim
150*> WORK is DOUBLE PRECISION array, dimension
151*> (NMAX*max(3,NSMAX))
152*> \endverbatim
153*>
154*> \param[out] RWORK
155*> \verbatim
156*> RWORK is DOUBLE PRECISION array, dimension
157*> (max(2*NMAX,2*NSMAX+NWORK))
158*> \endverbatim
159*>
160*> \param[out] IWORK
161*> \verbatim
162*> IWORK is INTEGER array, dimension (2*NMAX)
163*> \endverbatim
164*>
165*> \param[in] NOUT
166*> \verbatim
167*> NOUT is INTEGER
168*> The unit number for output.
169*> \endverbatim
170*
171* Authors:
172* ========
173*
174*> \author Univ. of Tennessee
175*> \author Univ. of California Berkeley
176*> \author Univ. of Colorado Denver
177*> \author NAG Ltd.
178*
179*> \ingroup double_lin
180*
181* =====================================================================
182 SUBROUTINE dchkge( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
183 $ NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B,
184 $ X, XACT, WORK, RWORK, IWORK, NOUT )
185*
186* -- LAPACK test routine --
187* -- LAPACK is a software package provided by Univ. of Tennessee, --
188* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189*
190* .. Scalar Arguments ..
191 LOGICAL TSTERR
192 INTEGER NM, NMAX, NN, NNB, NNS, NOUT
193 DOUBLE PRECISION THRESH
194* ..
195* .. Array Arguments ..
196 LOGICAL DOTYPE( * )
197 INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
198 $ nval( * )
199 DOUBLE PRECISION A( * ), AFAC( * ), AINV( * ), B( * ),
200 $ RWORK( * ), WORK( * ), X( * ), XACT( * )
201* ..
202*
203* =====================================================================
204*
205* .. Parameters ..
206 DOUBLE PRECISION ONE, ZERO
207 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
208 INTEGER NTYPES
209 parameter( ntypes = 11 )
210 INTEGER NTESTS
211 parameter( ntests = 8 )
212 INTEGER NTRAN
213 parameter( ntran = 3 )
214* ..
215* .. Local Scalars ..
216 LOGICAL TRFCON, ZEROT
217 CHARACTER DIST, NORM, TRANS, TYPE, XTYPE
218 CHARACTER*3 PATH
219 INTEGER I, IM, IMAT, IN, INB, INFO, IOFF, IRHS, ITRAN,
220 $ izero, k, kl, ku, lda, lwork, m, mode, n, nb,
221 $ nerrs, nfail, nimat, nrhs, nrun, nt
222 DOUBLE PRECISION AINVNM, ANORM, ANORMI, ANORMO, CNDNUM, DUMMY,
223 $ RCOND, RCONDC, RCONDI, RCONDO
224* ..
225* .. Local Arrays ..
226 CHARACTER TRANSS( NTRAN )
227 INTEGER ISEED( 4 ), ISEEDY( 4 )
228 DOUBLE PRECISION RESULT( NTESTS )
229* ..
230* .. External Functions ..
231 DOUBLE PRECISION DGET06, DLANGE
232 EXTERNAL DGET06, DLANGE
233* ..
234* .. External Subroutines ..
235 EXTERNAL alaerh, alahd, alasum, derrge, dgecon, dgerfs,
238 $ dlatms, xlaenv
239* ..
240* .. Intrinsic Functions ..
241 INTRINSIC max, min
242* ..
243* .. Scalars in Common ..
244 LOGICAL LERR, OK
245 CHARACTER*32 SRNAMT
246 INTEGER INFOT, NUNIT
247* ..
248* .. Common blocks ..
249 COMMON / infoc / infot, nunit, ok, lerr
250 COMMON / srnamc / srnamt
251* ..
252* .. Data statements ..
253 DATA iseedy / 1988, 1989, 1990, 1991 / ,
254 $ transs / 'N', 'T', 'C' /
255* ..
256* .. Executable Statements ..
257*
258* Initialize constants and the random number seed.
259*
260 path( 1: 1 ) = 'Double precision'
261 path( 2: 3 ) = 'GE'
262 nrun = 0
263 nfail = 0
264 nerrs = 0
265 DO 10 i = 1, 4
266 iseed( i ) = iseedy( i )
267 10 CONTINUE
268*
269* Test the error exits
270*
271 CALL xlaenv( 1, 1 )
272 IF( tsterr )
273 $ CALL derrge( path, nout )
274 infot = 0
275 CALL xlaenv( 2, 2 )
276*
277* Do for each value of M in MVAL
278*
279 DO 120 im = 1, nm
280 m = mval( im )
281 lda = max( 1, m )
282*
283* Do for each value of N in NVAL
284*
285 DO 110 in = 1, nn
286 n = nval( in )
287 xtype = 'N'
288 nimat = ntypes
289 IF( m.LE.0 .OR. n.LE.0 )
290 $ nimat = 1
291*
292 DO 100 imat = 1, nimat
293*
294* Do the tests only if DOTYPE( IMAT ) is true.
295*
296 IF( .NOT.dotype( imat ) )
297 $ GO TO 100
298*
299* Skip types 5, 6, or 7 if the matrix size is too small.
300*
301 zerot = imat.GE.5 .AND. imat.LE.7
302 IF( zerot .AND. n.LT.imat-4 )
303 $ GO TO 100
304*
305* Set up parameters with DLATB4 and generate a test matrix
306* with DLATMS.
307*
308 CALL dlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
309 $ cndnum, dist )
310*
311 srnamt = 'DLATMS'
312 CALL dlatms( m, n, dist, iseed, TYPE, rwork, mode,
313 $ cndnum, anorm, kl, ku, 'No packing', a, lda,
314 $ work, info )
315*
316* Check error code from DLATMS.
317*
318 IF( info.NE.0 ) THEN
319 CALL alaerh( path, 'DLATMS', info, 0, ' ', m, n, -1,
320 $ -1, -1, imat, nfail, nerrs, nout )
321 GO TO 100
322 END IF
323*
324* For types 5-7, zero one or more columns of the matrix to
325* test that INFO is returned correctly.
326*
327 IF( zerot ) THEN
328 IF( imat.EQ.5 ) THEN
329 izero = 1
330 ELSE IF( imat.EQ.6 ) THEN
331 izero = min( m, n )
332 ELSE
333 izero = min( m, n ) / 2 + 1
334 END IF
335 ioff = ( izero-1 )*lda
336 IF( imat.LT.7 ) THEN
337 DO 20 i = 1, m
338 a( ioff+i ) = zero
339 20 CONTINUE
340 ELSE
341 CALL dlaset( 'Full', m, n-izero+1, zero, zero,
342 $ a( ioff+1 ), lda )
343 END IF
344 ELSE
345 izero = 0
346 END IF
347*
348* These lines, if used in place of the calls in the DO 60
349* loop, cause the code to bomb on a Sun SPARCstation.
350*
351* ANORMO = DLANGE( 'O', M, N, A, LDA, RWORK )
352* ANORMI = DLANGE( 'I', M, N, A, LDA, RWORK )
353*
354* Do for each blocksize in NBVAL
355*
356 DO 90 inb = 1, nnb
357 nb = nbval( inb )
358 CALL xlaenv( 1, nb )
359*
360* Compute the LU factorization of the matrix.
361*
362 CALL dlacpy( 'Full', m, n, a, lda, afac, lda )
363 srnamt = 'DGETRF'
364 CALL dgetrf( m, n, afac, lda, iwork, info )
365*
366* Check error code from DGETRF.
367*
368 IF( info.NE.izero )
369 $ CALL alaerh( path, 'DGETRF', info, izero, ' ', m,
370 $ n, -1, -1, nb, imat, nfail, nerrs,
371 $ nout )
372 trfcon = .false.
373*
374*+ TEST 1
375* Reconstruct matrix from factors and compute residual.
376*
377 CALL dlacpy( 'Full', m, n, afac, lda, ainv, lda )
378 CALL dget01( m, n, a, lda, ainv, lda, iwork, rwork,
379 $ result( 1 ) )
380 nt = 1
381*
382*+ TEST 2
383* Form the inverse if the factorization was successful
384* and compute the residual.
385*
386 IF( m.EQ.n .AND. info.EQ.0 ) THEN
387 CALL dlacpy( 'Full', n, n, afac, lda, ainv, lda )
388 srnamt = 'DGETRI'
389 nrhs = nsval( 1 )
390 lwork = nmax*max( 3, nrhs )
391 CALL dgetri( n, ainv, lda, iwork, work, lwork,
392 $ info )
393*
394* Check error code from DGETRI.
395*
396 IF( info.NE.0 )
397 $ CALL alaerh( path, 'DGETRI', info, 0, ' ', n, n,
398 $ -1, -1, nb, imat, nfail, nerrs,
399 $ nout )
400*
401* Compute the residual for the matrix times its
402* inverse. Also compute the 1-norm condition number
403* of A.
404*
405 CALL dget03( n, a, lda, ainv, lda, work, lda,
406 $ rwork, rcondo, result( 2 ) )
407 anormo = dlange( 'O', m, n, a, lda, rwork )
408*
409* Compute the infinity-norm condition number of A.
410*
411 anormi = dlange( 'I', m, n, a, lda, rwork )
412 ainvnm = dlange( 'I', n, n, ainv, lda, rwork )
413 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
414 rcondi = one
415 ELSE
416 rcondi = ( one / anormi ) / ainvnm
417 END IF
418 nt = 2
419 ELSE
420*
421* Do only the condition estimate if INFO > 0.
422*
423 trfcon = .true.
424 anormo = dlange( 'O', m, n, a, lda, rwork )
425 anormi = dlange( 'I', m, n, a, lda, rwork )
426 rcondo = zero
427 rcondi = zero
428 END IF
429*
430* Print information about the tests so far that did not
431* pass the threshold.
432*
433 DO 30 k = 1, nt
434 IF( result( k ).GE.thresh ) THEN
435 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
436 $ CALL alahd( nout, path )
437 WRITE( nout, fmt = 9999 )m, n, nb, imat, k,
438 $ result( k )
439 nfail = nfail + 1
440 END IF
441 30 CONTINUE
442 nrun = nrun + nt
443*
444* Skip the remaining tests if this is not the first
445* block size or if M .ne. N. Skip the solve tests if
446* the matrix is singular.
447*
448 IF( inb.GT.1 .OR. m.NE.n )
449 $ GO TO 90
450 IF( trfcon )
451 $ GO TO 70
452*
453 DO 60 irhs = 1, nns
454 nrhs = nsval( irhs )
455 xtype = 'N'
456*
457 DO 50 itran = 1, ntran
458 trans = transs( itran )
459 IF( itran.EQ.1 ) THEN
460 rcondc = rcondo
461 ELSE
462 rcondc = rcondi
463 END IF
464*
465*+ TEST 3
466* Solve and compute residual for A * X = B.
467*
468 srnamt = 'DLARHS'
469 CALL dlarhs( path, xtype, ' ', trans, n, n, kl,
470 $ ku, nrhs, a, lda, xact, lda, b,
471 $ lda, iseed, info )
472 xtype = 'C'
473*
474 CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
475 srnamt = 'DGETRS'
476 CALL dgetrs( trans, n, nrhs, afac, lda, iwork,
477 $ x, lda, info )
478*
479* Check error code from DGETRS.
480*
481 IF( info.NE.0 )
482 $ CALL alaerh( path, 'DGETRS', info, 0, trans,
483 $ n, n, -1, -1, nrhs, imat, nfail,
484 $ nerrs, nout )
485*
486 CALL dlacpy( 'Full', n, nrhs, b, lda, work,
487 $ lda )
488 CALL dget02( trans, n, n, nrhs, a, lda, x, lda,
489 $ work, lda, rwork, result( 3 ) )
490*
491*+ TEST 4
492* Check solution from generated exact solution.
493*
494 CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
495 $ result( 4 ) )
496*
497*+ TESTS 5, 6, and 7
498* Use iterative refinement to improve the
499* solution.
500*
501 srnamt = 'DGERFS'
502 CALL dgerfs( trans, n, nrhs, a, lda, afac, lda,
503 $ iwork, b, lda, x, lda, rwork,
504 $ rwork( nrhs+1 ), work,
505 $ iwork( n+1 ), info )
506*
507* Check error code from DGERFS.
508*
509 IF( info.NE.0 )
510 $ CALL alaerh( path, 'DGERFS', info, 0, trans,
511 $ n, n, -1, -1, nrhs, imat, nfail,
512 $ nerrs, nout )
513*
514 CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
515 $ result( 5 ) )
516 CALL dget07( trans, n, nrhs, a, lda, b, lda, x,
517 $ lda, xact, lda, rwork, .true.,
518 $ rwork( nrhs+1 ), result( 6 ) )
519*
520* Print information about the tests that did not
521* pass the threshold.
522*
523 DO 40 k = 3, 7
524 IF( result( k ).GE.thresh ) THEN
525 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
526 $ CALL alahd( nout, path )
527 WRITE( nout, fmt = 9998 )trans, n, nrhs,
528 $ imat, k, result( k )
529 nfail = nfail + 1
530 END IF
531 40 CONTINUE
532 nrun = nrun + 5
533 50 CONTINUE
534 60 CONTINUE
535*
536*+ TEST 8
537* Get an estimate of RCOND = 1/CNDNUM.
538*
539 70 CONTINUE
540 DO 80 itran = 1, 2
541 IF( itran.EQ.1 ) THEN
542 anorm = anormo
543 rcondc = rcondo
544 norm = 'O'
545 ELSE
546 anorm = anormi
547 rcondc = rcondi
548 norm = 'I'
549 END IF
550 srnamt = 'DGECON'
551 CALL dgecon( norm, n, afac, lda, anorm, rcond,
552 $ work, iwork( n+1 ), info )
553*
554* Check error code from DGECON.
555*
556 IF( info.NE.0 )
557 $ CALL alaerh( path, 'DGECON', info, 0, norm, n,
558 $ n, -1, -1, -1, imat, nfail, nerrs,
559 $ nout )
560*
561* This line is needed on a Sun SPARCstation.
562*
563 dummy = rcond
564*
565 result( 8 ) = dget06( rcond, rcondc )
566*
567* Print information about the tests that did not pass
568* the threshold.
569*
570 IF( result( 8 ).GE.thresh ) THEN
571 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
572 $ CALL alahd( nout, path )
573 WRITE( nout, fmt = 9997 )norm, n, imat, 8,
574 $ result( 8 )
575 nfail = nfail + 1
576 END IF
577 nrun = nrun + 1
578 80 CONTINUE
579 90 CONTINUE
580 100 CONTINUE
581 110 CONTINUE
582 120 CONTINUE
583*
584* Print a summary of the results.
585*
586 CALL alasum( path, nout, nfail, nrun, nerrs )
587*
588 9999 FORMAT( ' M = ', i5, ', N =', i5, ', NB =', i4, ', type ', i2,
589 $ ', test(', i2, ') =', g12.5 )
590 9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
591 $ i2, ', test(', i2, ') =', g12.5 )
592 9997 FORMAT( ' NORM =''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
593 $ ', test(', i2, ') =', g12.5 )
594 RETURN
595*
596* End of DCHKGE
597*
598 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.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 alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine dchkge(dotype, nm, mval, nn, nval, nnb, nbval, nns, nsval, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
DCHKGE
Definition dchkge.f:185
subroutine derrge(path, nunit)
DERRGE
Definition derrge.f:55
subroutine dget01(m, n, a, lda, afac, ldafac, ipiv, rwork, resid)
DGET01
Definition dget01.f:107
subroutine dget03(n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
DGET03
Definition dget03.f:109
subroutine dget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
DGET04
Definition dget04.f:102
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 dgecon(norm, n, a, lda, anorm, rcond, work, iwork, info)
DGECON
Definition dgecon.f:132
subroutine dgerfs(trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DGERFS
Definition dgerfs.f:185
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 dgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
DGETRS
Definition dgetrs.f:121
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