LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2011

Definition at line 153 of file ddrvab.f.

153 *
154 * -- LAPACK test routine (version 3.4.0) --
155 * -- LAPACK is a software package provided by Univ. of Tennessee, --
156 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
157 * November 2011
158 *
159 * .. Scalar Arguments ..
160  INTEGER nm, nmax, nns, nout
161  DOUBLE PRECISION thresh
162 * ..
163 * .. Array Arguments ..
164  LOGICAL dotype( * )
165  INTEGER mval( * ), nsval( * ), iwork( * )
166  REAL swork(*)
167  DOUBLE PRECISION a( * ), afac( * ), b( * ),
168  $ rwork( * ), work( * ), x( * )
169 * ..
170 *
171 * =====================================================================
172 *
173 * .. Parameters ..
174  DOUBLE PRECISION zero
175  parameter ( zero = 0.0d+0 )
176  INTEGER ntypes
177  parameter ( ntypes = 11 )
178  INTEGER ntests
179  parameter ( ntests = 1 )
180 * ..
181 * .. Local Scalars ..
182  LOGICAL zerot
183  CHARACTER dist, trans, TYPE, xtype
184  CHARACTER*3 path
185  INTEGER i, im, imat, info, ioff, irhs,
186  $ izero, kl, ku, lda, m, mode, n,
187  $ nerrs, nfail, nimat, nrhs, nrun
188  DOUBLE PRECISION anorm, cndnum
189 * ..
190 * .. Local Arrays ..
191  INTEGER iseed( 4 ), iseedy( 4 )
192  DOUBLE PRECISION result( ntests )
193 * ..
194 * .. Local Variables ..
195  INTEGER iter, kase
196 * ..
197 * .. External Subroutines ..
198  EXTERNAL alaerh, alahd, dget08, dlacpy, dlarhs, dlaset,
199  $ dlatb4, dlatms
200 * ..
201 * .. Intrinsic Functions ..
202  INTRINSIC dble, max, min, sqrt
203 * ..
204 * .. Scalars in Common ..
205  LOGICAL lerr, ok
206  CHARACTER*32 srnamt
207  INTEGER infot, nunit
208 * ..
209 * .. Common blocks ..
210  COMMON / infoc / infot, nunit, ok, lerr
211  COMMON / srnamc / srnamt
212 * ..
213 * .. Data statements ..
214  DATA iseedy / 2006, 2007, 2008, 2009 /
215 * ..
216 * .. Executable Statements ..
217 *
218 * Initialize constants and the random number seed.
219 *
220  kase = 0
221  path( 1: 1 ) = 'Double precision'
222  path( 2: 3 ) = 'GE'
223  nrun = 0
224  nfail = 0
225  nerrs = 0
226  DO 10 i = 1, 4
227  iseed( i ) = iseedy( i )
228  10 CONTINUE
229 *
230  infot = 0
231 *
232 * Do for each value of M in MVAL
233 *
234  DO 120 im = 1, nm
235  m = mval( im )
236  lda = max( 1, m )
237 *
238  n = m
239  nimat = ntypes
240  IF( m.LE.0 .OR. n.LE.0 )
241  $ nimat = 1
242 *
243  DO 100 imat = 1, nimat
244 *
245 * Do the tests only if DOTYPE( IMAT ) is true.
246 *
247  IF( .NOT.dotype( imat ) )
248  $ GO TO 100
249 *
250 * Skip types 5, 6, or 7 if the matrix size is too small.
251 *
252  zerot = imat.GE.5 .AND. imat.LE.7
253  IF( zerot .AND. n.LT.imat-4 )
254  $ GO TO 100
255 *
256 * Set up parameters with DLATB4 and generate a test matrix
257 * with DLATMS.
258 *
259  CALL dlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
260  $ cndnum, dist )
261 *
262  srnamt = 'DLATMS'
263  CALL dlatms( m, n, dist, iseed, TYPE, rwork, mode,
264  $ cndnum, anorm, kl, ku, 'No packing', a, lda,
265  $ work, info )
266 *
267 * Check error code from DLATMS.
268 *
269  IF( info.NE.0 ) THEN
270  CALL alaerh( path, 'DLATMS', info, 0, ' ', m, n, -1,
271  $ -1, -1, imat, nfail, nerrs, nout )
272  GO TO 100
273  END IF
274 *
275 * For types 5-7, zero one or more columns of the matrix to
276 * test that INFO is returned correctly.
277 *
278  IF( zerot ) THEN
279  IF( imat.EQ.5 ) THEN
280  izero = 1
281  ELSE IF( imat.EQ.6 ) THEN
282  izero = min( m, n )
283  ELSE
284  izero = min( m, n ) / 2 + 1
285  END IF
286  ioff = ( izero-1 )*lda
287  IF( imat.LT.7 ) THEN
288  DO 20 i = 1, m
289  a( ioff+i ) = zero
290  20 CONTINUE
291  ELSE
292  CALL dlaset( 'Full', m, n-izero+1, zero, zero,
293  $ a( ioff+1 ), lda )
294  END IF
295  ELSE
296  izero = 0
297  END IF
298 *
299  DO 60 irhs = 1, nns
300  nrhs = nsval( irhs )
301  xtype = 'N'
302  trans = 'N'
303 *
304  srnamt = 'DLARHS'
305  CALL dlarhs( path, xtype, ' ', trans, n, n, kl,
306  $ ku, nrhs, a, lda, x, lda, b,
307  $ lda, iseed, info )
308 *
309  srnamt = 'DSGESV'
310 *
311  kase = kase + 1
312 *
313  CALL dlacpy( 'Full', m, n, a, lda, afac, lda )
314 *
315  CALL dsgesv( n, nrhs, a, lda, iwork, b, lda, x, lda,
316  $ work, swork, iter, info)
317 *
318  IF (iter.LT.0) THEN
319  CALL dlacpy( 'Full', m, n, afac, lda, a, lda )
320  ENDIF
321 *
322 * Check error code from DSGESV. This should be the same as
323 * the one of DGETRF.
324 *
325  IF( info.NE.izero ) THEN
326 *
327  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
328  $ CALL alahd( nout, path )
329  nerrs = nerrs + 1
330 *
331  IF( info.NE.izero .AND. izero.NE.0 ) THEN
332  WRITE( nout, fmt = 9988 )'DSGESV',info,
333  $ izero,m,imat
334  ELSE
335  WRITE( nout, fmt = 9975 )'DSGESV',info,
336  $ m, imat
337  END IF
338  END IF
339 *
340 * Skip the remaining test if the matrix is singular.
341 *
342  IF( info.NE.0 )
343  $ GO TO 100
344 *
345 * Check the quality of the solution
346 *
347  CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
348 *
349  CALL dget08( trans, n, n, nrhs, a, lda, x, lda, work,
350  $ lda, rwork, result( 1 ) )
351 *
352 * Check if the test passes the tesing.
353 * Print information about the tests that did not
354 * pass the testing.
355 *
356 * If iterative refinement has been used and claimed to
357 * be successful (ITER>0), we want
358 * NORMI(B - A*X)/(NORMI(A)*NORMI(X)*EPS*SRQT(N)) < 1
359 *
360 * If double precision has been used (ITER<0), we want
361 * NORMI(B - A*X)/(NORMI(A)*NORMI(X)*EPS) < THRES
362 * (Cf. the linear solver testing routines)
363 *
364  IF ((thresh.LE.0.0e+00)
365  $ .OR.((iter.GE.0).AND.(n.GT.0)
366  $ .AND.(result(1).GE.sqrt(dble(n))))
367  $ .OR.((iter.LT.0).AND.(result(1).GE.thresh))) THEN
368 *
369  IF( nfail.EQ.0 .AND. nerrs.EQ.0 ) THEN
370  WRITE( nout, fmt = 8999 )'DGE'
371  WRITE( nout, fmt = '( '' Matrix types:'' )' )
372  WRITE( nout, fmt = 8979 )
373  WRITE( nout, fmt = '( '' Test ratios:'' )' )
374  WRITE( nout, fmt = 8960 )1
375  WRITE( nout, fmt = '( '' Messages:'' )' )
376  END IF
377 *
378  WRITE( nout, fmt = 9998 )trans, n, nrhs,
379  $ imat, 1, result( 1 )
380  nfail = nfail + 1
381  END IF
382  nrun = nrun + 1
383  60 CONTINUE
384  100 CONTINUE
385  120 CONTINUE
386 *
387 * Print a summary of the results.
388 *
389  IF( nfail.GT.0 ) THEN
390  WRITE( nout, fmt = 9996 )'DSGESV', nfail, nrun
391  ELSE
392  WRITE( nout, fmt = 9995 )'DSGESV', nrun
393  END IF
394  IF( nerrs.GT.0 ) THEN
395  WRITE( nout, fmt = 9994 )nerrs
396  END IF
397 *
398  9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
399  $ i2, ', test(', i2, ') =', g12.5 )
400  9996 FORMAT( 1x, a6, ': ', i6, ' out of ', i6,
401  $ ' tests failed to pass the threshold' )
402  9995 FORMAT( /1x, 'All tests for ', a6,
403  $ ' routines passed the threshold ( ', i6, ' tests run)' )
404  9994 FORMAT( 6x, i6, ' error messages recorded' )
405 *
406 * SUBNAM, INFO, INFOE, M, IMAT
407 *
408  9988 FORMAT( ' *** ', a6, ' returned with INFO =', i5, ' instead of ',
409  $ i5, / ' ==> M =', i5, ', type ',
410  $ i2 )
411 *
412 * SUBNAM, INFO, M, IMAT
413 *
414  9975 FORMAT( ' *** Error code from ', a6, '=', i5, ' for M=', i5,
415  $ ', type ', i2 )
416  8999 FORMAT( / 1x, a3, ': General dense matrices' )
417  8979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
418  $ '2. Upper triangular', 16x,
419  $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
420  $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
421  $ / 4x, '4. Random, CNDNUM = 2', 13x,
422  $ '10. Scaled near underflow', / 4x, '5. First column zero',
423  $ 14x, '11. Scaled near overflow', / 4x,
424  $ '6. Last column zero' )
425  8960 FORMAT( 3x, i2, ': norm_1( B - A * X ) / ',
426  $ '( norm_1(A) * norm_1(X) * EPS * SQRT(N) ) > 1 if ITERREF',
427  $ / 4x, 'or norm_1( B - A * X ) / ',
428  $ '( norm_1(A) * norm_1(X) * EPS ) > THRES if DGETRF' )
429  RETURN
430 *
431 * End of DDRVAB
432 *
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:112
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
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 precisio...
Definition: dsgesv.f:197
subroutine dget08(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DGET08
Definition: dget08.f:135
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323

Here is the call graph for this function:

Here is the caller graph for this function: