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

◆ ddrvab()

subroutine ddrvab ( logical, dimension( * )  dotype,
integer  nm,
integer, dimension( * )  mval,
integer  nns,
integer, dimension( * )  nsval,
double precision  thresh,
integer  nmax,
double precision, dimension( * )  a,
double precision, dimension( * )  afac,
double precision, dimension( * )  b,
double precision, dimension( * )  x,
double precision, dimension( * )  work,
double precision, dimension( * )  rwork,
real, dimension(*)  swork,
integer, dimension( * )  iwork,
integer  nout 
)

DDRVAB

Purpose:
 DDRVAB tests DSGESV
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]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[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]NMAX
          NMAX is INTEGER
          The maximum value permitted for M or 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]B
          B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(2*NMAX,2*NSMAX+NWORK))
[out]SWORK
          SWORK is REAL array, dimension
                      (NMAX*(NSMAX+NMAX))
[out]IWORK
          IWORK is INTEGER array, dimension
                      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 148 of file ddrvab.f.

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*
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 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
Here is the call graph for this function:
Here is the caller graph for this function: