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