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

◆ zchkhe()

subroutine zchkhe ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
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 
)

ZCHKHE

Purpose:
 ZCHKHE tests ZHETRF, -TRI2, -TRS, -TRS2, -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]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the blocksize NB.
[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)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
[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(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(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 168 of file zchkhe.f.

171*
172* -- LAPACK test routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 LOGICAL TSTERR
178 INTEGER NMAX, NN, NNB, NNS, NOUT
179 DOUBLE PRECISION THRESH
180* ..
181* .. Array Arguments ..
182 LOGICAL DOTYPE( * )
183 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
184 DOUBLE PRECISION RWORK( * )
185 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
186 $ WORK( * ), X( * ), XACT( * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 DOUBLE PRECISION ZERO
193 parameter( zero = 0.0d+0 )
194 COMPLEX*16 CZERO
195 parameter( czero = ( 0.0d+0, 0.0d+0 ) )
196 INTEGER NTYPES
197 parameter( ntypes = 10 )
198 INTEGER NTESTS
199 parameter( ntests = 9 )
200* ..
201* .. Local Scalars ..
202 LOGICAL TRFCON, ZEROT
203 CHARACTER DIST, TYPE, UPLO, XTYPE
204 CHARACTER*3 PATH
205 INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
206 $ IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
207 $ N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
208 DOUBLE PRECISION ANORM, CNDNUM, RCOND, RCONDC
209* ..
210* .. Local Arrays ..
211 CHARACTER UPLOS( 2 )
212 INTEGER ISEED( 4 ), ISEEDY( 4 )
213 DOUBLE PRECISION RESULT( NTESTS )
214* ..
215* .. External Functions ..
216 DOUBLE PRECISION DGET06, ZLANHE
217 EXTERNAL dget06, zlanhe
218* ..
219* .. External Subroutines ..
220 EXTERNAL alaerh, alahd, alasum, xlaenv, zerrhe, zget04,
224* ..
225* .. Intrinsic Functions ..
226 INTRINSIC max, min
227* ..
228* .. Scalars in Common ..
229 LOGICAL LERR, OK
230 CHARACTER*32 SRNAMT
231 INTEGER INFOT, NUNIT
232* ..
233* .. Common blocks ..
234 COMMON / infoc / infot, nunit, ok, lerr
235 COMMON / srnamc / srnamt
236* ..
237* .. Data statements ..
238 DATA iseedy / 1988, 1989, 1990, 1991 /
239 DATA uplos / 'U', 'L' /
240* ..
241* .. Executable Statements ..
242*
243* Initialize constants and the random number seed.
244*
245 path( 1: 1 ) = 'Zomplex precision'
246 path( 2: 3 ) = 'HE'
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 zerrhe( path, nout )
258 infot = 0
259*
260* Set the minimum block size for which the block routine should
261* be used, which will be later returned by ILAENV
262*
263 CALL xlaenv( 2, 2 )
264*
265* Do for each value of N in NVAL
266*
267 DO 180 in = 1, nn
268 n = nval( in )
269 lda = max( n, 1 )
270 xtype = 'N'
271 nimat = ntypes
272 IF( n.LE.0 )
273 $ nimat = 1
274*
275 izero = 0
276 DO 170 imat = 1, nimat
277*
278* Do the tests only if DOTYPE( IMAT ) is true.
279*
280 IF( .NOT.dotype( imat ) )
281 $ GO TO 170
282*
283* Skip types 3, 4, 5, or 6 if the matrix size is too small.
284*
285 zerot = imat.GE.3 .AND. imat.LE.6
286 IF( zerot .AND. n.LT.imat-2 )
287 $ GO TO 170
288*
289* Do first for UPLO = 'U', then for UPLO = 'L'
290*
291 DO 160 iuplo = 1, 2
292 uplo = uplos( iuplo )
293*
294* Set up parameters with ZLATB4 for the matrix generator
295* based on the type of matrix to be generated.
296*
297 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
298 $ CNDNUM, DIST )
299*
300* Generate a matrix with ZLATMS.
301*
302 srnamt = 'ZLATMS'
303 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
304 $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
305 $ INFO )
306*
307* Check error code from ZLATMS and handle error.
308*
309 IF( info.NE.0 ) THEN
310 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
311 $ -1, -1, imat, nfail, nerrs, nout )
312*
313* Skip all tests for this generated matrix
314*
315 GO TO 160
316 END IF
317*
318* For types 3-6, zero one or more rows and columns of
319* the matrix to test that INFO is returned correctly.
320*
321 IF( zerot ) THEN
322 IF( imat.EQ.3 ) THEN
323 izero = 1
324 ELSE IF( imat.EQ.4 ) THEN
325 izero = n
326 ELSE
327 izero = n / 2 + 1
328 END IF
329*
330 IF( imat.LT.6 ) THEN
331*
332* Set row and column IZERO to zero.
333*
334 IF( iuplo.EQ.1 ) THEN
335 ioff = ( izero-1 )*lda
336 DO 20 i = 1, izero - 1
337 a( ioff+i ) = czero
338 20 CONTINUE
339 ioff = ioff + izero
340 DO 30 i = izero, n
341 a( ioff ) = czero
342 ioff = ioff + lda
343 30 CONTINUE
344 ELSE
345 ioff = izero
346 DO 40 i = 1, izero - 1
347 a( ioff ) = czero
348 ioff = ioff + lda
349 40 CONTINUE
350 ioff = ioff - izero
351 DO 50 i = izero, n
352 a( ioff+i ) = czero
353 50 CONTINUE
354 END IF
355 ELSE
356 IF( iuplo.EQ.1 ) THEN
357*
358* Set the first IZERO rows and columns to zero.
359*
360 ioff = 0
361 DO 70 j = 1, n
362 i2 = min( j, izero )
363 DO 60 i = 1, i2
364 a( ioff+i ) = czero
365 60 CONTINUE
366 ioff = ioff + lda
367 70 CONTINUE
368 ELSE
369*
370* Set the last IZERO rows and columns to zero.
371*
372 ioff = 0
373 DO 90 j = 1, n
374 i1 = max( j, izero )
375 DO 80 i = i1, n
376 a( ioff+i ) = czero
377 80 CONTINUE
378 ioff = ioff + lda
379 90 CONTINUE
380 END IF
381 END IF
382 ELSE
383 izero = 0
384 END IF
385*
386* End generate test matrix A.
387*
388*
389* Set the imaginary part of the diagonals.
390*
391 CALL zlaipd( n, a, lda+1, 0 )
392*
393* Do for each value of NB in NBVAL
394*
395 DO 150 inb = 1, nnb
396*
397* Set the optimal blocksize, which will be later
398* returned by ILAENV.
399*
400 nb = nbval( inb )
401 CALL xlaenv( 1, nb )
402*
403* Copy the test matrix A into matrix AFAC which
404* will be factorized in place. This is needed to
405* preserve the test matrix A for subsequent tests.
406*
407 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
408*
409* Compute the L*D*L**T or U*D*U**T factorization of the
410* matrix. IWORK stores details of the interchanges and
411* the block structure of D. AINV is a work array for
412* block factorization, LWORK is the length of AINV.
413*
414 lwork = max( 2, nb )*lda
415 srnamt = 'ZHETRF'
416 CALL zhetrf( uplo, n, afac, lda, iwork, ainv, lwork,
417 $ info )
418*
419* Adjust the expected value of INFO to account for
420* pivoting.
421*
422 k = izero
423 IF( k.GT.0 ) THEN
424 100 CONTINUE
425 IF( iwork( k ).LT.0 ) THEN
426 IF( iwork( k ).NE.-k ) THEN
427 k = -iwork( k )
428 GO TO 100
429 END IF
430 ELSE IF( iwork( k ).NE.k ) THEN
431 k = iwork( k )
432 GO TO 100
433 END IF
434 END IF
435*
436* Check error code from ZHETRF and handle error.
437*
438 IF( info.NE.k )
439 $ CALL alaerh( path, 'ZHETRF', info, k, uplo, n, n,
440 $ -1, -1, nb, imat, nfail, nerrs, nout )
441*
442* Set the condition estimate flag if the INFO is not 0.
443*
444 IF( info.NE.0 ) THEN
445 trfcon = .true.
446 ELSE
447 trfcon = .false.
448 END IF
449*
450*+ TEST 1
451* Reconstruct matrix from factors and compute residual.
452*
453 CALL zhet01( uplo, n, a, lda, afac, lda, iwork, ainv,
454 $ lda, rwork, result( 1 ) )
455 nt = 1
456*
457*+ TEST 2
458* Form the inverse and compute the residual.
459*
460 IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
461 CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
462 srnamt = 'ZHETRI2'
463 lwork = (n+nb+1)*(nb+3)
464 CALL zhetri2( uplo, n, ainv, lda, iwork, work,
465 $ lwork, info )
466*
467* Check error code from ZHETRI and handle error.
468*
469 IF( info.NE.0 )
470 $ CALL alaerh( path, 'ZHETRI', info, -1, uplo, n,
471 $ n, -1, -1, -1, imat, nfail, nerrs,
472 $ nout )
473*
474* Compute the residual for a symmetric matrix times
475* its inverse.
476*
477 CALL zpot03( uplo, n, a, lda, ainv, lda, work, lda,
478 $ rwork, rcondc, result( 2 ) )
479 nt = 2
480 END IF
481*
482* Print information about the tests that did not pass
483* the threshold.
484*
485 DO 110 k = 1, nt
486 IF( result( k ).GE.thresh ) THEN
487 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
488 $ CALL alahd( nout, path )
489 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
490 $ result( k )
491 nfail = nfail + 1
492 END IF
493 110 CONTINUE
494 nrun = nrun + nt
495*
496* Skip the other tests if this is not the first block
497* size.
498*
499 IF( inb.GT.1 )
500 $ GO TO 150
501*
502* Do only the condition estimate if INFO is not 0.
503*
504 IF( trfcon ) THEN
505 rcondc = zero
506 GO TO 140
507 END IF
508*
509* Do for each value of NRHS in NSVAL.
510*
511 DO 130 irhs = 1, nns
512 nrhs = nsval( irhs )
513*
514*+ TEST 3 (Using TRS)
515* Solve and compute residual for A * X = B.
516*
517* Choose a set of NRHS random solution vectors
518* stored in XACT and set up the right hand side B
519*
520 srnamt = 'ZLARHS'
521 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
522 $ nrhs, a, lda, xact, lda, b, lda,
523 $ iseed, info )
524 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
525*
526 srnamt = 'ZHETRS'
527 CALL zhetrs( uplo, n, nrhs, afac, lda, iwork, x,
528 $ lda, info )
529*
530* Check error code from ZHETRS and handle error.
531*
532 IF( info.NE.0 )
533 $ CALL alaerh( path, 'ZHETRS', info, 0, uplo, n,
534 $ n, -1, -1, nrhs, imat, nfail,
535 $ nerrs, nout )
536*
537 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
538*
539* Compute the residual for the solution
540*
541 CALL zpot02( uplo, n, nrhs, a, lda, x, lda, work,
542 $ lda, rwork, result( 3 ) )
543*
544*+ TEST 4 (Using TRS2)
545* Solve and compute residual for A * X = B.
546*
547* Choose a set of NRHS random solution vectors
548* stored in XACT and set up the right hand side B
549*
550 srnamt = 'ZLARHS'
551 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
552 $ nrhs, a, lda, xact, lda, b, lda,
553 $ iseed, info )
554 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
555*
556 srnamt = 'ZHETRS2'
557 CALL zhetrs2( uplo, n, nrhs, afac, lda, iwork, x,
558 $ lda, work, info )
559*
560* Check error code from ZHETRS2 and handle error.
561*
562 IF( info.NE.0 )
563 $ CALL alaerh( path, 'ZHETRS2', info, 0, uplo, n,
564 $ n, -1, -1, nrhs, imat, nfail,
565 $ nerrs, nout )
566*
567 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
568*
569* Compute the residual for the solution
570*
571 CALL zpot02( uplo, n, nrhs, a, lda, x, lda, work,
572 $ lda, rwork, result( 4 ) )
573*
574*+ TEST 5
575* Check solution from generated exact solution.
576*
577 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
578 $ result( 5 ) )
579*
580*+ TESTS 6, 7, and 8
581* Use iterative refinement to improve the solution.
582*
583 srnamt = 'ZHERFS'
584 CALL zherfs( uplo, n, nrhs, a, lda, afac, lda,
585 $ iwork, b, lda, x, lda, rwork,
586 $ rwork( nrhs+1 ), work,
587 $ rwork( 2*nrhs+1 ), info )
588*
589* Check error code from ZHERFS.
590*
591 IF( info.NE.0 )
592 $ CALL alaerh( path, 'ZHERFS', info, 0, uplo, n,
593 $ n, -1, -1, nrhs, imat, nfail,
594 $ nerrs, nout )
595*
596 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
597 $ result( 6 ) )
598 CALL zpot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
599 $ xact, lda, rwork, rwork( nrhs+1 ),
600 $ result( 7 ) )
601*
602* Print information about the tests that did not pass
603* the threshold.
604*
605 DO 120 k = 3, 8
606 IF( result( k ).GE.thresh ) THEN
607 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
608 $ CALL alahd( nout, path )
609 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
610 $ imat, k, result( k )
611 nfail = nfail + 1
612 END IF
613 120 CONTINUE
614 nrun = nrun + 6
615*
616* End do for each value of NRHS in NSVAL.
617*
618 130 CONTINUE
619*
620*+ TEST 9
621* Get an estimate of RCOND = 1/CNDNUM.
622*
623 140 CONTINUE
624 anorm = zlanhe( '1', uplo, n, a, lda, rwork )
625 srnamt = 'ZHECON'
626 CALL zhecon( uplo, n, afac, lda, iwork, anorm, rcond,
627 $ work, info )
628*
629* Check error code from ZHECON and handle error.
630*
631 IF( info.NE.0 )
632 $ CALL alaerh( path, 'ZHECON', info, 0, uplo, n, n,
633 $ -1, -1, -1, imat, nfail, nerrs, nout )
634*
635 result( 9 ) = dget06( rcond, rcondc )
636*
637* Print information about the tests that did not pass
638* the threshold.
639*
640 IF( result( 9 ).GE.thresh ) THEN
641 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
642 $ CALL alahd( nout, path )
643 WRITE( nout, fmt = 9997 )uplo, n, imat, 9,
644 $ result( 9 )
645 nfail = nfail + 1
646 END IF
647 nrun = nrun + 1
648 150 CONTINUE
649 160 CONTINUE
650 170 CONTINUE
651 180 CONTINUE
652*
653* Print a summary of the results.
654*
655 CALL alasum( path, nout, nfail, nrun, nerrs )
656*
657 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
658 $ i2, ', test ', i2, ', ratio =', g12.5 )
659 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
660 $ i2, ', test(', i2, ') =', g12.5 )
661 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
662 $ ', test(', i2, ') =', g12.5 )
663 RETURN
664*
665* End of ZCHKHE
666*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.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 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 zhecon(uplo, n, a, lda, ipiv, anorm, rcond, work, info)
ZHECON
Definition zhecon.f:125
subroutine zherfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZHERFS
Definition zherfs.f:192
subroutine zhetrf(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRF
Definition zhetrf.f:177
subroutine zhetri2(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRI2
Definition zhetri2.f:127
subroutine zhetrs2(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, info)
ZHETRS2
Definition zhetrs2.f:127
subroutine zhetrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
ZHETRS
Definition zhetrs.f:120
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 zlanhe(norm, uplo, n, a, lda, work)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhe.f:124
subroutine zerrhe(path, nunit)
ZERRHE
Definition zerrhe.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zhet01(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
ZHET01
Definition zhet01.f:126
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 zpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZPOT02
Definition zpot02.f:127
subroutine zpot03(uplo, n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
ZPOT03
Definition zpot03.f:126
subroutine zpot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPOT05
Definition zpot05.f:165
Here is the call graph for this function:
Here is the caller graph for this function: