LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine dchksy_rook ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
double precision, dimension( * )  A,
double precision, dimension( * )  AFAC,
double precision, dimension( * )  AINV,
double precision, dimension( * )  B,
double precision, dimension( * )  X,
double precision, dimension( * )  XACT,
double precision, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

DCHKSY_ROOK

Purpose:
 DCHKSY_ROOK tests DSYTRF_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 DOUBLE PRECISION
          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 DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AINV
          AINV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION 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 dchksy_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  DOUBLE PRECISION thresh
183 * ..
184 * .. Array Arguments ..
185  LOGICAL dotype( * )
186  INTEGER iwork( * ), nbval( * ), nsval( * ), nval( * )
187  DOUBLE PRECISION a( * ), afac( * ), ainv( * ), b( * ),
188  $ rwork( * ), work( * ), x( * ), xact( * )
189 * ..
190 *
191 * =====================================================================
192 *
193 * .. Parameters ..
194  DOUBLE PRECISION zero, one
195  parameter ( zero = 0.0d+0, one = 1.0d+0 )
196  DOUBLE PRECISION 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  DOUBLE PRECISION alpha, anorm, cndnum, const, dtemp, sing_max,
211  $ sing_min, rcond, rcondc
212 * ..
213 * .. Local Arrays ..
214  CHARACTER uplos( 2 )
215  INTEGER iseed( 4 ), iseedy( 4 )
216  DOUBLE PRECISION block( 2, 2 ), ddummy( 1 ), result( ntests )
217 * ..
218 * .. External Functions ..
219  DOUBLE PRECISION dget06, dlange, dlansy
220  EXTERNAL dget06, dlange, dlansy
221 * ..
222 * .. External Subroutines ..
223  EXTERNAL alaerh, alahd, alasum, derrsy, dget04, dlacpy,
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 ) = 'Double precision'
253  path( 2: 3 ) = 'SR'
254 *
255 * Path to generate matrices
256 *
257  matpath( 1: 1 ) = 'Double 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 derrsy( 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 DLATB4 for the matrix generator
313 * based on the type of matrix to be generated.
314 *
315  CALL dlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
316  $ mode, cndnum, dist )
317 *
318 * Generate a matrix with DLATMS.
319 *
320  srnamt = 'DLATMS'
321  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
322  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
323  $ info )
324 *
325 * Check error code from DLATMS and handle error.
326 *
327  IF( info.NE.0 ) THEN
328  CALL alaerh( path, 'DLATMS', 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 dlacpy( 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 = 'DSYTRF_ROOK'
431  CALL dsytrf_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 DSYTRF_ROOK and handle error.
452 *
453  IF( info.NE.k)
454  $ CALL alaerh( path, 'DSYTRF_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 dsyt01_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 dlacpy( uplo, n, n, afac, lda, ainv, lda )
481  srnamt = 'DSYTRI_ROOK'
482  CALL dsytri_rook( uplo, n, ainv, lda, iwork, work,
483  $ info )
484 *
485 * Check error code from DSYTRI_ROOK and handle error.
486 *
487  IF( info.NE.0 )
488  $ CALL alaerh( path, 'DSYTRI_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 dpot03( 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  dtemp = 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  dtemp = dlange( '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  dtemp = dlange( 'M', k-2, 2,
544  $ afac( ( k-2 )*lda+1 ), lda, rwork )
545  k = k - 1
546 *
547  END IF
548 *
549 * DTEMP should be bounded by CONST
550 *
551  dtemp = dtemp - const + thresh
552  IF( dtemp.GT.result( 3 ) )
553  $ result( 3 ) = dtemp
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  dtemp = dlange( '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  dtemp = dlange( 'M', n-k-1, 2,
582  $ afac( ( k-1 )*lda+k+2 ), lda, rwork )
583  k = k + 1
584 *
585  END IF
586 *
587 * DTEMP should be bounded by CONST
588 *
589  dtemp = dtemp - const + thresh
590  IF( dtemp.GT.result( 3 ) )
591  $ result( 3 ) = dtemp
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  dtemp = zero
606 *
607  const = ( one+alpha ) / ( one-alpha )
608  CALL dlacpy( 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 dgesvd( 'N', 'N', 2, 2, block, 2, rwork,
631  $ ddummy, 1, ddummy, 1,
632  $ work, 10, info )
633 *
634  sing_max = rwork( 1 )
635  sing_min = rwork( 2 )
636 *
637  dtemp = sing_max / sing_min
638 *
639 * DTEMP should be bounded by CONST
640 *
641  dtemp = dtemp - const + thresh
642  IF( dtemp.GT.result( 4 ) )
643  $ result( 4 ) = dtemp
644  k = k - 1
645 *
646  END IF
647 *
648  k = k - 1
649 *
650  GO TO 160
651  170 CONTINUE
652 *
653  ELSE
654 *
655 * Loop forward for UPLO = 'L'
656 *
657  k = 1
658  180 CONTINUE
659  IF( k.GE.n )
660  $ GO TO 190
661 *
662  IF( iwork( k ).LT.zero ) THEN
663 *
664 * Get the two singular values
665 * (real and non-negative) of a 2-by-2 block,
666 * store them in RWORK array
667 *
668  block( 1, 1 ) = afac( ( k-1 )*lda+k )
669  block( 2, 1 ) = afac( ( k-1 )*lda+k+1 )
670  block( 1, 2 ) = block( 2, 1 )
671  block( 2, 2 ) = afac( k*lda+k+1 )
672 *
673  CALL dgesvd( 'N', 'N', 2, 2, block, 2, rwork,
674  $ ddummy, 1, ddummy, 1,
675  $ work, 10, info )
676 *
677 *
678  sing_max = rwork( 1 )
679  sing_min = rwork( 2 )
680 *
681  dtemp = sing_max / sing_min
682 *
683 * DTEMP should be bounded by CONST
684 *
685  dtemp = dtemp - const + thresh
686  IF( dtemp.GT.result( 4 ) )
687  $ result( 4 ) = dtemp
688  k = k + 1
689 *
690  END IF
691 *
692  k = k + 1
693 *
694  GO TO 180
695  190 CONTINUE
696  END IF
697 *
698 * Print information about the tests that did not pass
699 * the threshold.
700 *
701  DO 200 k = 3, 4
702  IF( result( k ).GE.thresh ) THEN
703  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
704  $ CALL alahd( nout, path )
705  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
706  $ result( k )
707  nfail = nfail + 1
708  END IF
709  200 CONTINUE
710  nrun = nrun + 2
711 *
712 * Skip the other tests if this is not the first block
713 * size.
714 *
715  IF( inb.GT.1 )
716  $ GO TO 240
717 *
718 * Do only the condition estimate if INFO is not 0.
719 *
720  IF( trfcon ) THEN
721  rcondc = zero
722  GO TO 230
723  END IF
724 *
725 * Do for each value of NRHS in NSVAL.
726 *
727  DO 220 irhs = 1, nns
728  nrhs = nsval( irhs )
729 *
730 *+ TEST 5 ( Using TRS_ROOK)
731 * Solve and compute residual for A * X = B.
732 *
733 * Choose a set of NRHS random solution vectors
734 * stored in XACT and set up the right hand side B
735 *
736  srnamt = 'DLARHS'
737  CALL dlarhs( matpath, xtype, uplo, ' ', n, n,
738  $ kl, ku, nrhs, a, lda, xact, lda,
739  $ b, lda, iseed, info )
740  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
741 *
742  srnamt = 'DSYTRS_ROOK'
743  CALL dsytrs_rook( uplo, n, nrhs, afac, lda, iwork,
744  $ x, lda, info )
745 *
746 * Check error code from DSYTRS_ROOK and handle error.
747 *
748  IF( info.NE.0 )
749  $ CALL alaerh( path, 'DSYTRS_ROOK', info, 0,
750  $ uplo, n, n, -1, -1, nrhs, imat,
751  $ nfail, nerrs, nout )
752 *
753  CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
754 *
755 * Compute the residual for the solution
756 *
757  CALL dpot02( uplo, n, nrhs, a, lda, x, lda, work,
758  $ lda, rwork, result( 5 ) )
759 *
760 *+ TEST 6
761 * Check solution from generated exact solution.
762 *
763  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
764  $ result( 6 ) )
765 *
766 * Print information about the tests that did not pass
767 * the threshold.
768 *
769  DO 210 k = 5, 6
770  IF( result( k ).GE.thresh ) THEN
771  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
772  $ CALL alahd( nout, path )
773  WRITE( nout, fmt = 9998 )uplo, n, nrhs,
774  $ imat, k, result( k )
775  nfail = nfail + 1
776  END IF
777  210 CONTINUE
778  nrun = nrun + 2
779 *
780 * End do for each value of NRHS in NSVAL.
781 *
782  220 CONTINUE
783 *
784 *+ TEST 7
785 * Get an estimate of RCOND = 1/CNDNUM.
786 *
787  230 CONTINUE
788  anorm = dlansy( '1', uplo, n, a, lda, rwork )
789  srnamt = 'DSYCON_ROOK'
790  CALL dsycon_rook( uplo, n, afac, lda, iwork, anorm,
791  $ rcond, work, iwork( n+1 ), info )
792 *
793 * Check error code from DSYCON_ROOK and handle error.
794 *
795  IF( info.NE.0 )
796  $ CALL alaerh( path, 'DSYCON_ROOK', info, 0,
797  $ uplo, n, n, -1, -1, -1, imat,
798  $ nfail, nerrs, nout )
799 *
800 * Compute the test ratio to compare to values of RCOND
801 *
802  result( 7 ) = dget06( rcond, rcondc )
803 *
804 * Print information about the tests that did not pass
805 * the threshold.
806 *
807  IF( result( 7 ).GE.thresh ) THEN
808  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
809  $ CALL alahd( nout, path )
810  WRITE( nout, fmt = 9997 )uplo, n, imat, 7,
811  $ result( 7 )
812  nfail = nfail + 1
813  END IF
814  nrun = nrun + 1
815  240 CONTINUE
816 *
817  250 CONTINUE
818  260 CONTINUE
819  270 CONTINUE
820 *
821 * Print a summary of the results.
822 *
823  CALL alasum( path, nout, nfail, nrun, nerrs )
824 *
825  9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
826  $ i2, ', test ', i2, ', ratio =', g12.5 )
827  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
828  $ i2, ', test(', i2, ') =', g12.5 )
829  9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
830  $ ', test(', i2, ') =', g12.5 )
831  RETURN
832 *
833 * End of DCHKSY_ROOK
834 *
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dsyt01_rook(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
DSYT01_ROOK
Definition: dsyt01_rook.f:126
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:107
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY 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: dlansy.f:124
subroutine dsytrs_rook(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
DSYTRS_ROOK
Definition: dsytrs_rook.f:138
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine dsytri_rook(UPLO, N, A, LDA, IPIV, WORK, INFO)
DSYTRI_ROOK
Definition: dsytri_rook.f:131
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine dsycon_rook(UPLO, N, A, LDA, IPIV, ANORM, RCOND, WORK, IWORK, INFO)
DSYCON_ROOK
Definition: dsycon_rook.f:146
double precision function dlange(NORM, M, N, A, LDA, WORK)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlange.f:116
subroutine dpot03(UPLO, N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
DPOT03
Definition: dpot03.f:127
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine derrsy(PATH, NUNIT)
DERRSY
Definition: derrsy.f:57
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dsytrf_rook(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
DSYTRF_ROOK
Definition: dsytrf_rook.f:210
subroutine dgesvd(JOBU, JOBVT, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, INFO)
DGESVD computes the singular value decomposition (SVD) for GE matrices
Definition: dgesvd.f:213
subroutine dpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DPOT02
Definition: dpot02.f:129
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: