LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ddrvpp ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
double precision, dimension( * )  A,
double precision, dimension( * )  AFAC,
double precision, dimension( * )  ASAV,
double precision, dimension( * )  B,
double precision, dimension( * )  BSAV,
double precision, dimension( * )  X,
double precision, dimension( * )  XACT,
double precision, dimension( * )  S,
double precision, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

DDRVPP

Purpose:
 DDRVPP tests the driver routines DPPSV and -SVX.
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]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[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]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[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+1)/2)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension
                      (NMAX*(NMAX+1)/2)
[out]ASAV
          ASAV is DOUBLE PRECISION array, dimension
                      (NMAX*(NMAX+1)/2)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]S
          S is DOUBLE PRECISION array, dimension (NMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[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 169 of file ddrvpp.f.

169 *
170 * -- LAPACK test routine (version 3.4.0) --
171 * -- LAPACK is a software package provided by Univ. of Tennessee, --
172 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
173 * November 2011
174 *
175 * .. Scalar Arguments ..
176  LOGICAL tsterr
177  INTEGER nmax, nn, nout, nrhs
178  DOUBLE PRECISION thresh
179 * ..
180 * .. Array Arguments ..
181  LOGICAL dotype( * )
182  INTEGER iwork( * ), nval( * )
183  DOUBLE PRECISION a( * ), afac( * ), asav( * ), b( * ),
184  $ bsav( * ), rwork( * ), s( * ), work( * ),
185  $ x( * ), xact( * )
186 * ..
187 *
188 * =====================================================================
189 *
190 * .. Parameters ..
191  DOUBLE PRECISION one, zero
192  parameter ( one = 1.0d+0, zero = 0.0d+0 )
193  INTEGER ntypes
194  parameter ( ntypes = 9 )
195  INTEGER ntests
196  parameter ( ntests = 6 )
197 * ..
198 * .. Local Scalars ..
199  LOGICAL equil, nofact, prefac, zerot
200  CHARACTER dist, equed, fact, packit, TYPE, uplo, xtype
201  CHARACTER*3 path
202  INTEGER i, iequed, ifact, imat, in, info, ioff, iuplo,
203  $ izero, k, k1, kl, ku, lda, mode, n, nerrs,
204  $ nfact, nfail, nimat, npp, nrun, nt
205  DOUBLE PRECISION ainvnm, amax, anorm, cndnum, rcond, rcondc,
206  $ roldc, scond
207 * ..
208 * .. Local Arrays ..
209  CHARACTER equeds( 2 ), facts( 3 ), packs( 2 ), uplos( 2 )
210  INTEGER iseed( 4 ), iseedy( 4 )
211  DOUBLE PRECISION result( ntests )
212 * ..
213 * .. External Functions ..
214  LOGICAL lsame
215  DOUBLE PRECISION dget06, dlansp
216  EXTERNAL lsame, dget06, dlansp
217 * ..
218 * .. External Subroutines ..
219  EXTERNAL aladhd, alaerh, alasvm, dcopy, derrvx, dget04,
222  $ dpptrf, dpptri
223 * ..
224 * .. Scalars in Common ..
225  LOGICAL lerr, ok
226  CHARACTER*32 srnamt
227  INTEGER infot, nunit
228 * ..
229 * .. Common blocks ..
230  COMMON / infoc / infot, nunit, ok, lerr
231  COMMON / srnamc / srnamt
232 * ..
233 * .. Intrinsic Functions ..
234  INTRINSIC max
235 * ..
236 * .. Data statements ..
237  DATA iseedy / 1988, 1989, 1990, 1991 /
238  DATA uplos / 'U', 'L' / , facts / 'F', 'N', 'E' / ,
239  $ packs / 'C', 'R' / , equeds / 'N', 'Y' /
240 * ..
241 * .. Executable Statements ..
242 *
243 * Initialize constants and the random number seed.
244 *
245  path( 1: 1 ) = 'Double precision'
246  path( 2: 3 ) = 'PP'
247  nrun = 0
248  nfail = 0
249  nerrs = 0
250  DO 10 i = 1, 4
251  iseed( i ) = iseedy( i )
252  10 CONTINUE
253 *
254 * Test the error exits
255 *
256  IF( tsterr )
257  $ CALL derrvx( path, nout )
258  infot = 0
259 *
260 * Do for each value of N in NVAL
261 *
262  DO 140 in = 1, nn
263  n = nval( in )
264  lda = max( n, 1 )
265  npp = n*( n+1 ) / 2
266  xtype = 'N'
267  nimat = ntypes
268  IF( n.LE.0 )
269  $ nimat = 1
270 *
271  DO 130 imat = 1, nimat
272 *
273 * Do the tests only if DOTYPE( IMAT ) is true.
274 *
275  IF( .NOT.dotype( imat ) )
276  $ GO TO 130
277 *
278 * Skip types 3, 4, or 5 if the matrix size is too small.
279 *
280  zerot = imat.GE.3 .AND. imat.LE.5
281  IF( zerot .AND. n.LT.imat-2 )
282  $ GO TO 130
283 *
284 * Do first for UPLO = 'U', then for UPLO = 'L'
285 *
286  DO 120 iuplo = 1, 2
287  uplo = uplos( iuplo )
288  packit = packs( iuplo )
289 *
290 * Set up parameters with DLATB4 and generate a test matrix
291 * with DLATMS.
292 *
293  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
294  $ cndnum, dist )
295  rcondc = one / cndnum
296 *
297  srnamt = 'DLATMS'
298  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
299  $ cndnum, anorm, kl, ku, packit, a, lda, work,
300  $ info )
301 *
302 * Check error code from DLATMS.
303 *
304  IF( info.NE.0 ) THEN
305  CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
306  $ -1, -1, imat, nfail, nerrs, nout )
307  GO TO 120
308  END IF
309 *
310 * For types 3-5, zero one row and column of the matrix to
311 * test that INFO is returned correctly.
312 *
313  IF( zerot ) THEN
314  IF( imat.EQ.3 ) THEN
315  izero = 1
316  ELSE IF( imat.EQ.4 ) THEN
317  izero = n
318  ELSE
319  izero = n / 2 + 1
320  END IF
321 *
322 * Set row and column IZERO of A to 0.
323 *
324  IF( iuplo.EQ.1 ) THEN
325  ioff = ( izero-1 )*izero / 2
326  DO 20 i = 1, izero - 1
327  a( ioff+i ) = zero
328  20 CONTINUE
329  ioff = ioff + izero
330  DO 30 i = izero, n
331  a( ioff ) = zero
332  ioff = ioff + i
333  30 CONTINUE
334  ELSE
335  ioff = izero
336  DO 40 i = 1, izero - 1
337  a( ioff ) = zero
338  ioff = ioff + n - i
339  40 CONTINUE
340  ioff = ioff - izero
341  DO 50 i = izero, n
342  a( ioff+i ) = zero
343  50 CONTINUE
344  END IF
345  ELSE
346  izero = 0
347  END IF
348 *
349 * Save a copy of the matrix A in ASAV.
350 *
351  CALL dcopy( npp, a, 1, asav, 1 )
352 *
353  DO 110 iequed = 1, 2
354  equed = equeds( iequed )
355  IF( iequed.EQ.1 ) THEN
356  nfact = 3
357  ELSE
358  nfact = 1
359  END IF
360 *
361  DO 100 ifact = 1, nfact
362  fact = facts( ifact )
363  prefac = lsame( fact, 'F' )
364  nofact = lsame( fact, 'N' )
365  equil = lsame( fact, 'E' )
366 *
367  IF( zerot ) THEN
368  IF( prefac )
369  $ GO TO 100
370  rcondc = zero
371 *
372  ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
373 *
374 * Compute the condition number for comparison with
375 * the value returned by DPPSVX (FACT = 'N' reuses
376 * the condition number from the previous iteration
377 * with FACT = 'F').
378 *
379  CALL dcopy( npp, asav, 1, afac, 1 )
380  IF( equil .OR. iequed.GT.1 ) THEN
381 *
382 * Compute row and column scale factors to
383 * equilibrate the matrix A.
384 *
385  CALL dppequ( uplo, n, afac, s, scond, amax,
386  $ info )
387  IF( info.EQ.0 .AND. n.GT.0 ) THEN
388  IF( iequed.GT.1 )
389  $ scond = zero
390 *
391 * Equilibrate the matrix.
392 *
393  CALL dlaqsp( uplo, n, afac, s, scond,
394  $ amax, equed )
395  END IF
396  END IF
397 *
398 * Save the condition number of the
399 * non-equilibrated system for use in DGET04.
400 *
401  IF( equil )
402  $ roldc = rcondc
403 *
404 * Compute the 1-norm of A.
405 *
406  anorm = dlansp( '1', uplo, n, afac, rwork )
407 *
408 * Factor the matrix A.
409 *
410  CALL dpptrf( uplo, n, afac, info )
411 *
412 * Form the inverse of A.
413 *
414  CALL dcopy( npp, afac, 1, a, 1 )
415  CALL dpptri( uplo, n, a, info )
416 *
417 * Compute the 1-norm condition number of A.
418 *
419  ainvnm = dlansp( '1', uplo, n, a, rwork )
420  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
421  rcondc = one
422  ELSE
423  rcondc = ( one / anorm ) / ainvnm
424  END IF
425  END IF
426 *
427 * Restore the matrix A.
428 *
429  CALL dcopy( npp, asav, 1, a, 1 )
430 *
431 * Form an exact solution and set the right hand side.
432 *
433  srnamt = 'DLARHS'
434  CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
435  $ nrhs, a, lda, xact, lda, b, lda,
436  $ iseed, info )
437  xtype = 'C'
438  CALL dlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
439 *
440  IF( nofact ) THEN
441 *
442 * --- Test DPPSV ---
443 *
444 * Compute the L*L' or U'*U factorization of the
445 * matrix and solve the system.
446 *
447  CALL dcopy( npp, a, 1, afac, 1 )
448  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
449 *
450  srnamt = 'DPPSV '
451  CALL dppsv( uplo, n, nrhs, afac, x, lda, info )
452 *
453 * Check error code from DPPSV .
454 *
455  IF( info.NE.izero ) THEN
456  CALL alaerh( path, 'DPPSV ', info, izero,
457  $ uplo, n, n, -1, -1, nrhs, imat,
458  $ nfail, nerrs, nout )
459  GO TO 70
460  ELSE IF( info.NE.0 ) THEN
461  GO TO 70
462  END IF
463 *
464 * Reconstruct matrix from factors and compute
465 * residual.
466 *
467  CALL dppt01( uplo, n, a, afac, rwork,
468  $ result( 1 ) )
469 *
470 * Compute residual of the computed solution.
471 *
472  CALL dlacpy( 'Full', n, nrhs, b, lda, work,
473  $ lda )
474  CALL dppt02( uplo, n, nrhs, a, x, lda, work,
475  $ lda, rwork, result( 2 ) )
476 *
477 * Check solution from generated exact solution.
478 *
479  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
480  $ result( 3 ) )
481  nt = 3
482 *
483 * Print information about the tests that did not
484 * pass the threshold.
485 *
486  DO 60 k = 1, nt
487  IF( result( k ).GE.thresh ) THEN
488  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
489  $ CALL aladhd( nout, path )
490  WRITE( nout, fmt = 9999 )'DPPSV ', uplo,
491  $ n, imat, k, result( k )
492  nfail = nfail + 1
493  END IF
494  60 CONTINUE
495  nrun = nrun + nt
496  70 CONTINUE
497  END IF
498 *
499 * --- Test DPPSVX ---
500 *
501  IF( .NOT.prefac .AND. npp.GT.0 )
502  $ CALL dlaset( 'Full', npp, 1, zero, zero, afac,
503  $ npp )
504  CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
505  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
506 *
507 * Equilibrate the matrix if FACT='F' and
508 * EQUED='Y'.
509 *
510  CALL dlaqsp( uplo, n, a, s, scond, amax, equed )
511  END IF
512 *
513 * Solve the system and compute the condition number
514 * and error bounds using DPPSVX.
515 *
516  srnamt = 'DPPSVX'
517  CALL dppsvx( fact, uplo, n, nrhs, a, afac, equed,
518  $ s, b, lda, x, lda, rcond, rwork,
519  $ rwork( nrhs+1 ), work, iwork, info )
520 *
521 * Check the error code from DPPSVX.
522 *
523  IF( info.NE.izero ) THEN
524  CALL alaerh( path, 'DPPSVX', info, izero,
525  $ fact // uplo, n, n, -1, -1, nrhs,
526  $ imat, nfail, nerrs, nout )
527  GO TO 90
528  END IF
529 *
530  IF( info.EQ.0 ) THEN
531  IF( .NOT.prefac ) THEN
532 *
533 * Reconstruct matrix from factors and compute
534 * residual.
535 *
536  CALL dppt01( uplo, n, a, afac,
537  $ rwork( 2*nrhs+1 ), result( 1 ) )
538  k1 = 1
539  ELSE
540  k1 = 2
541  END IF
542 *
543 * Compute residual of the computed solution.
544 *
545  CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
546  $ lda )
547  CALL dppt02( uplo, n, nrhs, asav, x, lda, work,
548  $ lda, rwork( 2*nrhs+1 ),
549  $ result( 2 ) )
550 *
551 * Check solution from generated exact solution.
552 *
553  IF( nofact .OR. ( prefac .AND. lsame( equed,
554  $ 'N' ) ) ) THEN
555  CALL dget04( n, nrhs, x, lda, xact, lda,
556  $ rcondc, result( 3 ) )
557  ELSE
558  CALL dget04( n, nrhs, x, lda, xact, lda,
559  $ roldc, result( 3 ) )
560  END IF
561 *
562 * Check the error bounds from iterative
563 * refinement.
564 *
565  CALL dppt05( uplo, n, nrhs, asav, b, lda, x,
566  $ lda, xact, lda, rwork,
567  $ rwork( nrhs+1 ), result( 4 ) )
568  ELSE
569  k1 = 6
570  END IF
571 *
572 * Compare RCOND from DPPSVX with the computed value
573 * in RCONDC.
574 *
575  result( 6 ) = dget06( rcond, rcondc )
576 *
577 * Print information about the tests that did not pass
578 * the threshold.
579 *
580  DO 80 k = k1, 6
581  IF( result( k ).GE.thresh ) THEN
582  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
583  $ CALL aladhd( nout, path )
584  IF( prefac ) THEN
585  WRITE( nout, fmt = 9997 )'DPPSVX', fact,
586  $ uplo, n, equed, imat, k, result( k )
587  ELSE
588  WRITE( nout, fmt = 9998 )'DPPSVX', fact,
589  $ uplo, n, imat, k, result( k )
590  END IF
591  nfail = nfail + 1
592  END IF
593  80 CONTINUE
594  nrun = nrun + 7 - k1
595  90 CONTINUE
596  100 CONTINUE
597  110 CONTINUE
598  120 CONTINUE
599  130 CONTINUE
600  140 CONTINUE
601 *
602 * Print a summary of the results.
603 *
604  CALL alasvm( path, nout, nfail, nrun, nerrs )
605 *
606  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
607  $ ', test(', i1, ')=', g12.5 )
608  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
609  $ ', type ', i1, ', test(', i1, ')=', g12.5 )
610  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
611  $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ')=',
612  $ g12.5 )
613  RETURN
614 *
615 * End of DDRVPP
616 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
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 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 dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
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 dppt02(UPLO, N, NRHS, A, X, LDX, B, LDB, RWORK, RESID)
DPPT02
Definition: dppt02.f:124
subroutine dpptrf(UPLO, N, AP, INFO)
DPPTRF
Definition: dpptrf.f:121
subroutine dppequ(UPLO, N, AP, S, SCOND, AMAX, INFO)
DPPEQU
Definition: dppequ.f:118
subroutine dppsv(UPLO, N, NRHS, AP, B, LDB, INFO)
DPPSV computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: dppsv.f:146
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine dppsvx(FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
DPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: dppsvx.f:314
subroutine derrvx(PATH, NUNIT)
DERRVX
Definition: derrvx.f:57
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dppt05(UPLO, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DPPT05
Definition: dppt05.f:158
subroutine dpptri(UPLO, N, AP, INFO)
DPPTRI
Definition: dpptri.f:95
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dppt01(UPLO, N, A, AFAC, RWORK, RESID)
DPPT01
Definition: dppt01.f:95
subroutine dlaqsp(UPLO, N, AP, S, SCOND, AMAX, EQUED)
DLAQSP scales a symmetric/Hermitian matrix in packed storage, using scaling factors computed by sppeq...
Definition: dlaqsp.f:127
double precision function dlansp(NORM, UPLO, N, AP, WORK)
DLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric matrix supplied in packed form.
Definition: dlansp.f:116
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: