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

ZDRVPP

Purpose:
 ZDRVPP tests the driver routines ZPPSV 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 COMPLEX*16 array, dimension (NMAX*(NMAX+1)/2)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*(NMAX+1)/2)
[out]ASAV
          ASAV is COMPLEX*16 array, dimension (NMAX*(NMAX+1)/2)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]S
          S is DOUBLE PRECISION array, dimension (NMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[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 161 of file zdrvpp.f.

161 *
162 * -- LAPACK test routine (version 3.4.0) --
163 * -- LAPACK is a software package provided by Univ. of Tennessee, --
164 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
165 * November 2011
166 *
167 * .. Scalar Arguments ..
168  LOGICAL tsterr
169  INTEGER nmax, nn, nout, nrhs
170  DOUBLE PRECISION thresh
171 * ..
172 * .. Array Arguments ..
173  LOGICAL dotype( * )
174  INTEGER nval( * )
175  DOUBLE PRECISION rwork( * ), s( * )
176  COMPLEX*16 a( * ), afac( * ), asav( * ), b( * ),
177  $ bsav( * ), work( * ), x( * ), xact( * )
178 * ..
179 *
180 * =====================================================================
181 *
182 * .. Parameters ..
183  DOUBLE PRECISION one, zero
184  parameter ( one = 1.0d+0, zero = 0.0d+0 )
185  INTEGER ntypes
186  parameter ( ntypes = 9 )
187  INTEGER ntests
188  parameter ( ntests = 6 )
189 * ..
190 * .. Local Scalars ..
191  LOGICAL equil, nofact, prefac, zerot
192  CHARACTER dist, equed, fact, packit, TYPE, uplo, xtype
193  CHARACTER*3 path
194  INTEGER i, iequed, ifact, imat, in, info, ioff, iuplo,
195  $ izero, k, k1, kl, ku, lda, mode, n, nerrs,
196  $ nfact, nfail, nimat, npp, nrun, nt
197  DOUBLE PRECISION ainvnm, amax, anorm, cndnum, rcond, rcondc,
198  $ roldc, scond
199 * ..
200 * .. Local Arrays ..
201  CHARACTER equeds( 2 ), facts( 3 ), packs( 2 ), uplos( 2 )
202  INTEGER iseed( 4 ), iseedy( 4 )
203  DOUBLE PRECISION result( ntests )
204 * ..
205 * .. External Functions ..
206  LOGICAL lsame
207  DOUBLE PRECISION dget06, zlanhp
208  EXTERNAL lsame, dget06, zlanhp
209 * ..
210 * .. External Subroutines ..
211  EXTERNAL aladhd, alaerh, alasvm, zcopy, zerrvx, zget04,
214  $ zppt05, zpptrf, zpptri
215 * ..
216 * .. Scalars in Common ..
217  LOGICAL lerr, ok
218  CHARACTER*32 srnamt
219  INTEGER infot, nunit
220 * ..
221 * .. Common blocks ..
222  COMMON / infoc / infot, nunit, ok, lerr
223  COMMON / srnamc / srnamt
224 * ..
225 * .. Intrinsic Functions ..
226  INTRINSIC dcmplx, max
227 * ..
228 * .. Data statements ..
229  DATA iseedy / 1988, 1989, 1990, 1991 /
230  DATA uplos / 'U', 'L' / , facts / 'F', 'N', 'E' / ,
231  $ packs / 'C', 'R' / , equeds / 'N', 'Y' /
232 * ..
233 * .. Executable Statements ..
234 *
235 * Initialize constants and the random number seed.
236 *
237  path( 1: 1 ) = 'Zomplex precision'
238  path( 2: 3 ) = 'PP'
239  nrun = 0
240  nfail = 0
241  nerrs = 0
242  DO 10 i = 1, 4
243  iseed( i ) = iseedy( i )
244  10 CONTINUE
245 *
246 * Test the error exits
247 *
248  IF( tsterr )
249  $ CALL zerrvx( path, nout )
250  infot = 0
251 *
252 * Do for each value of N in NVAL
253 *
254  DO 140 in = 1, nn
255  n = nval( in )
256  lda = max( n, 1 )
257  npp = n*( n+1 ) / 2
258  xtype = 'N'
259  nimat = ntypes
260  IF( n.LE.0 )
261  $ nimat = 1
262 *
263  DO 130 imat = 1, nimat
264 *
265 * Do the tests only if DOTYPE( IMAT ) is true.
266 *
267  IF( .NOT.dotype( imat ) )
268  $ GO TO 130
269 *
270 * Skip types 3, 4, or 5 if the matrix size is too small.
271 *
272  zerot = imat.GE.3 .AND. imat.LE.5
273  IF( zerot .AND. n.LT.imat-2 )
274  $ GO TO 130
275 *
276 * Do first for UPLO = 'U', then for UPLO = 'L'
277 *
278  DO 120 iuplo = 1, 2
279  uplo = uplos( iuplo )
280  packit = packs( iuplo )
281 *
282 * Set up parameters with ZLATB4 and generate a test matrix
283 * with ZLATMS.
284 *
285  CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
286  $ cndnum, dist )
287  rcondc = one / cndnum
288 *
289  srnamt = 'ZLATMS'
290  CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
291  $ cndnum, anorm, kl, ku, packit, a, lda, work,
292  $ info )
293 *
294 * Check error code from ZLATMS.
295 *
296  IF( info.NE.0 ) THEN
297  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
298  $ -1, -1, imat, nfail, nerrs, nout )
299  GO TO 120
300  END IF
301 *
302 * For types 3-5, zero one row and column of the matrix to
303 * test that INFO is returned correctly.
304 *
305  IF( zerot ) THEN
306  IF( imat.EQ.3 ) THEN
307  izero = 1
308  ELSE IF( imat.EQ.4 ) THEN
309  izero = n
310  ELSE
311  izero = n / 2 + 1
312  END IF
313 *
314 * Set row and column IZERO of A to 0.
315 *
316  IF( iuplo.EQ.1 ) THEN
317  ioff = ( izero-1 )*izero / 2
318  DO 20 i = 1, izero - 1
319  a( ioff+i ) = zero
320  20 CONTINUE
321  ioff = ioff + izero
322  DO 30 i = izero, n
323  a( ioff ) = zero
324  ioff = ioff + i
325  30 CONTINUE
326  ELSE
327  ioff = izero
328  DO 40 i = 1, izero - 1
329  a( ioff ) = zero
330  ioff = ioff + n - i
331  40 CONTINUE
332  ioff = ioff - izero
333  DO 50 i = izero, n
334  a( ioff+i ) = zero
335  50 CONTINUE
336  END IF
337  ELSE
338  izero = 0
339  END IF
340 *
341 * Set the imaginary part of the diagonals.
342 *
343  IF( iuplo.EQ.1 ) THEN
344  CALL zlaipd( n, a, 2, 1 )
345  ELSE
346  CALL zlaipd( n, a, n, -1 )
347  END IF
348 *
349 * Save a copy of the matrix A in ASAV.
350 *
351  CALL zcopy( 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 ZPPSVX (FACT = 'N' reuses
376 * the condition number from the previous iteration
377 * with FACT = 'F').
378 *
379  CALL zcopy( 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 zppequ( 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 zlaqhp( 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 ZGET04.
400 *
401  IF( equil )
402  $ roldc = rcondc
403 *
404 * Compute the 1-norm of A.
405 *
406  anorm = zlanhp( '1', uplo, n, afac, rwork )
407 *
408 * Factor the matrix A.
409 *
410  CALL zpptrf( uplo, n, afac, info )
411 *
412 * Form the inverse of A.
413 *
414  CALL zcopy( npp, afac, 1, a, 1 )
415  CALL zpptri( uplo, n, a, info )
416 *
417 * Compute the 1-norm condition number of A.
418 *
419  ainvnm = zlanhp( '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 zcopy( npp, asav, 1, a, 1 )
430 *
431 * Form an exact solution and set the right hand side.
432 *
433  srnamt = 'ZLARHS'
434  CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
435  $ nrhs, a, lda, xact, lda, b, lda,
436  $ iseed, info )
437  xtype = 'C'
438  CALL zlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
439 *
440  IF( nofact ) THEN
441 *
442 * --- Test ZPPSV ---
443 *
444 * Compute the L*L' or U'*U factorization of the
445 * matrix and solve the system.
446 *
447  CALL zcopy( npp, a, 1, afac, 1 )
448  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
449 *
450  srnamt = 'ZPPSV '
451  CALL zppsv( uplo, n, nrhs, afac, x, lda, info )
452 *
453 * Check error code from ZPPSV .
454 *
455  IF( info.NE.izero ) THEN
456  CALL alaerh( path, 'ZPPSV ', 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 zppt01( uplo, n, a, afac, rwork,
468  $ result( 1 ) )
469 *
470 * Compute residual of the computed solution.
471 *
472  CALL zlacpy( 'Full', n, nrhs, b, lda, work,
473  $ lda )
474  CALL zppt02( uplo, n, nrhs, a, x, lda, work,
475  $ lda, rwork, result( 2 ) )
476 *
477 * Check solution from generated exact solution.
478 *
479  CALL zget04( 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 )'ZPPSV ', 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 ZPPSVX ---
500 *
501  IF( .NOT.prefac .AND. npp.GT.0 )
502  $ CALL zlaset( 'Full', npp, 1, dcmplx( zero ),
503  $ dcmplx( zero ), afac, npp )
504  CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
505  $ dcmplx( zero ), x, lda )
506  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
507 *
508 * Equilibrate the matrix if FACT='F' and
509 * EQUED='Y'.
510 *
511  CALL zlaqhp( uplo, n, a, s, scond, amax, equed )
512  END IF
513 *
514 * Solve the system and compute the condition number
515 * and error bounds using ZPPSVX.
516 *
517  srnamt = 'ZPPSVX'
518  CALL zppsvx( fact, uplo, n, nrhs, a, afac, equed,
519  $ s, b, lda, x, lda, rcond, rwork,
520  $ rwork( nrhs+1 ), work,
521  $ rwork( 2*nrhs+1 ), info )
522 *
523 * Check the error code from ZPPSVX.
524 *
525  IF( info.NE.izero ) THEN
526  CALL alaerh( path, 'ZPPSVX', info, izero,
527  $ fact // uplo, n, n, -1, -1, nrhs,
528  $ imat, nfail, nerrs, nout )
529  GO TO 90
530  END IF
531 *
532  IF( info.EQ.0 ) THEN
533  IF( .NOT.prefac ) THEN
534 *
535 * Reconstruct matrix from factors and compute
536 * residual.
537 *
538  CALL zppt01( uplo, n, a, afac,
539  $ rwork( 2*nrhs+1 ), result( 1 ) )
540  k1 = 1
541  ELSE
542  k1 = 2
543  END IF
544 *
545 * Compute residual of the computed solution.
546 *
547  CALL zlacpy( 'Full', n, nrhs, bsav, lda, work,
548  $ lda )
549  CALL zppt02( uplo, n, nrhs, asav, x, lda, work,
550  $ lda, rwork( 2*nrhs+1 ),
551  $ result( 2 ) )
552 *
553 * Check solution from generated exact solution.
554 *
555  IF( nofact .OR. ( prefac .AND. lsame( equed,
556  $ 'N' ) ) ) THEN
557  CALL zget04( n, nrhs, x, lda, xact, lda,
558  $ rcondc, result( 3 ) )
559  ELSE
560  CALL zget04( n, nrhs, x, lda, xact, lda,
561  $ roldc, result( 3 ) )
562  END IF
563 *
564 * Check the error bounds from iterative
565 * refinement.
566 *
567  CALL zppt05( uplo, n, nrhs, asav, b, lda, x,
568  $ lda, xact, lda, rwork,
569  $ rwork( nrhs+1 ), result( 4 ) )
570  ELSE
571  k1 = 6
572  END IF
573 *
574 * Compare RCOND from ZPPSVX with the computed value
575 * in RCONDC.
576 *
577  result( 6 ) = dget06( rcond, rcondc )
578 *
579 * Print information about the tests that did not pass
580 * the threshold.
581 *
582  DO 80 k = k1, 6
583  IF( result( k ).GE.thresh ) THEN
584  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
585  $ CALL aladhd( nout, path )
586  IF( prefac ) THEN
587  WRITE( nout, fmt = 9997 )'ZPPSVX', fact,
588  $ uplo, n, equed, imat, k, result( k )
589  ELSE
590  WRITE( nout, fmt = 9998 )'ZPPSVX', fact,
591  $ uplo, n, imat, k, result( k )
592  END IF
593  nfail = nfail + 1
594  END IF
595  80 CONTINUE
596  nrun = nrun + 7 - k1
597  90 CONTINUE
598  100 CONTINUE
599  110 CONTINUE
600  120 CONTINUE
601  130 CONTINUE
602  140 CONTINUE
603 *
604 * Print a summary of the results.
605 *
606  CALL alasvm( path, nout, nfail, nrun, nerrs )
607 *
608  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
609  $ ', test(', i1, ')=', g12.5 )
610  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
611  $ ', type ', i1, ', test(', i1, ')=', g12.5 )
612  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
613  $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ')=',
614  $ g12.5 )
615  RETURN
616 *
617 * End of ZDRVPP
618 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
double precision function zlanhp(NORM, UPLO, N, AP, WORK)
ZLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix supplied in packed form.
Definition: zlanhp.f:119
subroutine zlaqhp(UPLO, N, AP, S, SCOND, AMAX, EQUED)
ZLAQHP scales a Hermitian matrix stored in packed form.
Definition: zlaqhp.f:128
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zlaipd(N, A, INDA, VINDA)
ZLAIPD
Definition: zlaipd.f:85
subroutine zppsv(UPLO, N, NRHS, AP, B, LDB, INFO)
ZPPSV computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: zppsv.f:146
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine zppt02(UPLO, N, NRHS, A, X, LDX, B, LDB, RWORK, RESID)
ZPPT02
Definition: zppt02.f:125
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine zppt05(UPLO, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
ZPPT05
Definition: zppt05.f:159
subroutine zerrvx(PATH, NUNIT)
ZERRVX
Definition: zerrvx.f:57
subroutine zppequ(UPLO, N, AP, S, SCOND, AMAX, INFO)
ZPPEQU
Definition: zppequ.f:119
subroutine zppsvx(FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO)
ZPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: zppsvx.f:313
subroutine zpptri(UPLO, N, AP, INFO)
ZPPTRI
Definition: zpptri.f:95
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zppt01(UPLO, N, A, AFAC, RWORK, RESID)
ZPPT01
Definition: zppt01.f:97
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zpptrf(UPLO, N, AP, INFO)
ZPPTRF
Definition: zpptrf.f:121

Here is the call graph for this function:

Here is the caller graph for this function: