LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchksy_rk.f
Go to the documentation of this file.
1*> \brief \b ZCHKSY_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 ZCHKSY_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* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
23* DOUBLE PRECISION RWORK( * )
24* COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ), E( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZCHKSY_RK tests ZSYTRF_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 DOUBLE PRECISION
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*16 array, dimension (NMAX*NMAX)
109*> \endverbatim
110
111*> \param[out] AFAC
112*> \verbatim
113*> AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
114*> \endverbatim
115*>
116*> \param[out] E
117*> \verbatim
118*> E is COMPLEX*16 array, dimension (NMAX)
119*> \endverbatim
120*>
121*> \param[out] AINV
122*> \verbatim
123*> AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
124*> \endverbatim
125*>
126*> \param[out] B
127*> \verbatim
128*> B is COMPLEX*16 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*16 array, dimension (NMAX*NSMAX)
135*> \endverbatim
136*>
137*> \param[out] XACT
138*> \verbatim
139*> XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
140*> \endverbatim
141*>
142*> \param[out] WORK
143*> \verbatim
144*> WORK is COMPLEX*16 array, dimension (NMAX*max(3,NSMAX))
145*> \endverbatim
146*>
147*> \param[out] RWORK
148*> \verbatim
149*> RWORK is DOUBLE PRECISION 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 complex16_lin
172*
173* =====================================================================
174 SUBROUTINE zchksy_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 DOUBLE PRECISION THRESH
186* ..
187* .. Array Arguments ..
188 LOGICAL DOTYPE( * )
189 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
190 DOUBLE PRECISION RWORK( * )
191 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ), E( * ),
192 $ work( * ), x( * ), xact( * )
193* ..
194*
195* =====================================================================
196*
197* .. Parameters ..
198 DOUBLE PRECISION ZERO, ONE
199 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
200 DOUBLE PRECISION ONEHALF
201 parameter( onehalf = 0.5d+0 )
202 DOUBLE PRECISION EIGHT, SEVTEN
203 parameter( eight = 8.0d+0, sevten = 17.0d+0 )
204 COMPLEX*16 CZERO
205 parameter( czero = ( 0.0d+0, 0.0d+0 ) )
206 INTEGER NTYPES
207 parameter( ntypes = 11 )
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 DOUBLE PRECISION ALPHA, ANORM, CNDNUM, CONST, DTEMP, SING_MAX,
220 $ SING_MIN, RCOND, RCONDC
221* ..
222* .. Local Arrays ..
223 CHARACTER UPLOS( 2 )
224 INTEGER ISEED( 4 ), ISEEDY( 4 )
225 DOUBLE PRECISION RESULT( NTESTS )
226 COMPLEX*16 BLOCK( 2, 2 ), ZDUMMY( 1 )
227* ..
228* .. External Functions ..
229 DOUBLE PRECISION DGET06, ZLANGE, ZLANSY
230 EXTERNAL DGET06, ZLANGE, ZLANSY
231* ..
232* .. External Subroutines ..
233 EXTERNAL alaerh, alahd, alasum, zerrsy, zgesvd, zget04,
237* ..
238* .. Intrinsic Functions ..
239 INTRINSIC 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 ) = 'Zomplex precision'
263 path( 2: 3 ) = 'SK'
264*
265* Path to generate matrices
266*
267 matpath( 1: 1 ) = 'Zomplex precision'
268 matpath( 2: 3 ) = 'SY'
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 zerrsy( 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 test matrix A.
321*
322 IF( imat.NE.ntypes ) THEN
323*
324* Set up parameters with ZLATB4 for the matrix generator
325* based on the type of matrix to be generated.
326*
327 CALL zlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
328 $ mode, cndnum, dist )
329*
330* Generate a matrix with ZLATMS.
331*
332 srnamt = 'ZLATMS'
333 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
334 $ cndnum, anorm, kl, ku, uplo, a, lda,
335 $ work, info )
336*
337* Check error code from ZLATMS and handle error.
338*
339 IF( info.NE.0 ) THEN
340 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
341 $ -1, -1, -1, imat, nfail, nerrs, nout )
342*
343* Skip all tests for this generated matrix
344*
345 GO TO 250
346 END IF
347*
348* For matrix types 3-6, zero one or more rows and
349* columns of the matrix to test that INFO is returned
350* correctly.
351*
352 IF( zerot ) THEN
353 IF( imat.EQ.3 ) THEN
354 izero = 1
355 ELSE IF( imat.EQ.4 ) THEN
356 izero = n
357 ELSE
358 izero = n / 2 + 1
359 END IF
360*
361 IF( imat.LT.6 ) THEN
362*
363* Set row and column IZERO to zero.
364*
365 IF( iuplo.EQ.1 ) THEN
366 ioff = ( izero-1 )*lda
367 DO 20 i = 1, izero - 1
368 a( ioff+i ) = czero
369 20 CONTINUE
370 ioff = ioff + izero
371 DO 30 i = izero, n
372 a( ioff ) = czero
373 ioff = ioff + lda
374 30 CONTINUE
375 ELSE
376 ioff = izero
377 DO 40 i = 1, izero - 1
378 a( ioff ) = czero
379 ioff = ioff + lda
380 40 CONTINUE
381 ioff = ioff - izero
382 DO 50 i = izero, n
383 a( ioff+i ) = czero
384 50 CONTINUE
385 END IF
386 ELSE
387 IF( iuplo.EQ.1 ) THEN
388*
389* Set the first IZERO rows and columns to zero.
390*
391 ioff = 0
392 DO 70 j = 1, n
393 i2 = min( j, izero )
394 DO 60 i = 1, i2
395 a( ioff+i ) = czero
396 60 CONTINUE
397 ioff = ioff + lda
398 70 CONTINUE
399 ELSE
400*
401* Set the last IZERO rows and columns to zero.
402*
403 ioff = 0
404 DO 90 j = 1, n
405 i1 = max( j, izero )
406 DO 80 i = i1, n
407 a( ioff+i ) = czero
408 80 CONTINUE
409 ioff = ioff + lda
410 90 CONTINUE
411 END IF
412 END IF
413 ELSE
414 izero = 0
415 END IF
416*
417 ELSE
418*
419* For matrix kind IMAT = 11, generate special block
420* diagonal matrix to test alternate code
421* for the 2 x 2 blocks.
422*
423 CALL zlatsy( uplo, n, a, lda, iseed )
424*
425 END IF
426*
427* End generate test matrix A.
428*
429*
430* Do for each value of NB in NBVAL
431*
432 DO 240 inb = 1, nnb
433*
434* Set the optimal blocksize, which will be later
435* returned by ILAENV.
436*
437 nb = nbval( inb )
438 CALL xlaenv( 1, nb )
439*
440* Copy the test matrix A into matrix AFAC which
441* will be factorized in place. This is needed to
442* preserve the test matrix A for subsequent tests.
443*
444 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
445*
446* Compute the L*D*L**T or U*D*U**T factorization of the
447* matrix. IWORK stores details of the interchanges and
448* the block structure of D. AINV is a work array for
449* block factorization, LWORK is the length of AINV.
450*
451 lwork = max( 2, nb )*lda
452 srnamt = 'ZSYTRF_RK'
453 CALL zsytrf_rk( uplo, n, afac, lda, e, iwork, ainv,
454 $ lwork, info )
455*
456* Adjust the expected value of INFO to account for
457* pivoting.
458*
459 k = izero
460 IF( k.GT.0 ) THEN
461 100 CONTINUE
462 IF( iwork( k ).LT.0 ) THEN
463 IF( iwork( k ).NE.-k ) THEN
464 k = -iwork( k )
465 GO TO 100
466 END IF
467 ELSE IF( iwork( k ).NE.k ) THEN
468 k = iwork( k )
469 GO TO 100
470 END IF
471 END IF
472*
473* Check error code from ZSYTRF_RK and handle error.
474*
475 IF( info.NE.k)
476 $ CALL alaerh( path, 'ZSYTRF_RK', info, k,
477 $ uplo, n, n, -1, -1, nb, imat,
478 $ nfail, nerrs, nout )
479*
480* Set the condition estimate flag if the INFO is not 0.
481*
482 IF( info.NE.0 ) THEN
483 trfcon = .true.
484 ELSE
485 trfcon = .false.
486 END IF
487*
488*+ TEST 1
489* Reconstruct matrix from factors and compute residual.
490*
491 CALL zsyt01_3( uplo, n, a, lda, afac, lda, e, iwork,
492 $ ainv, lda, rwork, result( 1 ) )
493 nt = 1
494*
495*+ TEST 2
496* Form the inverse and compute the residual,
497* if the factorization was competed without INFO > 0
498* (i.e. there is no zero rows and columns).
499* Do it only for the first block size.
500*
501 IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
502 CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
503 srnamt = 'ZSYTRI_3'
504*
505* Another reason that we need to compute the inverse
506* is that ZSYT03 produces RCONDC which is used later
507* in TEST6 and TEST7.
508*
509 lwork = (n+nb+1)*(nb+3)
510 CALL zsytri_3( uplo, n, ainv, lda, e, iwork, work,
511 $ lwork, info )
512*
513* Check error code from ZSYTRI_3 and handle error.
514*
515 IF( info.NE.0 )
516 $ CALL alaerh( path, 'ZSYTRI_3', info, -1,
517 $ uplo, n, n, -1, -1, -1, imat,
518 $ nfail, nerrs, nout )
519*
520* Compute the residual for a symmetric matrix times
521* its inverse.
522*
523 CALL zsyt03( uplo, n, a, lda, ainv, lda, work, lda,
524 $ rwork, rcondc, result( 2 ) )
525 nt = 2
526 END IF
527*
528* Print information about the tests that did not pass
529* the threshold.
530*
531 DO 110 k = 1, nt
532 IF( result( k ).GE.thresh ) THEN
533 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
534 $ CALL alahd( nout, path )
535 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
536 $ result( k )
537 nfail = nfail + 1
538 END IF
539 110 CONTINUE
540 nrun = nrun + nt
541*
542*+ TEST 3
543* Compute largest element in U or L
544*
545 result( 3 ) = zero
546 dtemp = zero
547*
548 const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) ) /
549 $ ( one-alpha )
550*
551 IF( iuplo.EQ.1 ) THEN
552*
553* Compute largest element in U
554*
555 k = n
556 120 CONTINUE
557 IF( k.LE.1 )
558 $ GO TO 130
559*
560 IF( iwork( k ).GT.zero ) THEN
561*
562* Get max absolute value from elements
563* in column k in in U
564*
565 dtemp = zlange( 'M', k-1, 1,
566 $ afac( ( k-1 )*lda+1 ), lda, rwork )
567 ELSE
568*
569* Get max absolute value from elements
570* in columns k and k-1 in U
571*
572 dtemp = zlange( 'M', k-2, 2,
573 $ afac( ( k-2 )*lda+1 ), lda, rwork )
574 k = k - 1
575*
576 END IF
577*
578* DTEMP should be bounded by CONST
579*
580 dtemp = dtemp - const + thresh
581 IF( dtemp.GT.result( 3 ) )
582 $ result( 3 ) = dtemp
583*
584 k = k - 1
585*
586 GO TO 120
587 130 CONTINUE
588*
589 ELSE
590*
591* Compute largest element in L
592*
593 k = 1
594 140 CONTINUE
595 IF( k.GE.n )
596 $ GO TO 150
597*
598 IF( iwork( k ).GT.zero ) THEN
599*
600* Get max absolute value from elements
601* in column k in in L
602*
603 dtemp = zlange( 'M', n-k, 1,
604 $ afac( ( k-1 )*lda+k+1 ), lda, rwork )
605 ELSE
606*
607* Get max absolute value from elements
608* in columns k and k+1 in L
609*
610 dtemp = zlange( 'M', n-k-1, 2,
611 $ afac( ( k-1 )*lda+k+2 ), lda, rwork )
612 k = k + 1
613*
614 END IF
615*
616* DTEMP should be bounded by CONST
617*
618 dtemp = dtemp - const + thresh
619 IF( dtemp.GT.result( 3 ) )
620 $ result( 3 ) = dtemp
621*
622 k = k + 1
623*
624 GO TO 140
625 150 CONTINUE
626 END IF
627*
628*
629*+ TEST 4
630* Compute largest 2-Norm (condition number)
631* of 2-by-2 diag blocks
632*
633 result( 4 ) = zero
634 dtemp = zero
635*
636 const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) )*
637 $ ( ( one + alpha ) / ( one - alpha ) )
638*
639 IF( iuplo.EQ.1 ) THEN
640*
641* Loop backward for UPLO = 'U'
642*
643 k = n
644 160 CONTINUE
645 IF( k.LE.1 )
646 $ GO TO 170
647*
648 IF( iwork( k ).LT.zero ) THEN
649*
650* Get the two singular values
651* (real and non-negative) of a 2-by-2 block,
652* store them in RWORK array
653*
654 block( 1, 1 ) = afac( ( k-2 )*lda+k-1 )
655 block( 1, 2 ) = e( k )
656 block( 2, 1 ) = block( 1, 2 )
657 block( 2, 2 ) = afac( (k-1)*lda+k )
658*
659 CALL zgesvd( 'N', 'N', 2, 2, block, 2, rwork,
660 $ zdummy, 1, zdummy, 1,
661 $ work, 6, rwork( 3 ), info )
662*
663*
664 sing_max = rwork( 1 )
665 sing_min = rwork( 2 )
666*
667 dtemp = sing_max / sing_min
668*
669* DTEMP should be bounded by CONST
670*
671 dtemp = dtemp - const + thresh
672 IF( dtemp.GT.result( 4 ) )
673 $ result( 4 ) = dtemp
674 k = k - 1
675*
676 END IF
677*
678 k = k - 1
679*
680 GO TO 160
681 170 CONTINUE
682*
683 ELSE
684*
685* Loop forward for UPLO = 'L'
686*
687 k = 1
688 180 CONTINUE
689 IF( k.GE.n )
690 $ GO TO 190
691*
692 IF( iwork( k ).LT.zero ) THEN
693*
694* Get the two singular values
695* (real and non-negative) of a 2-by-2 block,
696* store them in RWORK array
697*
698 block( 1, 1 ) = afac( ( k-1 )*lda+k )
699 block( 2, 1 ) = e( k )
700 block( 1, 2 ) = block( 2, 1 )
701 block( 2, 2 ) = afac( k*lda+k+1 )
702*
703 CALL zgesvd( 'N', 'N', 2, 2, block, 2, rwork,
704 $ zdummy, 1, zdummy, 1,
705 $ work, 6, rwork(3), info )
706*
707 sing_max = rwork( 1 )
708 sing_min = rwork( 2 )
709*
710 dtemp = sing_max / sing_min
711*
712* DTEMP should be bounded by CONST
713*
714 dtemp = dtemp - const + thresh
715 IF( dtemp.GT.result( 4 ) )
716 $ result( 4 ) = dtemp
717 k = k + 1
718*
719 END IF
720*
721 k = k + 1
722*
723 GO TO 180
724 190 CONTINUE
725 END IF
726*
727* Print information about the tests that did not pass
728* the threshold.
729*
730 DO 200 k = 3, 4
731 IF( result( k ).GE.thresh ) THEN
732 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
733 $ CALL alahd( nout, path )
734 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
735 $ result( k )
736 nfail = nfail + 1
737 END IF
738 200 CONTINUE
739 nrun = nrun + 2
740*
741* Skip the other tests if this is not the first block
742* size.
743*
744 IF( inb.GT.1 )
745 $ GO TO 240
746*
747* Do only the condition estimate if INFO is not 0.
748*
749 IF( trfcon ) THEN
750 rcondc = zero
751 GO TO 230
752 END IF
753*
754* Do for each value of NRHS in NSVAL.
755*
756 DO 220 irhs = 1, nns
757 nrhs = nsval( irhs )
758*
759*+ TEST 5 ( Using TRS_3)
760* Solve and compute residual for A * X = B.
761*
762* Choose a set of NRHS random solution vectors
763* stored in XACT and set up the right hand side B
764*
765 srnamt = 'ZLARHS'
766 CALL zlarhs( matpath, xtype, uplo, ' ', n, n,
767 $ kl, ku, nrhs, a, lda, xact, lda,
768 $ b, lda, iseed, info )
769 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
770*
771 srnamt = 'ZSYTRS_3'
772 CALL zsytrs_3( uplo, n, nrhs, afac, lda, e, iwork,
773 $ x, lda, info )
774*
775* Check error code from ZSYTRS_3 and handle error.
776*
777 IF( info.NE.0 )
778 $ CALL alaerh( path, 'ZSYTRS_3', info, 0,
779 $ uplo, n, n, -1, -1, nrhs, imat,
780 $ nfail, nerrs, nout )
781*
782 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
783*
784* Compute the residual for the solution
785*
786 CALL zsyt02( uplo, n, nrhs, a, lda, x, lda, work,
787 $ lda, rwork, result( 5 ) )
788*
789*+ TEST 6
790* Check solution from generated exact solution.
791*
792 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
793 $ result( 6 ) )
794*
795* Print information about the tests that did not pass
796* the threshold.
797*
798 DO 210 k = 5, 6
799 IF( result( k ).GE.thresh ) THEN
800 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
801 $ CALL alahd( nout, path )
802 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
803 $ imat, k, result( k )
804 nfail = nfail + 1
805 END IF
806 210 CONTINUE
807 nrun = nrun + 2
808*
809* End do for each value of NRHS in NSVAL.
810*
811 220 CONTINUE
812*
813*+ TEST 7
814* Get an estimate of RCOND = 1/CNDNUM.
815*
816 230 CONTINUE
817 anorm = zlansy( '1', uplo, n, a, lda, rwork )
818 srnamt = 'ZSYCON_3'
819 CALL zsycon_3( uplo, n, afac, lda, e, iwork, anorm,
820 $ rcond, work, info )
821*
822* Check error code from ZSYCON_3 and handle error.
823*
824 IF( info.NE.0 )
825 $ CALL alaerh( path, 'ZSYCON_3', info, 0,
826 $ uplo, n, n, -1, -1, -1, imat,
827 $ nfail, nerrs, nout )
828*
829* Compute the test ratio to compare values of RCOND
830*
831 result( 7 ) = dget06( rcond, rcondc )
832*
833* Print information about the tests that did not pass
834* the threshold.
835*
836 IF( result( 7 ).GE.thresh ) THEN
837 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
838 $ CALL alahd( nout, path )
839 WRITE( nout, fmt = 9997 )uplo, n, imat, 7,
840 $ result( 7 )
841 nfail = nfail + 1
842 END IF
843 nrun = nrun + 1
844 240 CONTINUE
845*
846 250 CONTINUE
847 260 CONTINUE
848 270 CONTINUE
849*
850* Print a summary of the results.
851*
852 CALL alasum( path, nout, nfail, nrun, nerrs )
853*
854 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
855 $ i2, ', test ', i2, ', ratio =', g12.5 )
856 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
857 $ i2, ', test(', i2, ') =', g12.5 )
858 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
859 $ ', test(', i2, ') =', g12.5 )
860 RETURN
861*
862* End of ZCHKSY_RK
863*
864 END
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
subroutine zgesvd(jobu, jobvt, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, info)
ZGESVD computes the singular value decomposition (SVD) for GE matrices
Definition zgesvd.f:214
subroutine zsycon_3(uplo, n, a, lda, e, ipiv, anorm, rcond, work, info)
ZSYCON_3
Definition zsycon_3.f:166
subroutine zsytrf_rk(uplo, n, a, lda, e, ipiv, work, lwork, info)
ZSYTRF_RK computes the factorization of a complex symmetric indefinite matrix using the bounded Bunch...
Definition zsytrf_rk.f:259
subroutine zsytri_3(uplo, n, a, lda, e, ipiv, work, lwork, info)
ZSYTRI_3
Definition zsytri_3.f:170
subroutine zsytrs_3(uplo, n, nrhs, a, lda, e, ipiv, b, ldb, info)
ZSYTRS_3
Definition zsytrs_3.f:165
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
subroutine zchksy_rk(dotype, nn, nval, nnb, nbval, nns, nsval, thresh, tsterr, nmax, a, afac, e, ainv, b, x, xact, work, rwork, iwork, nout)
ZCHKSY_RK
Definition zchksy_rk.f:177
subroutine zerrsy(path, nunit)
ZERRSY
Definition zerrsy.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
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 zlatsy(uplo, n, x, ldx, iseed)
ZLATSY
Definition zlatsy.f:89
subroutine zsyt01_3(uplo, n, a, lda, afac, ldafac, e, ipiv, c, ldc, rwork, resid)
ZSYT01_3
Definition zsyt01_3.f:141
subroutine zsyt02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZSYT02
Definition zsyt02.f:127
subroutine zsyt03(uplo, n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
ZSYT03
Definition zsyt03.f:126