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

SCHKSY_ROOK

Purpose:
 SCHKSY_ROOK tests SSYTRF_ROOK, -TRI_ROOK, -TRS_ROOK,
 and -CON_ROOK.
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]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)
[out]XACT
          XACT is REAL array, dimension (NMAX*NSMAX)
[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
November 2015

Definition at line 173 of file schksy_rook.f.

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