LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ddrvab.f
Go to the documentation of this file.
1*> \brief \b DDRVAB
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 DDRVAB( DOTYPE, NM, MVAL, NNS,
12* NSVAL, THRESH, NMAX, A, AFAC, B,
13* X, WORK, RWORK, SWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* INTEGER NM, NMAX, NNS, NOUT
17* DOUBLE PRECISION THRESH
18* ..
19* .. Array Arguments ..
20* LOGICAL DOTYPE( * )
21* INTEGER MVAL( * ), NSVAL( * ), IWORK( * )
22* REAL SWORK(*)
23* DOUBLE PRECISION A( * ), AFAC( * ), B( * ),
24* $ RWORK( * ), WORK( * ), X( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> DDRVAB tests DSGESV
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] DOTYPE
40*> \verbatim
41*> DOTYPE is LOGICAL array, dimension (NTYPES)
42*> The matrix types to be used for testing. Matrices of type j
43*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
44*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
45*> \endverbatim
46*>
47*> \param[in] NM
48*> \verbatim
49*> NM is INTEGER
50*> The number of values of M contained in the vector MVAL.
51*> \endverbatim
52*>
53*> \param[in] MVAL
54*> \verbatim
55*> MVAL is INTEGER array, dimension (NM)
56*> The values of the matrix row dimension M.
57*> \endverbatim
58*>
59*> \param[in] NNS
60*> \verbatim
61*> NNS is INTEGER
62*> The number of values of NRHS contained in the vector NSVAL.
63*> \endverbatim
64*>
65*> \param[in] NSVAL
66*> \verbatim
67*> NSVAL is INTEGER array, dimension (NNS)
68*> The values of the number of right hand sides NRHS.
69*> \endverbatim
70*>
71*> \param[in] THRESH
72*> \verbatim
73*> THRESH is DOUBLE PRECISION
74*> The threshold value for the test ratios. A result is
75*> included in the output file if RESULT >= THRESH. To have
76*> every test ratio printed, use THRESH = 0.
77*> \endverbatim
78*>
79*> \param[in] NMAX
80*> \verbatim
81*> NMAX is INTEGER
82*> The maximum value permitted for M or N, used in dimensioning
83*> the work arrays.
84*> \endverbatim
85*>
86*> \param[out] A
87*> \verbatim
88*> A is DOUBLE PRECISION array, dimension (NMAX*NMAX)
89*> \endverbatim
90*>
91*> \param[out] AFAC
92*> \verbatim
93*> AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
94*> \endverbatim
95*>
96*> \param[out] B
97*> \verbatim
98*> B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
99*> where NSMAX is the largest entry in NSVAL.
100*> \endverbatim
101*>
102*> \param[out] X
103*> \verbatim
104*> X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
105*> \endverbatim
106*>
107*> \param[out] WORK
108*> \verbatim
109*> WORK is DOUBLE PRECISION array, dimension
110*> (NMAX*max(3,NSMAX))
111*> \endverbatim
112*>
113*> \param[out] RWORK
114*> \verbatim
115*> RWORK is DOUBLE PRECISION array, dimension
116*> (max(2*NMAX,2*NSMAX+NWORK))
117*> \endverbatim
118*>
119*> \param[out] SWORK
120*> \verbatim
121*> SWORK is REAL array, dimension
122*> (NMAX*(NSMAX+NMAX))
123*> \endverbatim
124*>
125*> \param[out] IWORK
126*> \verbatim
127*> IWORK is INTEGER array, dimension
128*> NMAX
129*> \endverbatim
130*>
131*> \param[in] NOUT
132*> \verbatim
133*> NOUT is INTEGER
134*> The unit number for output.
135*> \endverbatim
136*
137* Authors:
138* ========
139*
140*> \author Univ. of Tennessee
141*> \author Univ. of California Berkeley
142*> \author Univ. of Colorado Denver
143*> \author NAG Ltd.
144*
145*> \ingroup double_lin
146*
147* =====================================================================
148 SUBROUTINE ddrvab( DOTYPE, NM, MVAL, NNS,
149 $ NSVAL, THRESH, NMAX, A, AFAC, B,
150 $ X, WORK, RWORK, SWORK, IWORK, NOUT )
151*
152* -- LAPACK test routine --
153* -- LAPACK is a software package provided by Univ. of Tennessee, --
154* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155*
156* .. Scalar Arguments ..
157 INTEGER NM, NMAX, NNS, NOUT
158 DOUBLE PRECISION THRESH
159* ..
160* .. Array Arguments ..
161 LOGICAL DOTYPE( * )
162 INTEGER MVAL( * ), NSVAL( * ), IWORK( * )
163 REAL SWORK(*)
164 DOUBLE PRECISION A( * ), AFAC( * ), B( * ),
165 $ rwork( * ), work( * ), x( * )
166* ..
167*
168* =====================================================================
169*
170* .. Parameters ..
171 DOUBLE PRECISION ZERO
172 PARAMETER ( ZERO = 0.0d+0 )
173 INTEGER NTYPES
174 parameter( ntypes = 11 )
175 INTEGER NTESTS
176 parameter( ntests = 1 )
177* ..
178* .. Local Scalars ..
179 LOGICAL ZEROT
180 CHARACTER DIST, TRANS, TYPE, XTYPE
181 CHARACTER*3 PATH
182 INTEGER I, IM, IMAT, INFO, IOFF, IRHS,
183 $ izero, kl, ku, lda, m, mode, n,
184 $ nerrs, nfail, nimat, nrhs, nrun
185 DOUBLE PRECISION ANORM, CNDNUM
186* ..
187* .. Local Arrays ..
188 INTEGER ISEED( 4 ), ISEEDY( 4 )
189 DOUBLE PRECISION RESULT( NTESTS )
190* ..
191* .. Local Variables ..
192 INTEGER ITER, KASE
193* ..
194* .. External Subroutines ..
195 EXTERNAL alaerh, alahd, dget08, dlacpy, dlarhs, dlaset,
196 $ dlatb4, dlatms
197* ..
198* .. Intrinsic Functions ..
199 INTRINSIC dble, max, min, sqrt
200* ..
201* .. Scalars in Common ..
202 LOGICAL LERR, OK
203 CHARACTER*32 SRNAMT
204 INTEGER INFOT, NUNIT
205* ..
206* .. Common blocks ..
207 COMMON / infoc / infot, nunit, ok, lerr
208 COMMON / srnamc / srnamt
209* ..
210* .. Data statements ..
211 DATA iseedy / 2006, 2007, 2008, 2009 /
212* ..
213* .. Executable Statements ..
214*
215* Initialize constants and the random number seed.
216*
217 kase = 0
218 path( 1: 1 ) = 'Double precision'
219 path( 2: 3 ) = 'GE'
220 nrun = 0
221 nfail = 0
222 nerrs = 0
223 DO 10 i = 1, 4
224 iseed( i ) = iseedy( i )
225 10 CONTINUE
226*
227 infot = 0
228*
229* Do for each value of M in MVAL
230*
231 DO 120 im = 1, nm
232 m = mval( im )
233 lda = max( 1, m )
234*
235 n = m
236 nimat = ntypes
237 IF( m.LE.0 .OR. n.LE.0 )
238 $ nimat = 1
239*
240 DO 100 imat = 1, nimat
241*
242* Do the tests only if DOTYPE( IMAT ) is true.
243*
244 IF( .NOT.dotype( imat ) )
245 $ GO TO 100
246*
247* Skip types 5, 6, or 7 if the matrix size is too small.
248*
249 zerot = imat.GE.5 .AND. imat.LE.7
250 IF( zerot .AND. n.LT.imat-4 )
251 $ GO TO 100
252*
253* Set up parameters with DLATB4 and generate a test matrix
254* with DLATMS.
255*
256 CALL dlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
257 $ cndnum, dist )
258*
259 srnamt = 'DLATMS'
260 CALL dlatms( m, n, dist, iseed, TYPE, rwork, mode,
261 $ cndnum, anorm, kl, ku, 'No packing', a, lda,
262 $ work, info )
263*
264* Check error code from DLATMS.
265*
266 IF( info.NE.0 ) THEN
267 CALL alaerh( path, 'DLATMS', info, 0, ' ', m, n, -1,
268 $ -1, -1, imat, nfail, nerrs, nout )
269 GO TO 100
270 END IF
271*
272* For types 5-7, zero one or more columns of the matrix to
273* test that INFO is returned correctly.
274*
275 IF( zerot ) THEN
276 IF( imat.EQ.5 ) THEN
277 izero = 1
278 ELSE IF( imat.EQ.6 ) THEN
279 izero = min( m, n )
280 ELSE
281 izero = min( m, n ) / 2 + 1
282 END IF
283 ioff = ( izero-1 )*lda
284 IF( imat.LT.7 ) THEN
285 DO 20 i = 1, m
286 a( ioff+i ) = zero
287 20 CONTINUE
288 ELSE
289 CALL dlaset( 'Full', m, n-izero+1, zero, zero,
290 $ a( ioff+1 ), lda )
291 END IF
292 ELSE
293 izero = 0
294 END IF
295*
296 DO 60 irhs = 1, nns
297 nrhs = nsval( irhs )
298 xtype = 'N'
299 trans = 'N'
300*
301 srnamt = 'DLARHS'
302 CALL dlarhs( path, xtype, ' ', trans, n, n, kl,
303 $ ku, nrhs, a, lda, x, lda, b,
304 $ lda, iseed, info )
305*
306 srnamt = 'DSGESV'
307*
308 kase = kase + 1
309*
310 CALL dlacpy( 'Full', m, n, a, lda, afac, lda )
311*
312 CALL dsgesv( n, nrhs, a, lda, iwork, b, lda, x, lda,
313 $ work, swork, iter, info)
314*
315 IF (iter.LT.0) THEN
316 CALL dlacpy( 'Full', m, n, afac, lda, a, lda )
317 ENDIF
318*
319* Check error code from DSGESV. This should be the same as
320* the one of DGETRF.
321*
322 IF( info.NE.izero ) THEN
323*
324 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
325 $ CALL alahd( nout, path )
326 nerrs = nerrs + 1
327*
328 IF( info.NE.izero .AND. izero.NE.0 ) THEN
329 WRITE( nout, fmt = 9988 )'DSGESV',info,
330 $ izero,m,imat
331 ELSE
332 WRITE( nout, fmt = 9975 )'DSGESV',info,
333 $ m, imat
334 END IF
335 END IF
336*
337* Skip the remaining test if the matrix is singular.
338*
339 IF( info.NE.0 )
340 $ GO TO 100
341*
342* Check the quality of the solution
343*
344 CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
345*
346 CALL dget08( trans, n, n, nrhs, a, lda, x, lda, work,
347 $ lda, rwork, result( 1 ) )
348*
349* Check if the test passes the testing.
350* Print information about the tests that did not
351* pass the testing.
352*
353* If iterative refinement has been used and claimed to
354* be successful (ITER>0), we want
355* NORMI(B - A*X)/(NORMI(A)*NORMI(X)*EPS*SRQT(N)) < 1
356*
357* If double precision has been used (ITER<0), we want
358* NORMI(B - A*X)/(NORMI(A)*NORMI(X)*EPS) < THRES
359* (Cf. the linear solver testing routines)
360*
361 IF ((thresh.LE.0.0e+00)
362 $ .OR.((iter.GE.0).AND.(n.GT.0)
363 $ .AND.(result(1).GE.sqrt(dble(n))))
364 $ .OR.((iter.LT.0).AND.(result(1).GE.thresh))) THEN
365*
366 IF( nfail.EQ.0 .AND. nerrs.EQ.0 ) THEN
367 WRITE( nout, fmt = 8999 )'DGE'
368 WRITE( nout, fmt = '( '' Matrix types:'' )' )
369 WRITE( nout, fmt = 8979 )
370 WRITE( nout, fmt = '( '' Test ratios:'' )' )
371 WRITE( nout, fmt = 8960 )1
372 WRITE( nout, fmt = '( '' Messages:'' )' )
373 END IF
374*
375 WRITE( nout, fmt = 9998 )trans, n, nrhs,
376 $ imat, 1, result( 1 )
377 nfail = nfail + 1
378 END IF
379 nrun = nrun + 1
380 60 CONTINUE
381 100 CONTINUE
382 120 CONTINUE
383*
384* Print a summary of the results.
385*
386 IF( nfail.GT.0 ) THEN
387 WRITE( nout, fmt = 9996 )'DSGESV', nfail, nrun
388 ELSE
389 WRITE( nout, fmt = 9995 )'DSGESV', nrun
390 END IF
391 IF( nerrs.GT.0 ) THEN
392 WRITE( nout, fmt = 9994 )nerrs
393 END IF
394*
395 9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
396 $ i2, ', test(', i2, ') =', g12.5 )
397 9996 FORMAT( 1x, a6, ': ', i6, ' out of ', i6,
398 $ ' tests failed to pass the threshold' )
399 9995 FORMAT( /1x, 'All tests for ', a6,
400 $ ' routines passed the threshold ( ', i6, ' tests run)' )
401 9994 FORMAT( 6x, i6, ' error messages recorded' )
402*
403* SUBNAM, INFO, INFOE, M, IMAT
404*
405 9988 FORMAT( ' *** ', a6, ' returned with INFO =', i5, ' instead of ',
406 $ i5, / ' ==> M =', i5, ', type ',
407 $ i2 )
408*
409* SUBNAM, INFO, M, IMAT
410*
411 9975 FORMAT( ' *** Error code from ', a6, '=', i5, ' for M=', i5,
412 $ ', type ', i2 )
413 8999 FORMAT( / 1x, a3, ': General dense matrices' )
414 8979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
415 $ '2. Upper triangular', 16x,
416 $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
417 $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
418 $ / 4x, '4. Random, CNDNUM = 2', 13x,
419 $ '10. Scaled near underflow', / 4x, '5. First column zero',
420 $ 14x, '11. Scaled near overflow', / 4x,
421 $ '6. Last column zero' )
422 8960 FORMAT( 3x, i2, ': norm_1( B - A * X ) / ',
423 $ '( norm_1(A) * norm_1(X) * EPS * SQRT(N) ) > 1 if ITERREF',
424 $ / 4x, 'or norm_1( B - A * X ) / ',
425 $ '( norm_1(A) * norm_1(X) * EPS ) > THRES if DGETRF' )
426 RETURN
427*
428* End of DDRVAB
429*
430 END
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 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 ddrvab(dotype, nm, mval, nns, nsval, thresh, nmax, a, afac, b, x, work, rwork, swork, iwork, nout)
DDRVAB
Definition ddrvab.f:151
subroutine dget08(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DGET08
Definition dget08.f:133
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 dsgesv(n, nrhs, a, lda, ipiv, b, ldb, x, ldx, work, swork, iter, info)
DSGESV computes the solution to system of linear equations A * X = B for GE matrices (mixed precision...
Definition dsgesv.f:195
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