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

◆ zchkhp()

subroutine zchkhp ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nns,
integer, dimension( * )  nsval,
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 
)

ZCHKHP

Purpose:
 ZCHKHP tests ZHPTRF, -TRI, -TRS, -RFS, and -CON
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]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]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*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(2,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array,
                                 dimension (NMAX+2*NSMAX)
[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 161 of file zchkhp.f.

164*
165* -- LAPACK test routine --
166* -- LAPACK is a software package provided by Univ. of Tennessee, --
167* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
168*
169* .. Scalar Arguments ..
170 LOGICAL TSTERR
171 INTEGER NMAX, NN, NNS, NOUT
172 DOUBLE PRECISION THRESH
173* ..
174* .. Array Arguments ..
175 LOGICAL DOTYPE( * )
176 INTEGER IWORK( * ), NSVAL( * ), NVAL( * )
177 DOUBLE PRECISION RWORK( * )
178 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
179 $ WORK( * ), X( * ), XACT( * )
180* ..
181*
182* =====================================================================
183*
184* .. Parameters ..
185 DOUBLE PRECISION ZERO
186 parameter( zero = 0.0d+0 )
187 INTEGER NTYPES
188 parameter( ntypes = 10 )
189 INTEGER NTESTS
190 parameter( ntests = 8 )
191* ..
192* .. Local Scalars ..
193 LOGICAL TRFCON, ZEROT
194 CHARACTER DIST, PACKIT, TYPE, UPLO, XTYPE
195 CHARACTER*3 PATH
196 INTEGER I, I1, I2, IMAT, IN, INFO, IOFF, IRHS, IUPLO,
197 $ IZERO, J, K, KL, KU, LDA, MODE, N, NERRS,
198 $ NFAIL, NIMAT, NPP, NRHS, NRUN, NT
199 DOUBLE PRECISION ANORM, CNDNUM, RCOND, RCONDC
200* ..
201* .. Local Arrays ..
202 CHARACTER UPLOS( 2 )
203 INTEGER ISEED( 4 ), ISEEDY( 4 )
204 DOUBLE PRECISION RESULT( NTESTS )
205* ..
206* .. External Functions ..
207 LOGICAL LSAME
208 DOUBLE PRECISION DGET06, ZLANHP
209 EXTERNAL lsame, dget06, zlanhp
210* ..
211* .. External Subroutines ..
212 EXTERNAL alaerh, alahd, alasum, zcopy, zerrsy, zget04,
215 $ zppt03, zppt05
216* ..
217* .. Intrinsic Functions ..
218 INTRINSIC max, min
219* ..
220* .. Scalars in Common ..
221 LOGICAL LERR, OK
222 CHARACTER*32 SRNAMT
223 INTEGER INFOT, NUNIT
224* ..
225* .. Common blocks ..
226 COMMON / infoc / infot, nunit, ok, lerr
227 COMMON / srnamc / srnamt
228* ..
229* .. Data statements ..
230 DATA iseedy / 1988, 1989, 1990, 1991 /
231 DATA uplos / 'U', 'L' /
232* ..
233* .. Executable Statements ..
234*
235* Initialize constants and the random number seed.
236*
237 path( 1: 1 ) = 'Zomplex precision'
238 path( 2: 3 ) = 'HP'
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 zerrsy( path, nout )
250 infot = 0
251*
252* Do for each value of N in NVAL
253*
254 DO 170 in = 1, nn
255 n = nval( in )
256 lda = max( n, 1 )
257 xtype = 'N'
258 nimat = ntypes
259 IF( n.LE.0 )
260 $ nimat = 1
261*
262 izero = 0
263 DO 160 imat = 1, nimat
264*
265* Do the tests only if DOTYPE( IMAT ) is true.
266*
267 IF( .NOT.dotype( imat ) )
268 $ GO TO 160
269*
270* Skip types 3, 4, 5, or 6 if the matrix size is too small.
271*
272 zerot = imat.GE.3 .AND. imat.LE.6
273 IF( zerot .AND. n.LT.imat-2 )
274 $ GO TO 160
275*
276* Do first for UPLO = 'U', then for UPLO = 'L'
277*
278 DO 150 iuplo = 1, 2
279 uplo = uplos( iuplo )
280 IF( lsame( uplo, 'U' ) ) THEN
281 packit = 'C'
282 ELSE
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 150
303 END IF
304*
305* For types 3-6, zero one or more rows and columns of
306* the 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* Compute the L*D*L' or U*D*U' factorization of the matrix.
381*
382 npp = n*( n+1 ) / 2
383 CALL zcopy( npp, a, 1, afac, 1 )
384 srnamt = 'ZHPTRF'
385 CALL zhptrf( uplo, n, afac, iwork, info )
386*
387* Adjust the expected value of INFO to account for
388* pivoting.
389*
390 k = izero
391 IF( k.GT.0 ) THEN
392 100 CONTINUE
393 IF( iwork( k ).LT.0 ) THEN
394 IF( iwork( k ).NE.-k ) THEN
395 k = -iwork( k )
396 GO TO 100
397 END IF
398 ELSE IF( iwork( k ).NE.k ) THEN
399 k = iwork( k )
400 GO TO 100
401 END IF
402 END IF
403*
404* Check error code from ZHPTRF.
405*
406 IF( info.NE.k )
407 $ CALL alaerh( path, 'ZHPTRF', info, k, uplo, n, n, -1,
408 $ -1, -1, imat, nfail, nerrs, nout )
409 IF( info.NE.0 ) THEN
410 trfcon = .true.
411 ELSE
412 trfcon = .false.
413 END IF
414*
415*+ TEST 1
416* Reconstruct matrix from factors and compute residual.
417*
418 CALL zhpt01( uplo, n, a, afac, iwork, ainv, lda, rwork,
419 $ result( 1 ) )
420 nt = 1
421*
422*+ TEST 2
423* Form the inverse and compute the residual.
424*
425 IF( .NOT.trfcon ) THEN
426 CALL zcopy( npp, afac, 1, ainv, 1 )
427 srnamt = 'ZHPTRI'
428 CALL zhptri( uplo, n, ainv, iwork, work, info )
429*
430* Check error code from ZHPTRI.
431*
432 IF( info.NE.0 )
433 $ CALL alaerh( path, 'ZHPTRI', info, 0, uplo, n, n,
434 $ -1, -1, -1, imat, nfail, nerrs, nout )
435*
436 CALL zppt03( uplo, n, a, ainv, work, lda, rwork,
437 $ rcondc, result( 2 ) )
438 nt = 2
439 END IF
440*
441* Print information about the tests that did not pass
442* the threshold.
443*
444 DO 110 k = 1, nt
445 IF( result( k ).GE.thresh ) THEN
446 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
447 $ CALL alahd( nout, path )
448 WRITE( nout, fmt = 9999 )uplo, n, imat, k,
449 $ result( k )
450 nfail = nfail + 1
451 END IF
452 110 CONTINUE
453 nrun = nrun + nt
454*
455* Do only the condition estimate if INFO is not 0.
456*
457 IF( trfcon ) THEN
458 rcondc = zero
459 GO TO 140
460 END IF
461*
462 DO 130 irhs = 1, nns
463 nrhs = nsval( irhs )
464*
465*+ TEST 3
466* Solve and compute residual for A * X = B.
467*
468 srnamt = 'ZLARHS'
469 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
470 $ nrhs, a, lda, xact, lda, b, lda, iseed,
471 $ info )
472 xtype = 'C'
473 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
474*
475 srnamt = 'ZHPTRS'
476 CALL zhptrs( uplo, n, nrhs, afac, iwork, x, lda,
477 $ info )
478*
479* Check error code from ZHPTRS.
480*
481 IF( info.NE.0 )
482 $ CALL alaerh( path, 'ZHPTRS', info, 0, uplo, n, n,
483 $ -1, -1, nrhs, imat, nfail, nerrs,
484 $ nout )
485*
486 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
487 CALL zppt02( uplo, n, nrhs, a, x, lda, work, lda,
488 $ rwork, result( 3 ) )
489*
490*+ TEST 4
491* Check solution from generated exact solution.
492*
493 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
494 $ result( 4 ) )
495*
496*+ TESTS 5, 6, and 7
497* Use iterative refinement to improve the solution.
498*
499 srnamt = 'ZHPRFS'
500 CALL zhprfs( uplo, n, nrhs, a, afac, iwork, b, lda, x,
501 $ lda, rwork, rwork( nrhs+1 ), work,
502 $ rwork( 2*nrhs+1 ), info )
503*
504* Check error code from ZHPRFS.
505*
506 IF( info.NE.0 )
507 $ CALL alaerh( path, 'ZHPRFS', info, 0, uplo, n, n,
508 $ -1, -1, nrhs, imat, nfail, nerrs,
509 $ nout )
510*
511 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
512 $ result( 5 ) )
513 CALL zppt05( uplo, n, nrhs, a, b, lda, x, lda, xact,
514 $ lda, rwork, rwork( nrhs+1 ),
515 $ result( 6 ) )
516*
517* Print information about the tests that did not pass
518* the threshold.
519*
520 DO 120 k = 3, 7
521 IF( result( k ).GE.thresh ) THEN
522 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
523 $ CALL alahd( nout, path )
524 WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat,
525 $ k, result( k )
526 nfail = nfail + 1
527 END IF
528 120 CONTINUE
529 nrun = nrun + 5
530 130 CONTINUE
531*
532*+ TEST 8
533* Get an estimate of RCOND = 1/CNDNUM.
534*
535 140 CONTINUE
536 anorm = zlanhp( '1', uplo, n, a, rwork )
537 srnamt = 'ZHPCON'
538 CALL zhpcon( uplo, n, afac, iwork, anorm, rcond, work,
539 $ info )
540*
541* Check error code from ZHPCON.
542*
543 IF( info.NE.0 )
544 $ CALL alaerh( path, 'ZHPCON', info, 0, uplo, n, n, -1,
545 $ -1, -1, imat, nfail, nerrs, nout )
546*
547 result( 8 ) = dget06( rcond, rcondc )
548*
549* Print the test ratio if it is .GE. THRESH.
550*
551 IF( result( 8 ).GE.thresh ) THEN
552 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
553 $ CALL alahd( nout, path )
554 WRITE( nout, fmt = 9999 )uplo, n, imat, 8,
555 $ result( 8 )
556 nfail = nfail + 1
557 END IF
558 nrun = nrun + 1
559 150 CONTINUE
560 160 CONTINUE
561 170 CONTINUE
562*
563* Print a summary of the results.
564*
565 CALL alasum( path, nout, nfail, nrun, nerrs )
566*
567 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', type ', i2, ', test ',
568 $ i2, ', ratio =', g12.5 )
569 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
570 $ i2, ', test(', i2, ') =', g12.5 )
571 RETURN
572*
573* End of ZCHKHP
574*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.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 alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zhpcon(uplo, n, ap, ipiv, anorm, rcond, work, info)
ZHPCON
Definition zhpcon.f:118
subroutine zhprfs(uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZHPRFS
Definition zhprfs.f:180
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 zhptrs(uplo, n, nrhs, ap, ipiv, b, ldb, info)
ZHPTRS
Definition zhptrs.f:115
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zerrsy(path, nunit)
ZERRSY
Definition zerrsy.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 zppt03(uplo, n, a, ainv, work, ldwork, rwork, rcond, resid)
ZPPT03
Definition zppt03.f:110
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: