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

SCHKSY_RK

Purpose:
 SCHKSY_RK tests SSYTRF_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 REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]E
          E is REAL array, dimension (NMAX)
[out]AINV
          AINV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NSMAX),
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is REAL array, dimension (NMAX*NSMAX),
          where NSMAX is the largest entry in NSVAL.
[out]XACT
          XACT is REAL array, dimension (NMAX*NSMAX),
          where NSMAX is the largest entry in NSVAL.
[out]WORK
          WORK is REAL 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 178 of file schksy_rk.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: