LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ schksy_rook()

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
December 2016

Definition at line 173 of file schksy_rook.f.

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