LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine cchksy_rk ( 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( * )  E,
complex, dimension( * )  AINV,
complex, dimension( * )  B,
complex, dimension( * )  X,
complex, dimension( * )  XACT,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

CCHKSY_RK

Purpose:
 CCHKSY_RK tests CSYTRF_RK, -TRI_3, -TRS_3,
 and -CON_3.
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 (NBVAL)
          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]E
          E is COMPLEX array, dimension (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 (2*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.
Date
December 2016

Definition at line 179 of file cchksy_rk.f.

179 *
180 * -- LAPACK test routine (version 3.7.0) --
181 * -- LAPACK is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * December 2016
184 *
185 * .. Scalar Arguments ..
186  LOGICAL tsterr
187  INTEGER nmax, nn, nnb, nns, nout
188  REAL thresh
189 * ..
190 * .. Array Arguments ..
191  LOGICAL dotype( * )
192  INTEGER iwork( * ), nbval( * ), nsval( * ), nval( * )
193  REAL rwork( * )
194  COMPLEX a( * ), afac( * ), ainv( * ), b( * ), e( * ),
195  $ work( * ), x( * ), xact( * )
196 * ..
197 *
198 * =====================================================================
199 *
200 * .. Parameters ..
201  REAL zero, one
202  parameter ( zero = 0.0e+0, one = 1.0e+0 )
203  REAL onehalf
204  parameter ( onehalf = 0.5e+0 )
205  REAL eight, sevten
206  parameter ( eight = 8.0e+0, sevten = 17.0e+0 )
207  COMPLEX czero
208  parameter ( czero = ( 0.0e+0, 0.0e+0 ) )
209  INTEGER ntypes
210  parameter ( ntypes = 11 )
211  INTEGER ntests
212  parameter ( ntests = 7 )
213 * ..
214 * .. Local Scalars ..
215  LOGICAL trfcon, zerot
216  CHARACTER dist, TYPE, uplo, xtype
217  CHARACTER*3 path, matpath
218  INTEGER i, i1, i2, imat, in, inb, info, ioff, irhs,
219  $ itemp, itemp2, iuplo, izero, j, k, kl, ku, lda,
220  $ lwork, mode, n, nb, nerrs, nfail, nimat, nrhs,
221  $ nrun, nt
222  REAL alpha, anorm, cndnum, const, sing_max,
223  $ sing_min, rcond, rcondc, stemp
224 * ..
225 * .. Local Arrays ..
226  CHARACTER uplos( 2 )
227  INTEGER iseed( 4 ), iseedy( 4 )
228  REAL result( ntests )
229  COMPLEX block( 2, 2 ), cdummy( 1 )
230 * ..
231 * .. External Functions ..
232  REAL clange, clansy, sget06
233  EXTERNAL clange, clansy, sget06
234 * ..
235 * .. External Subroutines ..
236  EXTERNAL alaerh, alahd, alasum, cerrsy, cgesvd, cget04,
240 * ..
241 * .. Intrinsic Functions ..
242  INTRINSIC max, min, sqrt
243 * ..
244 * .. Scalars in Common ..
245  LOGICAL lerr, ok
246  CHARACTER*32 srnamt
247  INTEGER infot, nunit
248 * ..
249 * .. Common blocks ..
250  COMMON / infoc / infot, nunit, ok, lerr
251  COMMON / srnamc / srnamt
252 * ..
253 * .. Data statements ..
254  DATA iseedy / 1988, 1989, 1990, 1991 /
255  DATA uplos / 'U', 'L' /
256 * ..
257 * .. Executable Statements ..
258 *
259 * Initialize constants and the random number seed.
260 *
261  alpha = ( one+sqrt( sevten ) ) / eight
262 *
263 * Test path
264 *
265  path( 1: 1 ) = 'Complex precision'
266  path( 2: 3 ) = 'SK'
267 *
268 * Path to generate matrices
269 *
270  matpath( 1: 1 ) = 'Complex precision'
271  matpath( 2: 3 ) = 'SY'
272 *
273  nrun = 0
274  nfail = 0
275  nerrs = 0
276  DO 10 i = 1, 4
277  iseed( i ) = iseedy( i )
278  10 CONTINUE
279 *
280 * Test the error exits
281 *
282  IF( tsterr )
283  $ CALL cerrsy( path, nout )
284  infot = 0
285 *
286 * Set the minimum block size for which the block routine should
287 * be used, which will be later returned by ILAENV
288 *
289  CALL xlaenv( 2, 2 )
290 *
291 * Do for each value of N in NVAL
292 *
293  DO 270 in = 1, nn
294  n = nval( in )
295  lda = max( n, 1 )
296  xtype = 'N'
297  nimat = ntypes
298  IF( n.LE.0 )
299  $ nimat = 1
300 *
301  izero = 0
302 *
303 * Do for each value of matrix type IMAT
304 *
305  DO 260 imat = 1, nimat
306 *
307 * Do the tests only if DOTYPE( IMAT ) is true.
308 *
309  IF( .NOT.dotype( imat ) )
310  $ GO TO 260
311 *
312 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
313 *
314  zerot = imat.GE.3 .AND. imat.LE.6
315  IF( zerot .AND. n.LT.imat-2 )
316  $ GO TO 260
317 *
318 * Do first for UPLO = 'U', then for UPLO = 'L'
319 *
320  DO 250 iuplo = 1, 2
321  uplo = uplos( iuplo )
322 *
323 * Begin generate test matrix A.
324 *
325  IF( imat.NE.ntypes ) THEN
326 *
327 * Set up parameters with CLATB4 for the matrix generator
328 * based on the type of matrix to be generated.
329 *
330  CALL clatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
331  $ mode, cndnum, dist )
332 *
333 * Generate a matrix with CLATMS.
334 *
335  srnamt = 'CLATMS'
336  CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
337  $ cndnum, anorm, kl, ku, uplo, a, lda,
338  $ work, info )
339 *
340 * Check error code from CLATMS and handle error.
341 *
342  IF( info.NE.0 ) THEN
343  CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n,
344  $ -1, -1, -1, imat, nfail, nerrs, nout )
345 *
346 * Skip all tests for this generated matrix
347 *
348  GO TO 250
349  END IF
350 *
351 * For matrix types 3-6, zero one or more rows and
352 * columns of the matrix to test that INFO is returned
353 * correctly.
354 *
355  IF( zerot ) THEN
356  IF( imat.EQ.3 ) THEN
357  izero = 1
358  ELSE IF( imat.EQ.4 ) THEN
359  izero = n
360  ELSE
361  izero = n / 2 + 1
362  END IF
363 *
364  IF( imat.LT.6 ) THEN
365 *
366 * Set row and column IZERO to zero.
367 *
368  IF( iuplo.EQ.1 ) THEN
369  ioff = ( izero-1 )*lda
370  DO 20 i = 1, izero - 1
371  a( ioff+i ) = czero
372  20 CONTINUE
373  ioff = ioff + izero
374  DO 30 i = izero, n
375  a( ioff ) = czero
376  ioff = ioff + lda
377  30 CONTINUE
378  ELSE
379  ioff = izero
380  DO 40 i = 1, izero - 1
381  a( ioff ) = czero
382  ioff = ioff + lda
383  40 CONTINUE
384  ioff = ioff - izero
385  DO 50 i = izero, n
386  a( ioff+i ) = czero
387  50 CONTINUE
388  END IF
389  ELSE
390  IF( iuplo.EQ.1 ) THEN
391 *
392 * Set the first IZERO rows and columns to zero.
393 *
394  ioff = 0
395  DO 70 j = 1, n
396  i2 = min( j, izero )
397  DO 60 i = 1, i2
398  a( ioff+i ) = czero
399  60 CONTINUE
400  ioff = ioff + lda
401  70 CONTINUE
402  ELSE
403 *
404 * Set the last IZERO rows and columns to zero.
405 *
406  ioff = 0
407  DO 90 j = 1, n
408  i1 = max( j, izero )
409  DO 80 i = i1, n
410  a( ioff+i ) = czero
411  80 CONTINUE
412  ioff = ioff + lda
413  90 CONTINUE
414  END IF
415  END IF
416  ELSE
417  izero = 0
418  END IF
419 *
420  ELSE
421 *
422 * For matrix kind IMAT = 11, generate special block
423 * diagonal matrix to test alternate code
424 * for the 2 x 2 blocks.
425 *
426  CALL clatsy( uplo, n, a, lda, iseed )
427 *
428  END IF
429 *
430 * End generate test matrix A.
431 *
432 *
433 * Do for each value of NB in NBVAL
434 *
435  DO 240 inb = 1, nnb
436 *
437 * Set the optimal blocksize, which will be later
438 * returned by ILAENV.
439 *
440  nb = nbval( inb )
441  CALL xlaenv( 1, nb )
442 *
443 * Copy the test matrix A into matrix AFAC which
444 * will be factorized in place. This is needed to
445 * preserve the test matrix A for subsequent tests.
446 *
447  CALL clacpy( uplo, n, n, a, lda, afac, lda )
448 *
449 * Compute the L*D*L**T or U*D*U**T factorization of the
450 * matrix. IWORK stores details of the interchanges and
451 * the block structure of D. AINV is a work array for
452 * block factorization, LWORK is the length of AINV.
453 *
454  lwork = max( 2, nb )*lda
455  srnamt = 'CSYTRF_RK'
456  CALL csytrf_rk( uplo, n, afac, lda, e, iwork, ainv,
457  $ lwork, info )
458 *
459 * Adjust the expected value of INFO to account for
460 * pivoting.
461 *
462  k = izero
463  IF( k.GT.0 ) THEN
464  100 CONTINUE
465  IF( iwork( k ).LT.0 ) THEN
466  IF( iwork( k ).NE.-k ) THEN
467  k = -iwork( k )
468  GO TO 100
469  END IF
470  ELSE IF( iwork( k ).NE.k ) THEN
471  k = iwork( k )
472  GO TO 100
473  END IF
474  END IF
475 *
476 * Check error code from CSYTRF_RK and handle error.
477 *
478  IF( info.NE.k)
479  $ CALL alaerh( path, 'CSYTRF_RK', info, k,
480  $ uplo, n, n, -1, -1, nb, imat,
481  $ nfail, nerrs, nout )
482 *
483 * Set the condition estimate flag if the INFO is not 0.
484 *
485  IF( info.NE.0 ) THEN
486  trfcon = .true.
487  ELSE
488  trfcon = .false.
489  END IF
490 *
491 *+ TEST 1
492 * Reconstruct matrix from factors and compute residual.
493 *
494  CALL csyt01_3( uplo, n, a, lda, afac, lda, e, iwork,
495  $ ainv, lda, rwork, result( 1 ) )
496  nt = 1
497 *
498 *+ TEST 2
499 * Form the inverse and compute the residual,
500 * if the factorization was competed without INFO > 0
501 * (i.e. there is no zero rows and columns).
502 * Do it only for the first block size.
503 *
504  IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
505  CALL clacpy( uplo, n, n, afac, lda, ainv, lda )
506  srnamt = 'CSYTRI_3'
507 *
508 * Another reason that we need to compute the invesrse
509 * is that CSYT03 produces RCONDC which is used later
510 * in TEST6 and TEST7.
511 *
512  lwork = (n+nb+1)*(nb+3)
513  CALL csytri_3( uplo, n, ainv, lda, e, iwork, work,
514  $ lwork, info )
515 *
516 * Check error code from CSYTRI_3 and handle error.
517 *
518  IF( info.NE.0 )
519  $ CALL alaerh( path, 'CSYTRI_3', info, -1,
520  $ uplo, n, n, -1, -1, -1, imat,
521  $ nfail, nerrs, nout )
522 *
523 * Compute the residual for a symmetric matrix times
524 * its inverse.
525 *
526  CALL csyt03( uplo, n, a, lda, ainv, lda, work, lda,
527  $ rwork, rcondc, result( 2 ) )
528  nt = 2
529  END IF
530 *
531 * Print information about the tests that did not pass
532 * the threshold.
533 *
534  DO 110 k = 1, nt
535  IF( result( k ).GE.thresh ) THEN
536  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
537  $ CALL alahd( nout, path )
538  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
539  $ result( k )
540  nfail = nfail + 1
541  END IF
542  110 CONTINUE
543  nrun = nrun + nt
544 *
545 *+ TEST 3
546 * Compute largest element in U or L
547 *
548  result( 3 ) = zero
549  stemp = zero
550 *
551  const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) ) /
552  $ ( one-alpha )
553 *
554  IF( iuplo.EQ.1 ) THEN
555 *
556 * Compute largest element in U
557 *
558  k = n
559  120 CONTINUE
560  IF( k.LE.1 )
561  $ GO TO 130
562 *
563  IF( iwork( k ).GT.zero ) THEN
564 *
565 * Get max absolute value from elements
566 * in column k in in U
567 *
568  stemp = clange( 'M', k-1, 1,
569  $ afac( ( k-1 )*lda+1 ), lda, rwork )
570  ELSE
571 *
572 * Get max absolute value from elements
573 * in columns k and k-1 in U
574 *
575  stemp = clange( 'M', k-2, 2,
576  $ afac( ( k-2 )*lda+1 ), lda, rwork )
577  k = k - 1
578 *
579  END IF
580 *
581 * STEMP should be bounded by CONST
582 *
583  stemp = stemp - const + thresh
584  IF( stemp.GT.result( 3 ) )
585  $ result( 3 ) = stemp
586 *
587  k = k - 1
588 *
589  GO TO 120
590  130 CONTINUE
591 *
592  ELSE
593 *
594 * Compute largest element in L
595 *
596  k = 1
597  140 CONTINUE
598  IF( k.GE.n )
599  $ GO TO 150
600 *
601  IF( iwork( k ).GT.zero ) THEN
602 *
603 * Get max absolute value from elements
604 * in column k in in L
605 *
606  stemp = clange( 'M', n-k, 1,
607  $ afac( ( k-1 )*lda+k+1 ), lda, rwork )
608  ELSE
609 *
610 * Get max absolute value from elements
611 * in columns k and k+1 in L
612 *
613  stemp = clange( 'M', n-k-1, 2,
614  $ afac( ( k-1 )*lda+k+2 ), lda, rwork )
615  k = k + 1
616 *
617  END IF
618 *
619 * STEMP should be bounded by CONST
620 *
621  stemp = stemp - const + thresh
622  IF( stemp.GT.result( 3 ) )
623  $ result( 3 ) = stemp
624 *
625  k = k + 1
626 *
627  GO TO 140
628  150 CONTINUE
629  END IF
630 *
631 *
632 *+ TEST 4
633 * Compute largest 2-Norm (condition number)
634 * of 2-by-2 diag blocks
635 *
636  result( 4 ) = zero
637  stemp = zero
638 *
639  const = ( ( alpha**2-one ) / ( alpha**2-onehalf ) )*
640  $ ( ( one + alpha ) / ( one - alpha ) )
641 *
642  IF( iuplo.EQ.1 ) THEN
643 *
644 * Loop backward for UPLO = 'U'
645 *
646  k = n
647  160 CONTINUE
648  IF( k.LE.1 )
649  $ GO TO 170
650 *
651  IF( iwork( k ).LT.zero ) THEN
652 *
653 * Get the two singular values
654 * (real and non-negative) of a 2-by-2 block,
655 * store them in RWORK array
656 *
657  block( 1, 1 ) = afac( ( k-2 )*lda+k-1 )
658  block( 1, 2 ) = e( k )
659  block( 2, 1 ) = block( 1, 2 )
660  block( 2, 2 ) = afac( (k-1)*lda+k )
661 *
662  CALL cgesvd( 'N', 'N', 2, 2, block, 2, rwork,
663  $ cdummy, 1, cdummy, 1,
664  $ work, 6, rwork( 3 ), info )
665 *
666 *
667  sing_max = rwork( 1 )
668  sing_min = rwork( 2 )
669 *
670  stemp = sing_max / sing_min
671 *
672 * STEMP should be bounded by CONST
673 *
674  stemp = stemp - const + thresh
675  IF( stemp.GT.result( 4 ) )
676  $ result( 4 ) = stemp
677  k = k - 1
678 *
679  END IF
680 *
681  k = k - 1
682 *
683  GO TO 160
684  170 CONTINUE
685 *
686  ELSE
687 *
688 * Loop forward for UPLO = 'L'
689 *
690  k = 1
691  180 CONTINUE
692  IF( k.GE.n )
693  $ GO TO 190
694 *
695  IF( iwork( k ).LT.zero ) THEN
696 *
697 * Get the two singular values
698 * (real and non-negative) of a 2-by-2 block,
699 * store them in RWORK array
700 *
701  block( 1, 1 ) = afac( ( k-1 )*lda+k )
702  block( 2, 1 ) = e( k )
703  block( 1, 2 ) = block( 2, 1 )
704  block( 2, 2 ) = afac( k*lda+k+1 )
705 *
706  CALL cgesvd( 'N', 'N', 2, 2, block, 2, rwork,
707  $ cdummy, 1, cdummy, 1,
708  $ work, 6, rwork(3), info )
709 *
710  sing_max = rwork( 1 )
711  sing_min = rwork( 2 )
712 *
713  stemp = sing_max / sing_min
714 *
715 * STEMP should be bounded by CONST
716 *
717  stemp = stemp - const + thresh
718  IF( stemp.GT.result( 4 ) )
719  $ result( 4 ) = stemp
720  k = k + 1
721 *
722  END IF
723 *
724  k = k + 1
725 *
726  GO TO 180
727  190 CONTINUE
728  END IF
729 *
730 * Print information about the tests that did not pass
731 * the threshold.
732 *
733  DO 200 k = 3, 4
734  IF( result( k ).GE.thresh ) THEN
735  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
736  $ CALL alahd( nout, path )
737  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
738  $ result( k )
739  nfail = nfail + 1
740  END IF
741  200 CONTINUE
742  nrun = nrun + 2
743 *
744 * Skip the other tests if this is not the first block
745 * size.
746 *
747  IF( inb.GT.1 )
748  $ GO TO 240
749 *
750 * Do only the condition estimate if INFO is not 0.
751 *
752  IF( trfcon ) THEN
753  rcondc = zero
754  GO TO 230
755  END IF
756 *
757 * Do for each value of NRHS in NSVAL.
758 *
759  DO 220 irhs = 1, nns
760  nrhs = nsval( irhs )
761 *
762 *+ TEST 5 ( Using TRS_3)
763 * Solve and compute residual for A * X = B.
764 *
765 * Choose a set of NRHS random solution vectors
766 * stored in XACT and set up the right hand side B
767 *
768  srnamt = 'CLARHS'
769  CALL clarhs( matpath, xtype, uplo, ' ', n, n,
770  $ kl, ku, nrhs, a, lda, xact, lda,
771  $ b, lda, iseed, info )
772  CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
773 *
774  srnamt = 'CSYTRS_3'
775  CALL csytrs_3( uplo, n, nrhs, afac, lda, e, iwork,
776  $ x, lda, info )
777 *
778 * Check error code from CSYTRS_3 and handle error.
779 *
780  IF( info.NE.0 )
781  $ CALL alaerh( path, 'CSYTRS_3', info, 0,
782  $ uplo, n, n, -1, -1, nrhs, imat,
783  $ nfail, nerrs, nout )
784 *
785  CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
786 *
787 * Compute the residual for the solution
788 *
789  CALL csyt02( uplo, n, nrhs, a, lda, x, lda, work,
790  $ lda, rwork, result( 5 ) )
791 *
792 *+ TEST 6
793 * Check solution from generated exact solution.
794 *
795  CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
796  $ result( 6 ) )
797 *
798 * Print information about the tests that did not pass
799 * the threshold.
800 *
801  DO 210 k = 5, 6
802  IF( result( k ).GE.thresh ) THEN
803  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
804  $ CALL alahd( nout, path )
805  WRITE( nout, fmt = 9998 )uplo, n, nrhs,
806  $ imat, k, result( k )
807  nfail = nfail + 1
808  END IF
809  210 CONTINUE
810  nrun = nrun + 2
811 *
812 * End do for each value of NRHS in NSVAL.
813 *
814  220 CONTINUE
815 *
816 *+ TEST 7
817 * Get an estimate of RCOND = 1/CNDNUM.
818 *
819  230 CONTINUE
820  anorm = clansy( '1', uplo, n, a, lda, rwork )
821  srnamt = 'CSYCON_3'
822  CALL csycon_3( uplo, n, afac, lda, e, iwork, anorm,
823  $ rcond, work, info )
824 *
825 * Check error code from CSYCON_3 and handle error.
826 *
827  IF( info.NE.0 )
828  $ CALL alaerh( path, 'CSYCON_3', info, 0,
829  $ uplo, n, n, -1, -1, -1, imat,
830  $ nfail, nerrs, nout )
831 *
832 * Compute the test ratio to compare values of RCOND
833 *
834  result( 7 ) = sget06( rcond, rcondc )
835 *
836 * Print information about the tests that did not pass
837 * the threshold.
838 *
839  IF( result( 7 ).GE.thresh ) THEN
840  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
841  $ CALL alahd( nout, path )
842  WRITE( nout, fmt = 9997 )uplo, n, imat, 7,
843  $ result( 7 )
844  nfail = nfail + 1
845  END IF
846  nrun = nrun + 1
847  240 CONTINUE
848 *
849  250 CONTINUE
850  260 CONTINUE
851  270 CONTINUE
852 *
853 * Print a summary of the results.
854 *
855  CALL alasum( path, nout, nfail, nrun, nerrs )
856 *
857  9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
858  $ i2, ', test ', i2, ', ratio =', g12.5 )
859  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
860  $ i2, ', test(', i2, ') =', g12.5 )
861  9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
862  $ ', test(', i2, ') =', g12.5 )
863  RETURN
864 *
865 * End of CCHKSY_RK
866 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:107
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine clatsy(UPLO, N, X, LDX, ISEED)
CLATSY
Definition: clatsy.f:91
subroutine clarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
CLARHS
Definition: clarhs.f:211
subroutine csyt02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
CSYT02
Definition: csyt02.f:129
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:216
subroutine cerrsy(PATH, NUNIT)
CERRSY
Definition: cerrsy.f:57
subroutine csycon_3(UPLO, N, A, LDA, E, IPIV, ANORM, RCOND, WORK, INFO)
CSYCON_3
Definition: csycon_3.f:173
subroutine csytrf_rk(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
CSYTRF_RK computes the factorization of a complex symmetric indefinite matrix using the bounded Bunch...
Definition: csytrf_rk.f:261
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
real function clange(NORM, M, N, A, LDA, WORK)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: clange.f:117
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine csyt03(UPLO, N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
CSYT03
Definition: csyt03.f:128
subroutine clatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
CLATMS
Definition: clatms.f:334
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine csytrs_3(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, INFO)
CSYTRS_3
Definition: csytrs_3.f:167
subroutine csyt01_3(UPLO, N, A, LDA, AFAC, LDAFAC, E, IPIV, C, LDC, RWORK, RESID)
CSYT01_3
Definition: csyt01_3.f:143
subroutine csytri_3(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
CSYTRI_3
Definition: csytri_3.f:172
real function clansy(NORM, UPLO, N, A, LDA, WORK)
CLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex symmetric matrix.
Definition: clansy.f:125
subroutine clatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
CLATB4
Definition: clatb4.f:123
subroutine cget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
CGET04
Definition: cget04.f:104
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75

Here is the call graph for this function:

Here is the caller graph for this function: