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

◆ zdrvpp()

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.

Definition at line 156 of file zdrvpp.f.

159*
160* -- LAPACK test routine --
161* -- LAPACK is a software package provided by Univ. of Tennessee, --
162* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163*
164* .. Scalar Arguments ..
165 LOGICAL TSTERR
166 INTEGER NMAX, NN, NOUT, NRHS
167 DOUBLE PRECISION THRESH
168* ..
169* .. Array Arguments ..
170 LOGICAL DOTYPE( * )
171 INTEGER NVAL( * )
172 DOUBLE PRECISION RWORK( * ), S( * )
173 COMPLEX*16 A( * ), AFAC( * ), ASAV( * ), B( * ),
174 $ BSAV( * ), WORK( * ), X( * ), XACT( * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 DOUBLE PRECISION ONE, ZERO
181 parameter( one = 1.0d+0, zero = 0.0d+0 )
182 INTEGER NTYPES
183 parameter( ntypes = 9 )
184 INTEGER NTESTS
185 parameter( ntests = 6 )
186* ..
187* .. Local Scalars ..
188 LOGICAL EQUIL, NOFACT, PREFAC, ZEROT
189 CHARACTER DIST, EQUED, FACT, PACKIT, TYPE, UPLO, XTYPE
190 CHARACTER*3 PATH
191 INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
192 $ IZERO, K, K1, KL, KU, LDA, MODE, N, NERRS,
193 $ NFACT, NFAIL, NIMAT, NPP, NRUN, NT
194 DOUBLE PRECISION AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
195 $ ROLDC, SCOND
196* ..
197* .. Local Arrays ..
198 CHARACTER EQUEDS( 2 ), FACTS( 3 ), PACKS( 2 ), UPLOS( 2 )
199 INTEGER ISEED( 4 ), ISEEDY( 4 )
200 DOUBLE PRECISION RESULT( NTESTS )
201* ..
202* .. External Functions ..
203 LOGICAL LSAME
204 DOUBLE PRECISION DGET06, ZLANHP
205 EXTERNAL lsame, dget06, zlanhp
206* ..
207* .. External Subroutines ..
208 EXTERNAL aladhd, alaerh, alasvm, zcopy, zerrvx, zget04,
212* ..
213* .. Scalars in Common ..
214 LOGICAL LERR, OK
215 CHARACTER*32 SRNAMT
216 INTEGER INFOT, NUNIT
217* ..
218* .. Common blocks ..
219 COMMON / infoc / infot, nunit, ok, lerr
220 COMMON / srnamc / srnamt
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC dcmplx, max
224* ..
225* .. Data statements ..
226 DATA iseedy / 1988, 1989, 1990, 1991 /
227 DATA uplos / 'U', 'L' / , facts / 'F', 'N', 'E' / ,
228 $ packs / 'C', 'R' / , equeds / 'N', 'Y' /
229* ..
230* .. Executable Statements ..
231*
232* Initialize constants and the random number seed.
233*
234 path( 1: 1 ) = 'Zomplex precision'
235 path( 2: 3 ) = 'PP'
236 nrun = 0
237 nfail = 0
238 nerrs = 0
239 DO 10 i = 1, 4
240 iseed( i ) = iseedy( i )
241 10 CONTINUE
242*
243* Test the error exits
244*
245 IF( tsterr )
246 $ CALL zerrvx( path, nout )
247 infot = 0
248*
249* Do for each value of N in NVAL
250*
251 DO 140 in = 1, nn
252 n = nval( in )
253 lda = max( n, 1 )
254 npp = n*( n+1 ) / 2
255 xtype = 'N'
256 nimat = ntypes
257 IF( n.LE.0 )
258 $ nimat = 1
259*
260 DO 130 imat = 1, nimat
261*
262* Do the tests only if DOTYPE( IMAT ) is true.
263*
264 IF( .NOT.dotype( imat ) )
265 $ GO TO 130
266*
267* Skip types 3, 4, or 5 if the matrix size is too small.
268*
269 zerot = imat.GE.3 .AND. imat.LE.5
270 IF( zerot .AND. n.LT.imat-2 )
271 $ GO TO 130
272*
273* Do first for UPLO = 'U', then for UPLO = 'L'
274*
275 DO 120 iuplo = 1, 2
276 uplo = uplos( iuplo )
277 packit = packs( iuplo )
278*
279* Set up parameters with ZLATB4 and generate a test matrix
280* with ZLATMS.
281*
282 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
283 $ CNDNUM, DIST )
284 rcondc = one / cndnum
285*
286 srnamt = 'ZLATMS'
287 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
288 $ CNDNUM, ANORM, KL, KU, PACKIT, A, LDA, WORK,
289 $ INFO )
290*
291* Check error code from ZLATMS.
292*
293 IF( info.NE.0 ) THEN
294 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
295 $ -1, -1, imat, nfail, nerrs, nout )
296 GO TO 120
297 END IF
298*
299* For types 3-5, zero one row and column of the matrix to
300* test that INFO is returned correctly.
301*
302 IF( zerot ) THEN
303 IF( imat.EQ.3 ) THEN
304 izero = 1
305 ELSE IF( imat.EQ.4 ) THEN
306 izero = n
307 ELSE
308 izero = n / 2 + 1
309 END IF
310*
311* Set row and column IZERO of A to 0.
312*
313 IF( iuplo.EQ.1 ) THEN
314 ioff = ( izero-1 )*izero / 2
315 DO 20 i = 1, izero - 1
316 a( ioff+i ) = zero
317 20 CONTINUE
318 ioff = ioff + izero
319 DO 30 i = izero, n
320 a( ioff ) = zero
321 ioff = ioff + i
322 30 CONTINUE
323 ELSE
324 ioff = izero
325 DO 40 i = 1, izero - 1
326 a( ioff ) = zero
327 ioff = ioff + n - i
328 40 CONTINUE
329 ioff = ioff - izero
330 DO 50 i = izero, n
331 a( ioff+i ) = zero
332 50 CONTINUE
333 END IF
334 ELSE
335 izero = 0
336 END IF
337*
338* Set the imaginary part of the diagonals.
339*
340 IF( iuplo.EQ.1 ) THEN
341 CALL zlaipd( n, a, 2, 1 )
342 ELSE
343 CALL zlaipd( n, a, n, -1 )
344 END IF
345*
346* Save a copy of the matrix A in ASAV.
347*
348 CALL zcopy( npp, a, 1, asav, 1 )
349*
350 DO 110 iequed = 1, 2
351 equed = equeds( iequed )
352 IF( iequed.EQ.1 ) THEN
353 nfact = 3
354 ELSE
355 nfact = 1
356 END IF
357*
358 DO 100 ifact = 1, nfact
359 fact = facts( ifact )
360 prefac = lsame( fact, 'F' )
361 nofact = lsame( fact, 'N' )
362 equil = lsame( fact, 'E' )
363*
364 IF( zerot ) THEN
365 IF( prefac )
366 $ GO TO 100
367 rcondc = zero
368*
369 ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
370*
371* Compute the condition number for comparison with
372* the value returned by ZPPSVX (FACT = 'N' reuses
373* the condition number from the previous iteration
374* with FACT = 'F').
375*
376 CALL zcopy( npp, asav, 1, afac, 1 )
377 IF( equil .OR. iequed.GT.1 ) THEN
378*
379* Compute row and column scale factors to
380* equilibrate the matrix A.
381*
382 CALL zppequ( uplo, n, afac, s, scond, amax,
383 $ info )
384 IF( info.EQ.0 .AND. n.GT.0 ) THEN
385 IF( iequed.GT.1 )
386 $ scond = zero
387*
388* Equilibrate the matrix.
389*
390 CALL zlaqhp( uplo, n, afac, s, scond,
391 $ amax, equed )
392 END IF
393 END IF
394*
395* Save the condition number of the
396* non-equilibrated system for use in ZGET04.
397*
398 IF( equil )
399 $ roldc = rcondc
400*
401* Compute the 1-norm of A.
402*
403 anorm = zlanhp( '1', uplo, n, afac, rwork )
404*
405* Factor the matrix A.
406*
407 CALL zpptrf( uplo, n, afac, info )
408*
409* Form the inverse of A.
410*
411 CALL zcopy( npp, afac, 1, a, 1 )
412 CALL zpptri( uplo, n, a, info )
413*
414* Compute the 1-norm condition number of A.
415*
416 ainvnm = zlanhp( '1', uplo, n, a, rwork )
417 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
418 rcondc = one
419 ELSE
420 rcondc = ( one / anorm ) / ainvnm
421 END IF
422 END IF
423*
424* Restore the matrix A.
425*
426 CALL zcopy( npp, asav, 1, a, 1 )
427*
428* Form an exact solution and set the right hand side.
429*
430 srnamt = 'ZLARHS'
431 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
432 $ nrhs, a, lda, xact, lda, b, lda,
433 $ iseed, info )
434 xtype = 'C'
435 CALL zlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
436*
437 IF( nofact ) THEN
438*
439* --- Test ZPPSV ---
440*
441* Compute the L*L' or U'*U factorization of the
442* matrix and solve the system.
443*
444 CALL zcopy( npp, a, 1, afac, 1 )
445 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
446*
447 srnamt = 'ZPPSV '
448 CALL zppsv( uplo, n, nrhs, afac, x, lda, info )
449*
450* Check error code from ZPPSV .
451*
452 IF( info.NE.izero ) THEN
453 CALL alaerh( path, 'ZPPSV ', info, izero,
454 $ uplo, n, n, -1, -1, nrhs, imat,
455 $ nfail, nerrs, nout )
456 GO TO 70
457 ELSE IF( info.NE.0 ) THEN
458 GO TO 70
459 END IF
460*
461* Reconstruct matrix from factors and compute
462* residual.
463*
464 CALL zppt01( uplo, n, a, afac, rwork,
465 $ result( 1 ) )
466*
467* Compute residual of the computed solution.
468*
469 CALL zlacpy( 'Full', n, nrhs, b, lda, work,
470 $ lda )
471 CALL zppt02( uplo, n, nrhs, a, x, lda, work,
472 $ lda, rwork, result( 2 ) )
473*
474* Check solution from generated exact solution.
475*
476 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
477 $ result( 3 ) )
478 nt = 3
479*
480* Print information about the tests that did not
481* pass the threshold.
482*
483 DO 60 k = 1, nt
484 IF( result( k ).GE.thresh ) THEN
485 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
486 $ CALL aladhd( nout, path )
487 WRITE( nout, fmt = 9999 )'ZPPSV ', uplo,
488 $ n, imat, k, result( k )
489 nfail = nfail + 1
490 END IF
491 60 CONTINUE
492 nrun = nrun + nt
493 70 CONTINUE
494 END IF
495*
496* --- Test ZPPSVX ---
497*
498 IF( .NOT.prefac .AND. npp.GT.0 )
499 $ CALL zlaset( 'Full', npp, 1, dcmplx( zero ),
500 $ dcmplx( zero ), afac, npp )
501 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
502 $ dcmplx( zero ), x, lda )
503 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
504*
505* Equilibrate the matrix if FACT='F' and
506* EQUED='Y'.
507*
508 CALL zlaqhp( uplo, n, a, s, scond, amax, equed )
509 END IF
510*
511* Solve the system and compute the condition number
512* and error bounds using ZPPSVX.
513*
514 srnamt = 'ZPPSVX'
515 CALL zppsvx( fact, uplo, n, nrhs, a, afac, equed,
516 $ s, b, lda, x, lda, rcond, rwork,
517 $ rwork( nrhs+1 ), work,
518 $ rwork( 2*nrhs+1 ), info )
519*
520* Check the error code from ZPPSVX.
521*
522 IF( info.NE.izero ) THEN
523 CALL alaerh( path, 'ZPPSVX', info, izero,
524 $ fact // uplo, n, n, -1, -1, nrhs,
525 $ imat, nfail, nerrs, nout )
526 GO TO 90
527 END IF
528*
529 IF( info.EQ.0 ) THEN
530 IF( .NOT.prefac ) THEN
531*
532* Reconstruct matrix from factors and compute
533* residual.
534*
535 CALL zppt01( uplo, n, a, afac,
536 $ rwork( 2*nrhs+1 ), result( 1 ) )
537 k1 = 1
538 ELSE
539 k1 = 2
540 END IF
541*
542* Compute residual of the computed solution.
543*
544 CALL zlacpy( 'Full', n, nrhs, bsav, lda, work,
545 $ lda )
546 CALL zppt02( uplo, n, nrhs, asav, x, lda, work,
547 $ lda, rwork( 2*nrhs+1 ),
548 $ result( 2 ) )
549*
550* Check solution from generated exact solution.
551*
552 IF( nofact .OR. ( prefac .AND. lsame( equed,
553 $ 'N' ) ) ) THEN
554 CALL zget04( n, nrhs, x, lda, xact, lda,
555 $ rcondc, result( 3 ) )
556 ELSE
557 CALL zget04( n, nrhs, x, lda, xact, lda,
558 $ roldc, result( 3 ) )
559 END IF
560*
561* Check the error bounds from iterative
562* refinement.
563*
564 CALL zppt05( uplo, n, nrhs, asav, b, lda, x,
565 $ lda, xact, lda, rwork,
566 $ rwork( nrhs+1 ), result( 4 ) )
567 ELSE
568 k1 = 6
569 END IF
570*
571* Compare RCOND from ZPPSVX with the computed value
572* in RCONDC.
573*
574 result( 6 ) = dget06( rcond, rcondc )
575*
576* Print information about the tests that did not pass
577* the threshold.
578*
579 DO 80 k = k1, 6
580 IF( result( k ).GE.thresh ) THEN
581 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
582 $ CALL aladhd( nout, path )
583 IF( prefac ) THEN
584 WRITE( nout, fmt = 9997 )'ZPPSVX', fact,
585 $ uplo, n, equed, imat, k, result( k )
586 ELSE
587 WRITE( nout, fmt = 9998 )'ZPPSVX', fact,
588 $ uplo, n, imat, k, result( k )
589 END IF
590 nfail = nfail + 1
591 END IF
592 80 CONTINUE
593 nrun = nrun + 7 - k1
594 90 CONTINUE
595 100 CONTINUE
596 110 CONTINUE
597 120 CONTINUE
598 130 CONTINUE
599 140 CONTINUE
600*
601* Print a summary of the results.
602*
603 CALL alasvm( path, nout, nfail, nrun, nerrs )
604*
605 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
606 $ ', test(', i1, ')=', g12.5 )
607 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
608 $ ', type ', i1, ', test(', i1, ')=', g12.5 )
609 9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
610 $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ')=',
611 $ g12.5 )
612 RETURN
613*
614* End of ZDRVPP
615*
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.f:208
subroutine aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
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,...
Definition zlanhp.f:117
subroutine zlaqhp(uplo, n, ap, s, scond, amax, equed)
ZLAQHP scales a Hermitian matrix stored in packed form.
Definition zlaqhp.f:126
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:106
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zppequ(uplo, n, ap, s, scond, amax, info)
ZPPEQU
Definition zppequ.f:117
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:144
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:311
subroutine zpptrf(uplo, n, ap, info)
ZPPTRF
Definition zpptrf.f:119
subroutine zpptri(uplo, n, ap, info)
ZPPTRI
Definition zpptri.f:93
subroutine zerrvx(path, nunit)
ZERRVX
Definition zerrvx.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlaipd(n, a, inda, vinda)
ZLAIPD
Definition zlaipd.f:83
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332
subroutine zppt01(uplo, n, a, afac, rwork, resid)
ZPPT01
Definition zppt01.f:95
subroutine zppt02(uplo, n, nrhs, a, x, ldx, b, ldb, rwork, resid)
ZPPT02
Definition zppt02.f:123
subroutine zppt05(uplo, n, nrhs, ap, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPPT05
Definition zppt05.f:157
Here is the call graph for this function:
Here is the caller graph for this function: