LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cchkhe_rk.f
Go to the documentation of this file.
1*> \brief \b CCHKHE_RK
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE CCHKHE_RK( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
12* THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B,
13* X, XACT, WORK, RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NNB, NNS, NOUT
18* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
23* REAL RWORK( * )
24* COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ), E( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> CCHKHE_RK tests CHETRF_RK, -TRI_3, -TRS_3,
35*> and -CON_3.
36*> \endverbatim
37*
38* Arguments:
39* ==========
40*
41*> \param[in] DOTYPE
42*> \verbatim
43*> DOTYPE is LOGICAL array, dimension (NTYPES)
44*> The matrix types to be used for testing. Matrices of type j
45*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
46*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
47*> \endverbatim
48*>
49*> \param[in] NN
50*> \verbatim
51*> NN is INTEGER
52*> The number of values of N contained in the vector NVAL.
53*> \endverbatim
54*>
55*> \param[in] NVAL
56*> \verbatim
57*> NVAL is INTEGER array, dimension (NN)
58*> The values of the matrix dimension N.
59*> \endverbatim
60*>
61*> \param[in] NNB
62*> \verbatim
63*> NNB is INTEGER
64*> The number of values of NB contained in the vector NBVAL.
65*> \endverbatim
66*>
67*> \param[in] NBVAL
68*> \verbatim
69*> NBVAL is INTEGER array, dimension (NNB)
70*> The values of the blocksize NB.
71*> \endverbatim
72*>
73*> \param[in] NNS
74*> \verbatim
75*> NNS is INTEGER
76*> The number of values of NRHS contained in the vector NSVAL.
77*> \endverbatim
78*>
79*> \param[in] NSVAL
80*> \verbatim
81*> NSVAL is INTEGER array, dimension (NNS)
82*> The values of the number of right hand sides NRHS.
83*> \endverbatim
84*>
85*> \param[in] THRESH
86*> \verbatim
87*> THRESH is REAL
88*> The threshold value for the test ratios. A result is
89*> included in the output file if RESULT >= THRESH. To have
90*> every test ratio printed, use THRESH = 0.
91*> \endverbatim
92*>
93*> \param[in] TSTERR
94*> \verbatim
95*> TSTERR is LOGICAL
96*> Flag that indicates whether error exits are to be tested.
97*> \endverbatim
98*>
99*> \param[in] NMAX
100*> \verbatim
101*> NMAX is INTEGER
102*> The maximum value permitted for N, used in dimensioning the
103*> work arrays.
104*> \endverbatim
105*>
106*> \param[out] A
107*> \verbatim
108*> A is COMPLEX array, dimension (NMAX*NMAX)
109*> \endverbatim
110*>
111*> \param[out] AFAC
112*> \verbatim
113*> AFAC is COMPLEX array, dimension (NMAX*NMAX)
114*> \endverbatim
115*>
116*> \param[out] E
117*> \verbatim
118*> E is COMPLEX array, dimension (NMAX)
119*> \endverbatim
120*>
121*> \param[out] AINV
122*> \verbatim
123*> AINV is COMPLEX array, dimension (NMAX*NMAX)
124*> \endverbatim
125*>
126*> \param[out] B
127*> \verbatim
128*> B is COMPLEX array, dimension (NMAX*NSMAX)
129*> where NSMAX is the largest entry in NSVAL.
130*> \endverbatim
131*>
132*> \param[out] X
133*> \verbatim
134*> X is COMPLEX array, dimension (NMAX*NSMAX)
135*> \endverbatim
136*>
137*> \param[out] XACT
138*> \verbatim
139*> XACT is COMPLEX array, dimension (NMAX*NSMAX)
140*> \endverbatim
141*>
142*> \param[out] WORK
143*> \verbatim
144*> WORK is COMPLEX array, dimension (NMAX*max(3,NSMAX))
145*> \endverbatim
146*>
147*> \param[out] RWORK
148*> \verbatim
149*> RWORK is REAL array, dimension (max(NMAX,2*NSMAX)
150*> \endverbatim
151*>
152*> \param[out] IWORK
153*> \verbatim
154*> IWORK is INTEGER array, dimension (2*NMAX)
155*> \endverbatim
156*>
157*> \param[in] NOUT
158*> \verbatim
159*> NOUT is INTEGER
160*> The unit number for output.
161*> \endverbatim
162*
163* Authors:
164* ========
165*
166*> \author Univ. of Tennessee
167*> \author Univ. of California Berkeley
168*> \author Univ. of Colorado Denver
169*> \author NAG Ltd.
170*
171*> \ingroup complex_lin
172*
173* =====================================================================
174 SUBROUTINE cchkhe_rk( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
175 $ THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B,
176 $ X, XACT, WORK, RWORK, IWORK, NOUT )
177*
178* -- LAPACK test routine --
179* -- LAPACK is a software package provided by Univ. of Tennessee, --
180* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181*
182* .. Scalar Arguments ..
183 LOGICAL TSTERR
184 INTEGER NMAX, NN, NNB, NNS, NOUT
185 REAL THRESH
186* ..
187* .. Array Arguments ..
188 LOGICAL DOTYPE( * )
189 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
190 REAL RWORK( * )
191 COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ), E( * ),
192 $ work( * ), x( * ), xact( * )
193* ..
194*
195* =====================================================================
196*
197* .. Parameters ..
198 REAL ZERO, ONE
199 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
200 REAL ONEHALF
201 parameter( onehalf = 0.5e+0 )
202 REAL EIGHT, SEVTEN
203 parameter( eight = 8.0e+0, sevten = 17.0e+0 )
204 COMPLEX CZERO
205 parameter( czero = ( 0.0e+0, 0.0e+0 ) )
206 INTEGER NTYPES
207 parameter( ntypes = 10 )
208 INTEGER NTESTS
209 parameter( ntests = 7 )
210* ..
211* .. Local Scalars ..
212 LOGICAL TRFCON, ZEROT
213 CHARACTER DIST, TYPE, UPLO, XTYPE
214 CHARACTER*3 PATH, MATPATH
215 INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
216 $ itemp, itemp2, iuplo, izero, j, k, kl, ku, lda,
217 $ lwork, mode, n, nb, nerrs, nfail, nimat, nrhs,
218 $ nrun, nt
219 REAL ALPHA, ANORM, CNDNUM, CONST, SING_MAX,
220 $ SING_MIN, RCOND, RCONDC, STEMP
221* ..
222* .. Local Arrays ..
223 CHARACTER UPLOS( 2 )
224 INTEGER ISEED( 4 ), ISEEDY( 4 ), IDUMMY( 1 )
225 REAL RESULT( NTESTS )
226 COMPLEX BLOCK( 2, 2 ), CDUMMY( 1 )
227* ..
228* .. External Functions ..
229 REAL CLANGE, CLANHE, SGET06
230 EXTERNAL CLANGE, CLANHE, SGET06
231* ..
232* .. External Subroutines ..
233 EXTERNAL alaerh, alahd, alasum, cerrhe, cgesvd, cget04,
237* ..
238* .. Intrinsic Functions ..
239 INTRINSIC conjg, max, min, sqrt
240* ..
241* .. Scalars in Common ..
242 LOGICAL LERR, OK
243 CHARACTER*32 SRNAMT
244 INTEGER INFOT, NUNIT
245* ..
246* .. Common blocks ..
247 COMMON / infoc / infot, nunit, ok, lerr
248 COMMON / srnamc / srnamt
249* ..
250* .. Data statements ..
251 DATA iseedy / 1988, 1989, 1990, 1991 /
252 DATA uplos / 'U', 'L' /
253* ..
254* .. Executable Statements ..
255*
256* Initialize constants and the random number seed.
257*
258 alpha = ( one+sqrt( sevten ) ) / eight
259*
260* Test path
261*
262 path( 1: 1 ) = 'Complex precision'
263 path( 2: 3 ) = 'HK'
264*
265* Path to generate matrices
266*
267 matpath( 1: 1 ) = 'Complex precision'
268 matpath( 2: 3 ) = 'HE'
269*
270 nrun = 0
271 nfail = 0
272 nerrs = 0
273 DO 10 i = 1, 4
274 iseed( i ) = iseedy( i )
275 10 CONTINUE
276*
277* Test the error exits
278*
279 IF( tsterr )
280 $ CALL cerrhe( path, nout )
281 infot = 0
282*
283* Set the minimum block size for which the block routine should
284* be used, which will be later returned by ILAENV
285*
286 CALL xlaenv( 2, 2 )
287*
288* Do for each value of N in NVAL
289*
290 DO 270 in = 1, nn
291 n = nval( in )
292 lda = max( n, 1 )
293 xtype = 'N'
294 nimat = ntypes
295 IF( n.LE.0 )
296 $ nimat = 1
297*
298 izero = 0
299*
300* Do for each value of matrix type IMAT
301*
302 DO 260 imat = 1, nimat
303*
304* Do the tests only if DOTYPE( IMAT ) is true.
305*
306 IF( .NOT.dotype( imat ) )
307 $ GO TO 260
308*
309* Skip types 3, 4, 5, or 6 if the matrix size is too small.
310*
311 zerot = imat.GE.3 .AND. imat.LE.6
312 IF( zerot .AND. n.LT.imat-2 )
313 $ GO TO 260
314*
315* Do first for UPLO = 'U', then for UPLO = 'L'
316*
317 DO 250 iuplo = 1, 2
318 uplo = uplos( iuplo )
319*
320* Begin generate the test matrix A.
321*
322* Set up parameters with CLATB4 for the matrix generator
323* based on the type of matrix to be generated.
324*
325 CALL clatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
326 $ mode, cndnum, dist )
327*
328* Generate a matrix with CLATMS.
329*
330 srnamt = 'CLATMS'
331 CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
332 $ cndnum, anorm, kl, ku, uplo, a, lda,
333 $ work, info )
334*
335* Check error code from CLATMS and handle error.
336*
337 IF( info.NE.0 ) THEN
338 CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n,
339 $ -1, -1, -1, imat, nfail, nerrs, nout )
340*
341* Skip all tests for this generated matrix
342*
343 GO TO 250
344 END IF
345*
346* For matrix types 3-6, zero one or more rows and
347* columns of the matrix to test that INFO is returned
348* correctly.
349*
350 IF( zerot ) THEN
351 IF( imat.EQ.3 ) THEN
352 izero = 1
353 ELSE IF( imat.EQ.4 ) THEN
354 izero = n
355 ELSE
356 izero = n / 2 + 1
357 END IF
358*
359 IF( imat.LT.6 ) THEN
360*
361* Set row and column IZERO to zero.
362*
363 IF( iuplo.EQ.1 ) THEN
364 ioff = ( izero-1 )*lda
365 DO 20 i = 1, izero - 1
366 a( ioff+i ) = czero
367 20 CONTINUE
368 ioff = ioff + izero
369 DO 30 i = izero, n
370 a( ioff ) = czero
371 ioff = ioff + lda
372 30 CONTINUE
373 ELSE
374 ioff = izero
375 DO 40 i = 1, izero - 1
376 a( ioff ) = czero
377 ioff = ioff + lda
378 40 CONTINUE
379 ioff = ioff - izero
380 DO 50 i = izero, n
381 a( ioff+i ) = czero
382 50 CONTINUE
383 END IF
384 ELSE
385 IF( iuplo.EQ.1 ) THEN
386*
387* Set the first IZERO rows and columns to zero.
388*
389 ioff = 0
390 DO 70 j = 1, n
391 i2 = min( j, izero )
392 DO 60 i = 1, i2
393 a( ioff+i ) = czero
394 60 CONTINUE
395 ioff = ioff + lda
396 70 CONTINUE
397 ELSE
398*
399* Set the last IZERO rows and columns to zero.
400*
401 ioff = 0
402 DO 90 j = 1, n
403 i1 = max( j, izero )
404 DO 80 i = i1, n
405 a( ioff+i ) = czero
406 80 CONTINUE
407 ioff = ioff + lda
408 90 CONTINUE
409 END IF
410 END IF
411 ELSE
412 izero = 0
413 END IF
414*
415* End generate the test matrix A.
416*
417*
418* Do for each value of NB in NBVAL
419*
420 DO 240 inb = 1, nnb
421*
422* Set the optimal blocksize, which will be later
423* returned by ILAENV.
424*
425 nb = nbval( inb )
426 CALL xlaenv( 1, nb )
427*
428* Copy the test matrix A into matrix AFAC which
429* will be factorized in place. This is needed to
430* preserve the test matrix A for subsequent tests.
431*
432 CALL clacpy( uplo, n, n, a, lda, afac, lda )
433*
434* Compute the L*D*L**T or U*D*U**T factorization of the
435* matrix. IWORK stores details of the interchanges and
436* the block structure of D. AINV is a work array for
437* block factorization, LWORK is the length of AINV.
438*
439 lwork = max( 2, nb )*lda
440 srnamt = 'CHETRF_RK'
441 CALL chetrf_rk( uplo, n, afac, lda, e, iwork, ainv,
442 $ lwork, info )
443*
444* Adjust the expected value of INFO to account for
445* pivoting.
446*
447 k = izero
448 IF( k.GT.0 ) THEN
449 100 CONTINUE
450 IF( iwork( k ).LT.0 ) THEN
451 IF( iwork( k ).NE.-k ) THEN
452 k = -iwork( k )
453 GO TO 100
454 END IF
455 ELSE IF( iwork( k ).NE.k ) THEN
456 k = iwork( k )
457 GO TO 100
458 END IF
459 END IF
460*
461* Check error code from CHETRF_RK and handle error.
462*
463 IF( info.NE.k)
464 $ CALL alaerh( path, 'CHETRF_RK', info, k,
465 $ uplo, n, n, -1, -1, nb, imat,
466 $ nfail, nerrs, nout )
467*
468* Set the condition estimate flag if the INFO is not 0.
469*
470 IF( info.NE.0 ) THEN
471 trfcon = .true.
472 ELSE
473 trfcon = .false.
474 END IF
475*
476*+ TEST 1
477* Reconstruct matrix from factors and compute residual.
478*
479 CALL chet01_3( uplo, n, a, lda, afac, lda, e, iwork,
480 $ ainv, lda, rwork, result( 1 ) )
481 nt = 1
482*
483*+ TEST 2
484* Form the inverse and compute the residual,
485* if the factorization was competed without INFO > 0
486* (i.e. there is no zero rows and columns).
487* Do it only for the first block size.
488*
489 IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
490 CALL clacpy( uplo, n, n, afac, lda, ainv, lda )
491 srnamt = 'CHETRI_3'
492*
493* Another reason that we need to compute the inverse
494* is that CPOT03 produces RCONDC which is used later
495* in TEST6 and TEST7.
496*
497 lwork = (n+nb+1)*(nb+3)
498 CALL chetri_3( uplo, n, ainv, lda, e, iwork, work,
499 $ lwork, info )
500*
501* Check error code from ZHETRI_3 and handle error.
502*
503 IF( info.NE.0 )
504 $ CALL alaerh( path, 'CHETRI_3', info, -1,
505 $ uplo, n, n, -1, -1, -1, imat,
506 $ nfail, nerrs, nout )
507*
508* Compute the residual for a Hermitian matrix times
509* its inverse.
510*
511 CALL cpot03( uplo, n, a, lda, ainv, lda, work, lda,
512 $ rwork, rcondc, result( 2 ) )
513 nt = 2
514 END IF
515*
516* Print information about the tests that did not pass
517* the threshold.
518*
519 DO 110 k = 1, nt
520 IF( result( k ).GE.thresh ) THEN
521 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
522 $ CALL alahd( nout, path )
523 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
524 $ result( k )
525 nfail = nfail + 1
526 END IF
527 110 CONTINUE
528 nrun = nrun + nt
529*
530*+ TEST 3
531* Compute largest element in U or L
532*
533 result( 3 ) = zero
534 stemp = zero
535*
536 const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) ) /
537 $ ( one-alpha )
538*
539 IF( iuplo.EQ.1 ) THEN
540*
541* Compute largest element in U
542*
543 k = n
544 120 CONTINUE
545 IF( k.LE.1 )
546 $ GO TO 130
547*
548 IF( iwork( k ).GT.zero ) THEN
549*
550* Get max absolute value from elements
551* in column k in U
552*
553 stemp = clange( 'M', k-1, 1,
554 $ afac( ( k-1 )*lda+1 ), lda, rwork )
555 ELSE
556*
557* Get max absolute value from elements
558* in columns k and k-1 in U
559*
560 stemp = clange( 'M', k-2, 2,
561 $ afac( ( k-2 )*lda+1 ), lda, rwork )
562 k = k - 1
563*
564 END IF
565*
566* STEMP should be bounded by CONST
567*
568 stemp = stemp - const + thresh
569 IF( stemp.GT.result( 3 ) )
570 $ result( 3 ) = stemp
571*
572 k = k - 1
573*
574 GO TO 120
575 130 CONTINUE
576*
577 ELSE
578*
579* Compute largest element in L
580*
581 k = 1
582 140 CONTINUE
583 IF( k.GE.n )
584 $ GO TO 150
585*
586 IF( iwork( k ).GT.zero ) THEN
587*
588* Get max absolute value from elements
589* in column k in L
590*
591 stemp = clange( 'M', n-k, 1,
592 $ afac( ( k-1 )*lda+k+1 ), lda, rwork )
593 ELSE
594*
595* Get max absolute value from elements
596* in columns k and k+1 in L
597*
598 stemp = clange( 'M', n-k-1, 2,
599 $ afac( ( k-1 )*lda+k+2 ), lda, rwork )
600 k = k + 1
601*
602 END IF
603*
604* STEMP should be bounded by CONST
605*
606 stemp = stemp - const + thresh
607 IF( stemp.GT.result( 3 ) )
608 $ result( 3 ) = stemp
609*
610 k = k + 1
611*
612 GO TO 140
613 150 CONTINUE
614 END IF
615*
616*
617*+ TEST 4
618* Compute largest 2-Norm (condition number)
619* of 2-by-2 diag blocks
620*
621 result( 4 ) = zero
622 stemp = zero
623*
624 const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) )*
625 $ ( ( one + alpha ) / ( one - alpha ) )
626 CALL clacpy( uplo, n, n, afac, lda, ainv, lda )
627*
628 IF( iuplo.EQ.1 ) THEN
629*
630* Loop backward for UPLO = 'U'
631*
632 k = n
633 160 CONTINUE
634 IF( k.LE.1 )
635 $ GO TO 170
636*
637 IF( iwork( k ).LT.zero ) THEN
638*
639* Get the two singular values
640* (real and non-negative) of a 2-by-2 block,
641* store them in RWORK array
642*
643 block( 1, 1 ) = afac( ( k-2 )*lda+k-1 )
644 block( 1, 2 ) = e( k )
645 block( 2, 1 ) = conjg( block( 1, 2 ) )
646 block( 2, 2 ) = afac( (k-1)*lda+k )
647*
648 CALL cgesvd( 'N', 'N', 2, 2, block, 2, rwork,
649 $ cdummy, 1, cdummy, 1,
650 $ work, 6, rwork( 3 ), info )
651*
652*
653 sing_max = rwork( 1 )
654 sing_min = rwork( 2 )
655*
656 stemp = sing_max / sing_min
657*
658* STEMP should be bounded by CONST
659*
660 stemp = stemp - const + thresh
661 IF( stemp.GT.result( 4 ) )
662 $ result( 4 ) = stemp
663 k = k - 1
664*
665 END IF
666*
667 k = k - 1
668*
669 GO TO 160
670 170 CONTINUE
671*
672 ELSE
673*
674* Loop forward for UPLO = 'L'
675*
676 k = 1
677 180 CONTINUE
678 IF( k.GE.n )
679 $ GO TO 190
680*
681 IF( iwork( k ).LT.zero ) THEN
682*
683* Get the two singular values
684* (real and non-negative) of a 2-by-2 block,
685* store them in RWORK array
686*
687 block( 1, 1 ) = afac( ( k-1 )*lda+k )
688 block( 2, 1 ) = e( k )
689 block( 1, 2 ) = conjg( block( 2, 1 ) )
690 block( 2, 2 ) = afac( k*lda+k+1 )
691*
692 CALL cgesvd( 'N', 'N', 2, 2, block, 2, rwork,
693 $ cdummy, 1, cdummy, 1,
694 $ work, 6, rwork(3), info )
695*
696 sing_max = rwork( 1 )
697 sing_min = rwork( 2 )
698*
699 stemp = sing_max / sing_min
700*
701* STEMP should be bounded by CONST
702*
703 stemp = stemp - const + thresh
704 IF( stemp.GT.result( 4 ) )
705 $ result( 4 ) = stemp
706 k = k + 1
707*
708 END IF
709*
710 k = k + 1
711*
712 GO TO 180
713 190 CONTINUE
714 END IF
715*
716* Print information about the tests that did not pass
717* the threshold.
718*
719 DO 200 k = 3, 4
720 IF( result( k ).GE.thresh ) THEN
721 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
722 $ CALL alahd( nout, path )
723 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
724 $ result( k )
725 nfail = nfail + 1
726 END IF
727 200 CONTINUE
728 nrun = nrun + 2
729*
730* Skip the other tests if this is not the first block
731* size.
732*
733 IF( inb.GT.1 )
734 $ GO TO 240
735*
736* Do only the condition estimate if INFO is not 0.
737*
738 IF( trfcon ) THEN
739 rcondc = zero
740 GO TO 230
741 END IF
742*
743* Do for each value of NRHS in NSVAL.
744*
745 DO 220 irhs = 1, nns
746 nrhs = nsval( irhs )
747*
748* Begin loop over NRHS values
749*
750*
751*+ TEST 5 ( Using TRS_3)
752* Solve and compute residual for A * X = B.
753*
754* Choose a set of NRHS random solution vectors
755* stored in XACT and set up the right hand side B
756*
757 srnamt = 'CLARHS'
758 CALL clarhs( matpath, xtype, uplo, ' ', n, n,
759 $ kl, ku, nrhs, a, lda, xact, lda,
760 $ b, lda, iseed, info )
761 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
762*
763 srnamt = 'CHETRS_3'
764 CALL chetrs_3( uplo, n, nrhs, afac, lda, e, iwork,
765 $ x, lda, info )
766*
767* Check error code from CHETRS_3 and handle error.
768*
769 IF( info.NE.0 )
770 $ CALL alaerh( path, 'CHETRS_3', info, 0,
771 $ uplo, n, n, -1, -1, nrhs, imat,
772 $ nfail, nerrs, nout )
773*
774 CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
775*
776* Compute the residual for the solution
777*
778 CALL cpot02( uplo, n, nrhs, a, lda, x, lda, work,
779 $ lda, rwork, result( 5 ) )
780*
781*+ TEST 6
782* Check solution from generated exact solution.
783*
784 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
785 $ result( 6 ) )
786*
787* Print information about the tests that did not pass
788* the threshold.
789*
790 DO 210 k = 5, 6
791 IF( result( k ).GE.thresh ) THEN
792 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
793 $ CALL alahd( nout, path )
794 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
795 $ imat, k, result( k )
796 nfail = nfail + 1
797 END IF
798 210 CONTINUE
799 nrun = nrun + 2
800*
801* End do for each value of NRHS in NSVAL.
802*
803 220 CONTINUE
804*
805*+ TEST 7
806* Get an estimate of RCOND = 1/CNDNUM.
807*
808 230 CONTINUE
809 anorm = clanhe( '1', uplo, n, a, lda, rwork )
810 srnamt = 'CHECON_3'
811 CALL checon_3( uplo, n, afac, lda, e, iwork, anorm,
812 $ rcond, work, info )
813*
814* Check error code from CHECON_3 and handle error.
815*
816 IF( info.NE.0 )
817 $ CALL alaerh( path, 'CHECON_3', info, 0,
818 $ uplo, n, n, -1, -1, -1, imat,
819 $ nfail, nerrs, nout )
820*
821* Compute the test ratio to compare values of RCOND
822*
823 result( 7 ) = sget06( rcond, rcondc )
824*
825* Print information about the tests that did not pass
826* the threshold.
827*
828 IF( result( 7 ).GE.thresh ) THEN
829 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
830 $ CALL alahd( nout, path )
831 WRITE( nout, fmt = 9997 )uplo, n, imat, 7,
832 $ result( 7 )
833 nfail = nfail + 1
834 END IF
835 nrun = nrun + 1
836 240 CONTINUE
837*
838 250 CONTINUE
839 260 CONTINUE
840 270 CONTINUE
841*
842* Print a summary of the results.
843*
844 CALL alasum( path, nout, nfail, nrun, nerrs )
845*
846 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
847 $ i2, ', test ', i2, ', ratio =', g12.5 )
848 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
849 $ i2, ', test ', i2, ', ratio =', g12.5 )
850 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
851 $ ', test ', i2, ', ratio =', g12.5 )
852 RETURN
853*
854* End of CCHKHE_RK
855*
856 END
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 cchkhe_rk(dotype, nn, nval, nnb, nbval, nns, nsval, thresh, tsterr, nmax, a, afac, e, ainv, b, x, xact, work, rwork, iwork, nout)
CCHKHE_RK
Definition cchkhe_rk.f:177
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_3(uplo, n, a, lda, afac, ldafac, e, ipiv, c, ldc, rwork, resid)
CHET01_3
Definition chet01_3.f:141
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 cgesvd(jobu, jobvt, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, info)
CGESVD computes the singular value decomposition (SVD) for GE matrices
Definition cgesvd.f:214
subroutine checon_3(uplo, n, a, lda, e, ipiv, anorm, rcond, work, info)
CHECON_3
Definition checon_3.f:166
subroutine chetrf_rk(uplo, n, a, lda, e, ipiv, work, lwork, info)
CHETRF_RK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bunch...
Definition chetrf_rk.f:259
subroutine chetri_3(uplo, n, a, lda, e, ipiv, work, lwork, info)
CHETRI_3
Definition chetri_3.f:170
subroutine chetrs_3(uplo, n, nrhs, a, lda, e, ipiv, b, ldb, info)
CHETRS_3
Definition chetrs_3.f:165
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