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