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

◆ zdrvhp()

subroutine zdrvhp ( 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( * )  ainv,
complex*16, dimension( * )  b,
complex*16, dimension( * )  x,
complex*16, dimension( * )  xact,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

ZDRVHP

Purpose:
 ZDRVHP tests the driver routines ZHPSV 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]AINV
          AINV is COMPLEX*16 array, dimension
                      (NMAX*(NMAX+1)/2)
[out]B
          B 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]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(2,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.

Definition at line 154 of file zdrvhp.f.

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