LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ddrvac ( 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  NOUT 
)

DDRVAC

Purpose:
 DDRVAC tests DSPOSV.
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 N contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix dimension N.
[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 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)
[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))
[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 146 of file ddrvac.f.

146 *
147 * -- LAPACK test routine (version 3.4.0) --
148 * -- LAPACK is a software package provided by Univ. of Tennessee, --
149 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150 * November 2011
151 *
152 * .. Scalar Arguments ..
153  INTEGER nmax, nm, nns, nout
154  DOUBLE PRECISION thresh
155 * ..
156 * .. Array Arguments ..
157  LOGICAL dotype( * )
158  INTEGER mval( * ), nsval( * )
159  REAL swork(*)
160  DOUBLE PRECISION a( * ), afac( * ), b( * ),
161  $ rwork( * ), work( * ), x( * )
162 * ..
163 *
164 * =====================================================================
165 *
166 * .. Parameters ..
167  DOUBLE PRECISION zero
168  parameter ( zero = 0.0d+0 )
169  INTEGER ntypes
170  parameter ( ntypes = 9 )
171  INTEGER ntests
172  parameter ( ntests = 1 )
173 * ..
174 * .. Local Scalars ..
175  LOGICAL zerot
176  CHARACTER dist, TYPE, uplo, xtype
177  CHARACTER*3 path
178  INTEGER i, im, imat, info, ioff, irhs, iuplo,
179  $ izero, kl, ku, lda, mode, n,
180  $ nerrs, nfail, nimat, nrhs, nrun
181  DOUBLE PRECISION anorm, cndnum
182 * ..
183 * .. Local Arrays ..
184  CHARACTER uplos( 2 )
185  INTEGER iseed( 4 ), iseedy( 4 )
186  DOUBLE PRECISION result( ntests )
187 * ..
188 * .. Local Variables ..
189  INTEGER iter, kase
190 * ..
191 * .. External Functions ..
192  LOGICAL lsame
193  EXTERNAL lsame
194 * ..
195 * .. External Subroutines ..
196  EXTERNAL alaerh, dlacpy,
197  $ dlarhs, dlaset, dlatb4, dlatms,
198  $ dpot06, dsposv
199 * ..
200 * .. Intrinsic Functions ..
201  INTRINSIC dble, max, sqrt
202 * ..
203 * .. Scalars in Common ..
204  LOGICAL lerr, ok
205  CHARACTER*32 srnamt
206  INTEGER infot, nunit
207 * ..
208 * .. Common blocks ..
209  COMMON / infoc / infot, nunit, ok, lerr
210  COMMON / srnamc / srnamt
211 * ..
212 * .. Data statements ..
213  DATA iseedy / 1988, 1989, 1990, 1991 /
214  DATA uplos / 'U', 'L' /
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 ) = 'PO'
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 N in MVAL
233 *
234  DO 120 im = 1, nm
235  n = mval( im )
236  lda = max( n, 1 )
237  nimat = ntypes
238  IF( n.LE.0 )
239  $ nimat = 1
240 *
241  DO 110 imat = 1, nimat
242 *
243 * Do the tests only if DOTYPE( IMAT ) is true.
244 *
245  IF( .NOT.dotype( imat ) )
246  $ GO TO 110
247 *
248 * Skip types 3, 4, or 5 if the matrix size is too small.
249 *
250  zerot = imat.GE.3 .AND. imat.LE.5
251  IF( zerot .AND. n.LT.imat-2 )
252  $ GO TO 110
253 *
254 * Do first for UPLO = 'U', then for UPLO = 'L'
255 *
256  DO 100 iuplo = 1, 2
257  uplo = uplos( iuplo )
258 *
259 * Set up parameters with DLATB4 and generate a test matrix
260 * with DLATMS.
261 *
262  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
263  $ cndnum, dist )
264 *
265  srnamt = 'DLATMS'
266  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
267  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
268  $ info )
269 *
270 * Check error code from DLATMS.
271 *
272  IF( info.NE.0 ) THEN
273  CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
274  $ -1, -1, imat, nfail, nerrs, nout )
275  GO TO 100
276  END IF
277 *
278 * For types 3-5, zero one row and column of the matrix to
279 * test that INFO is returned correctly.
280 *
281  IF( zerot ) THEN
282  IF( imat.EQ.3 ) THEN
283  izero = 1
284  ELSE IF( imat.EQ.4 ) THEN
285  izero = n
286  ELSE
287  izero = n / 2 + 1
288  END IF
289  ioff = ( izero-1 )*lda
290 *
291 * Set row and column IZERO of A to 0.
292 *
293  IF( iuplo.EQ.1 ) THEN
294  DO 20 i = 1, izero - 1
295  a( ioff+i ) = zero
296  20 CONTINUE
297  ioff = ioff + izero
298  DO 30 i = izero, n
299  a( ioff ) = zero
300  ioff = ioff + lda
301  30 CONTINUE
302  ELSE
303  ioff = izero
304  DO 40 i = 1, izero - 1
305  a( ioff ) = zero
306  ioff = ioff + lda
307  40 CONTINUE
308  ioff = ioff - izero
309  DO 50 i = izero, n
310  a( ioff+i ) = zero
311  50 CONTINUE
312  END IF
313  ELSE
314  izero = 0
315  END IF
316 *
317  DO 60 irhs = 1, nns
318  nrhs = nsval( irhs )
319  xtype = 'N'
320 *
321 * Form an exact solution and set the right hand side.
322 *
323  srnamt = 'DLARHS'
324  CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
325  $ nrhs, a, lda, x, lda, b, lda,
326  $ iseed, info )
327 *
328 * Compute the L*L' or U'*U factorization of the
329 * matrix and solve the system.
330 *
331  srnamt = 'DSPOSV '
332  kase = kase + 1
333 *
334  CALL dlacpy( 'All', n, n, a, lda, afac, lda)
335 *
336  CALL dsposv( uplo, n, nrhs, afac, lda, b, lda, x, lda,
337  $ work, swork, iter, info )
338 
339  IF (iter.LT.0) THEN
340  CALL dlacpy( 'All', n, n, a, lda, afac, lda )
341  ENDIF
342 *
343 * Check error code from DSPOSV .
344 *
345  IF( info.NE.izero ) THEN
346 *
347  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
348  $ CALL alahd( nout, path )
349  nerrs = nerrs + 1
350 *
351  IF( info.NE.izero .AND. izero.NE.0 ) THEN
352  WRITE( nout, fmt = 9988 )'DSPOSV',info,izero,n,
353  $ imat
354  ELSE
355  WRITE( nout, fmt = 9975 )'DSPOSV',info,n,imat
356  END IF
357  END IF
358 *
359 * Skip the remaining test if the matrix is singular.
360 *
361  IF( info.NE.0 )
362  $ GO TO 110
363 *
364 * Check the quality of the solution
365 *
366  CALL dlacpy( 'All', n, nrhs, b, lda, work, lda )
367 *
368  CALL dpot06( uplo, n, nrhs, a, lda, x, lda, work,
369  $ lda, rwork, result( 1 ) )
370 *
371 * Check if the test passes the tesing.
372 * Print information about the tests that did not
373 * pass the testing.
374 *
375 * If iterative refinement has been used and claimed to
376 * be successful (ITER>0), we want
377 * NORM1(B - A*X)/(NORM1(A)*NORM1(X)*EPS*SRQT(N)) < 1
378 *
379 * If double precision has been used (ITER<0), we want
380 * NORM1(B - A*X)/(NORM1(A)*NORM1(X)*EPS) < THRES
381 * (Cf. the linear solver testing routines)
382 *
383  IF ((thresh.LE.0.0e+00)
384  $ .OR.((iter.GE.0).AND.(n.GT.0)
385  $ .AND.(result(1).GE.sqrt(dble(n))))
386  $ .OR.((iter.LT.0).AND.(result(1).GE.thresh))) THEN
387 *
388  IF( nfail.EQ.0 .AND. nerrs.EQ.0 ) THEN
389  WRITE( nout, fmt = 8999 )'DPO'
390  WRITE( nout, fmt = '( '' Matrix types:'' )' )
391  WRITE( nout, fmt = 8979 )
392  WRITE( nout, fmt = '( '' Test ratios:'' )' )
393  WRITE( nout, fmt = 8960 )1
394  WRITE( nout, fmt = '( '' Messages:'' )' )
395  END IF
396 *
397  WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat, 1,
398  $ result( 1 )
399 *
400  nfail = nfail + 1
401 *
402  END IF
403 *
404  nrun = nrun + 1
405 *
406  60 CONTINUE
407  100 CONTINUE
408  110 CONTINUE
409  120 CONTINUE
410 *
411 * Print a summary of the results.
412 *
413  IF( nfail.GT.0 ) THEN
414  WRITE( nout, fmt = 9996 )'DSPOSV', nfail, nrun
415  ELSE
416  WRITE( nout, fmt = 9995 )'DSPOSV', nrun
417  END IF
418  IF( nerrs.GT.0 ) THEN
419  WRITE( nout, fmt = 9994 )nerrs
420  END IF
421 *
422  9998 FORMAT( ' UPLO=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
423  $ i2, ', test(', i2, ') =', g12.5 )
424  9996 FORMAT( 1x, a6, ': ', i6, ' out of ', i6,
425  $ ' tests failed to pass the threshold' )
426  9995 FORMAT( /1x, 'All tests for ', a6,
427  $ ' routines passed the threshold ( ', i6, ' tests run)' )
428  9994 FORMAT( 6x, i6, ' error messages recorded' )
429 *
430 * SUBNAM, INFO, INFOE, N, IMAT
431 *
432  9988 FORMAT( ' *** ', a6, ' returned with INFO =', i5, ' instead of ',
433  $ i5, / ' ==> N =', i5, ', type ',
434  $ i2 )
435 *
436 * SUBNAM, INFO, N, IMAT
437 *
438  9975 FORMAT( ' *** Error code from ', a6, '=', i5, ' for M=', i5,
439  $ ', type ', i2 )
440  8999 FORMAT( / 1x, a3, ': positive definite dense matrices' )
441  8979 FORMAT( 4x, '1. Diagonal', 24x, '7. Last n/2 columns zero', / 4x,
442  $ '2. Upper triangular', 16x,
443  $ '8. Random, CNDNUM = sqrt(0.1/EPS)', / 4x,
444  $ '3. Lower triangular', 16x, '9. Random, CNDNUM = 0.1/EPS',
445  $ / 4x, '4. Random, CNDNUM = 2', 13x,
446  $ '10. Scaled near underflow', / 4x, '5. First column zero',
447  $ 14x, '11. Scaled near overflow', / 4x,
448  $ '6. Last column zero' )
449  8960 FORMAT( 3x, i2, ': norm_1( B - A * X ) / ',
450  $ '( norm_1(A) * norm_1(X) * EPS * SQRT(N) ) > 1 if ITERREF',
451  $ / 4x, 'or norm_1( B - A * X ) / ',
452  $ '( norm_1(A) * norm_1(X) * EPS ) > THRES if DPOTRF' )
453 
454  RETURN
455 *
456 * End of DDRVAC
457 *
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 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 dsposv(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, WORK, SWORK, ITER, INFO)
DSPOSV computes the solution to system of linear equations A * X = B for PO matrices ...
Definition: dsposv.f:201
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dpot06(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DPOT06
Definition: dpot06.f:129
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: