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

◆ cchkhe()

subroutine cchkhe ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
integer  nns,
integer, dimension( * )  nsval,
real  thresh,
logical  tsterr,
integer  nmax,
complex, dimension( * )  a,
complex, dimension( * )  afac,
complex, dimension( * )  ainv,
complex, dimension( * )  b,
complex, dimension( * )  x,
complex, dimension( * )  xact,
complex, dimension( * )  work,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

CCHKHE

Purpose:
 CCHKHE tests CHETRF, -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 REAL
          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 array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX array, dimension (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL 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 cchkhe.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 REAL THRESH
180* ..
181* .. Array Arguments ..
182 LOGICAL DOTYPE( * )
183 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
184 REAL RWORK( * )
185 COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ),
186 $ WORK( * ), X( * ), XACT( * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 REAL ZERO
193 parameter( zero = 0.0e+0 )
194 COMPLEX CZERO
195 parameter( czero = ( 0.0e+0, 0.0e+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 REAL ANORM, CNDNUM, RCOND, RCONDC
209* ..
210* .. Local Arrays ..
211 CHARACTER UPLOS( 2 )
212 INTEGER ISEED( 4 ), ISEEDY( 4 )
213 REAL RESULT( NTESTS )
214* ..
215* .. External Functions ..
216 REAL CLANHE, SGET06
217 EXTERNAL clanhe, sget06
218* ..
219* .. External Subroutines ..
220 EXTERNAL alaerh, alahd, alasum, cerrhe, cget04, checon,
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 ) = 'Complex 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 cerrhe( 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*
277* Do for each value of matrix type IMAT
278*
279 DO 170 imat = 1, nimat
280*
281* Do the tests only if DOTYPE( IMAT ) is true.
282*
283 IF( .NOT.dotype( imat ) )
284 $ GO TO 170
285*
286* Skip types 3, 4, 5, or 6 if the matrix size is too small.
287*
288 zerot = imat.GE.3 .AND. imat.LE.6
289 IF( zerot .AND. n.LT.imat-2 )
290 $ GO TO 170
291*
292* Do first for UPLO = 'U', then for UPLO = 'L'
293*
294 DO 160 iuplo = 1, 2
295 uplo = uplos( iuplo )
296*
297* Begin generate test matrix A.
298*
299*
300* Set up parameters with CLATB4 for the matrix generator
301* based on the type of matrix to be generated.
302*
303 CALL clatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
304 $ CNDNUM, DIST )
305*
306* Generate a matrix with CLATMS.
307*
308 srnamt = 'CLATMS'
309 CALL clatms( n, n, dist, iseed, TYPE, RWORK, MODE,
310 $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
311 $ INFO )
312*
313* Check error code from CLATMS and handle error.
314*
315 IF( info.NE.0 ) THEN
316 CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n, -1,
317 $ -1, -1, imat, nfail, nerrs, nout )
318*
319* Skip all tests for this generated matrix
320*
321 GO TO 160
322 END IF
323*
324* For matrix types 3-6, zero one or more rows and
325* columns of the matrix to test that INFO is returned
326* correctly.
327*
328 IF( zerot ) THEN
329 IF( imat.EQ.3 ) THEN
330 izero = 1
331 ELSE IF( imat.EQ.4 ) THEN
332 izero = n
333 ELSE
334 izero = n / 2 + 1
335 END IF
336*
337 IF( imat.LT.6 ) THEN
338*
339* Set row and column IZERO to zero.
340*
341 IF( iuplo.EQ.1 ) THEN
342 ioff = ( izero-1 )*lda
343 DO 20 i = 1, izero - 1
344 a( ioff+i ) = czero
345 20 CONTINUE
346 ioff = ioff + izero
347 DO 30 i = izero, n
348 a( ioff ) = czero
349 ioff = ioff + lda
350 30 CONTINUE
351 ELSE
352 ioff = izero
353 DO 40 i = 1, izero - 1
354 a( ioff ) = czero
355 ioff = ioff + lda
356 40 CONTINUE
357 ioff = ioff - izero
358 DO 50 i = izero, n
359 a( ioff+i ) = czero
360 50 CONTINUE
361 END IF
362 ELSE
363 IF( iuplo.EQ.1 ) THEN
364*
365* Set the first IZERO rows and columns to zero.
366*
367 ioff = 0
368 DO 70 j = 1, n
369 i2 = min( j, izero )
370 DO 60 i = 1, i2
371 a( ioff+i ) = czero
372 60 CONTINUE
373 ioff = ioff + lda
374 70 CONTINUE
375 ELSE
376*
377* Set the last IZERO rows and columns to zero.
378*
379 ioff = 0
380 DO 90 j = 1, n
381 i1 = max( j, izero )
382 DO 80 i = i1, n
383 a( ioff+i ) = czero
384 80 CONTINUE
385 ioff = ioff + lda
386 90 CONTINUE
387 END IF
388 END IF
389 ELSE
390 izero = 0
391 END IF
392*
393* Set the imaginary part of the diagonals.
394*
395 CALL claipd( n, a, lda+1, 0 )
396*
397* End generate test matrix A.
398*
399*
400* Do for each value of NB in NBVAL
401*
402 DO 150 inb = 1, nnb
403*
404* Set the optimal blocksize, which will be later
405* returned by ILAENV.
406*
407 nb = nbval( inb )
408 CALL xlaenv( 1, nb )
409*
410* Copy the test matrix A into matrix AFAC which
411* will be factorized in place. This is needed to
412* preserve the test matrix A for subsequent tests.
413*
414 CALL clacpy( uplo, n, n, a, lda, afac, lda )
415*
416* Compute the L*D*L**T or U*D*U**T factorization of the
417* matrix. IWORK stores details of the interchanges and
418* the block structure of D. AINV is a work array for
419* block factorization, LWORK is the length of AINV.
420*
421 lwork = max( 2, nb )*lda
422 srnamt = 'CHETRF'
423 CALL chetrf( uplo, n, afac, lda, iwork, ainv, lwork,
424 $ info )
425*
426* Adjust the expected value of INFO to account for
427* pivoting.
428*
429 k = izero
430 IF( k.GT.0 ) THEN
431 100 CONTINUE
432 IF( iwork( k ).LT.0 ) THEN
433 IF( iwork( k ).NE.-k ) THEN
434 k = -iwork( k )
435 GO TO 100
436 END IF
437 ELSE IF( iwork( k ).NE.k ) THEN
438 k = iwork( k )
439 GO TO 100
440 END IF
441 END IF
442*
443* Check error code from CHETRF and handle error.
444*
445 IF( info.NE.k )
446 $ CALL alaerh( path, 'CHETRF', info, k, uplo, n, n,
447 $ -1, -1, nb, imat, nfail, nerrs, nout )
448*
449* Set the condition estimate flag if the INFO is not 0.
450*
451 IF( info.NE.0 ) THEN
452 trfcon = .true.
453 ELSE
454 trfcon = .false.
455 END IF
456*
457*+ TEST 1
458* Reconstruct matrix from factors and compute residual.
459*
460 CALL chet01( uplo, n, a, lda, afac, lda, iwork, ainv,
461 $ lda, rwork, result( 1 ) )
462 nt = 1
463*
464*+ TEST 2
465* Form the inverse and compute the residual,
466* if the factorization was competed without INFO > 0
467* (i.e. there is no zero rows and columns).
468* Do it only for the first block size.
469*
470 IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
471 CALL clacpy( uplo, n, n, afac, lda, ainv, lda )
472 srnamt = 'CHETRI2'
473 lwork = (n+nb+1)*(nb+3)
474 CALL chetri2( uplo, n, ainv, lda, iwork, work,
475 $ lwork, info )
476*
477* Check error code from CHETRI2 and handle error.
478*
479 IF( info.NE.0 )
480 $ CALL alaerh( path, 'CHETRI2', info, -1, uplo, n,
481 $ n, -1, -1, -1, imat, nfail, nerrs,
482 $ nout )
483*
484* Compute the residual for a symmetric matrix times
485* its inverse.
486*
487 CALL cpot03( uplo, n, a, lda, ainv, lda, work, lda,
488 $ rwork, rcondc, result( 2 ) )
489 nt = 2
490 END IF
491*
492* Print information about the tests that did not pass
493* the threshold.
494*
495 DO 110 k = 1, nt
496 IF( result( k ).GE.thresh ) THEN
497 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
498 $ CALL alahd( nout, path )
499 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
500 $ result( k )
501 nfail = nfail + 1
502 END IF
503 110 CONTINUE
504 nrun = nrun + nt
505*
506* Skip the other tests if this is not the first block
507* size.
508*
509 IF( inb.GT.1 )
510 $ GO TO 150
511*
512* Do only the condition estimate if INFO is not 0.
513*
514 IF( trfcon ) THEN
515 rcondc = zero
516 GO TO 140
517 END IF
518*
519* Do for each value of NRHS in NSVAL.
520*
521 DO 130 irhs = 1, nns
522 nrhs = nsval( irhs )
523*
524*+ TEST 3 (Using TRS)
525* Solve and compute residual for A * X = B.
526*
527* Choose a set of NRHS random solution vectors
528* stored in XACT and set up the right hand side B
529*
530 srnamt = 'CLARHS'
531 CALL clarhs( path, xtype, uplo, ' ', n, n, kl, ku,
532 $ nrhs, a, lda, xact, lda, b, lda,
533 $ iseed, info )
534 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
535*
536 srnamt = 'CHETRS'
537 CALL chetrs( uplo, n, nrhs, afac, lda, iwork, x,
538 $ lda, info )
539*
540* Check error code from CHETRS and handle error.
541*
542 IF( info.NE.0 )
543 $ CALL alaerh( path, 'CHETRS', info, 0, uplo, n,
544 $ n, -1, -1, nrhs, imat, nfail,
545 $ nerrs, nout )
546*
547 CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
548*
549* Compute the residual for the solution
550*
551 CALL cpot02( uplo, n, nrhs, a, lda, x, lda, work,
552 $ lda, rwork, result( 3 ) )
553*
554*+ TEST 4 (Using TRS2)
555* Solve and compute residual for A * X = B.
556*
557* Choose a set of NRHS random solution vectors
558* stored in XACT and set up the right hand side B
559*
560 srnamt = 'CLARHS'
561 CALL clarhs( path, xtype, uplo, ' ', n, n, kl, ku,
562 $ nrhs, a, lda, xact, lda, b, lda,
563 $ iseed, info )
564 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
565*
566 srnamt = 'CHETRS2'
567 CALL chetrs2( uplo, n, nrhs, afac, lda, iwork, x,
568 $ lda, work, info )
569*
570* Check error code from CHETRS2 and handle error.
571*
572 IF( info.NE.0 )
573 $ CALL alaerh( path, 'CHETRS2', info, 0, uplo, n,
574 $ n, -1, -1, nrhs, imat, nfail,
575 $ nerrs, nout )
576*
577 CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
578*
579* Compute the residual for the solution
580*
581 CALL cpot02( uplo, n, nrhs, a, lda, x, lda, work,
582 $ lda, rwork, result( 4 ) )
583*
584*+ TEST 5
585* Check solution from generated exact solution.
586*
587 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
588 $ result( 5 ) )
589*
590*+ TESTS 6, 7, and 8
591* Use iterative refinement to improve the solution.
592*
593 srnamt = 'CHERFS'
594 CALL cherfs( uplo, n, nrhs, a, lda, afac, lda,
595 $ iwork, b, lda, x, lda, rwork,
596 $ rwork( nrhs+1 ), work,
597 $ rwork( 2*nrhs+1 ), info )
598*
599* Check error code from CHERFS and handle error.
600*
601 IF( info.NE.0 )
602 $ CALL alaerh( path, 'CHERFS', info, 0, uplo, n,
603 $ n, -1, -1, nrhs, imat, nfail,
604 $ nerrs, nout )
605*
606 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
607 $ result( 6 ) )
608 CALL cpot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
609 $ xact, lda, rwork, rwork( nrhs+1 ),
610 $ result( 7 ) )
611*
612* Print information about the tests that did not pass
613* the threshold.
614*
615 DO 120 k = 3, 8
616 IF( result( k ).GE.thresh ) THEN
617 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
618 $ CALL alahd( nout, path )
619 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
620 $ imat, k, result( k )
621 nfail = nfail + 1
622 END IF
623 120 CONTINUE
624 nrun = nrun + 6
625*
626* End do for each value of NRHS in NSVAL.
627*
628 130 CONTINUE
629*
630*+ TEST 9
631* Get an estimate of RCOND = 1/CNDNUM.
632*
633 140 CONTINUE
634 anorm = clanhe( '1', uplo, n, a, lda, rwork )
635 srnamt = 'CHECON'
636 CALL checon( uplo, n, afac, lda, iwork, anorm, rcond,
637 $ work, info )
638*
639* Check error code from CHECON and handle error.
640*
641 IF( info.NE.0 )
642 $ CALL alaerh( path, 'CHECON', info, 0, uplo, n, n,
643 $ -1, -1, -1, imat, nfail, nerrs, nout )
644*
645* Compute the test ratio to compare values of RCOND
646*
647 result( 9 ) = sget06( rcond, rcondc )
648*
649* Print information about the tests that did not pass
650* the threshold.
651*
652 IF( result( 9 ).GE.thresh ) THEN
653 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
654 $ CALL alahd( nout, path )
655 WRITE( nout, fmt = 9997 )uplo, n, imat, 8,
656 $ result( 9 )
657 nfail = nfail + 1
658 END IF
659 nrun = nrun + 1
660 150 CONTINUE
661 160 CONTINUE
662 170 CONTINUE
663 180 CONTINUE
664*
665* Print a summary of the results.
666*
667 CALL alasum( path, nout, nfail, nrun, nerrs )
668*
669 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
670 $ i2, ', test ', i2, ', ratio =', g12.5 )
671 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
672 $ i2, ', test(', i2, ') =', g12.5 )
673 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
674 $ ', test(', i2, ') =', g12.5 )
675 RETURN
676*
677* End of CCHKHE
678*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine clarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
CLARHS
Definition clarhs.f:208
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
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
subroutine cerrhe(path, nunit)
CERRHE
Definition cerrhe.f:55
subroutine cget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
CGET04
Definition cget04.f:102
subroutine chet01(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
CHET01
Definition chet01.f:126
subroutine claipd(n, a, inda, vinda)
CLAIPD
Definition claipd.f:83
subroutine clatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
CLATB4
Definition clatb4.f:121
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine cpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
CPOT02
Definition cpot02.f:127
subroutine cpot03(uplo, n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
CPOT03
Definition cpot03.f:126
subroutine cpot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
CPOT05
Definition cpot05.f:165
subroutine checon(uplo, n, a, lda, ipiv, anorm, rcond, work, info)
CHECON
Definition checon.f:125
subroutine cherfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CHERFS
Definition cherfs.f:192
subroutine chetrf(uplo, n, a, lda, ipiv, work, lwork, info)
CHETRF
Definition chetrf.f:177
subroutine chetri2(uplo, n, a, lda, ipiv, work, lwork, info)
CHETRI2
Definition chetri2.f:127
subroutine chetrs2(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, info)
CHETRS2
Definition chetrs2.f:127
subroutine chetrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
CHETRS
Definition chetrs.f:120
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
real function clanhe(norm, uplo, n, a, lda, work)
CLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clanhe.f:124
real function sget06(rcond, rcondc)
SGET06
Definition sget06.f:55
Here is the call graph for this function:
Here is the caller graph for this function: