LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine schkge ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
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 
)

SCHKGE

Purpose:
 SCHKGE tests SGETRF, -TRI, -TRS, -RFS, and -CON.
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]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[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 column 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 M or 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(2*NMAX,2*NSMAX+NWORK))
[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
April 2012

Definition at line 187 of file schkge.f.

187 *
188 * -- LAPACK test routine (version 3.4.1) --
189 * -- LAPACK is a software package provided by Univ. of Tennessee, --
190 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191 * April 2012
192 *
193 * .. Scalar Arguments ..
194  LOGICAL tsterr
195  INTEGER nm, nmax, nn, nnb, nns, nout
196  REAL thresh
197 * ..
198 * .. Array Arguments ..
199  LOGICAL dotype( * )
200  INTEGER iwork( * ), mval( * ), nbval( * ), nsval( * ),
201  $ nval( * )
202  REAL a( * ), afac( * ), ainv( * ), b( * ),
203  $ rwork( * ), work( * ), x( * ), xact( * )
204 * ..
205 *
206 * =====================================================================
207 *
208 * .. Parameters ..
209  REAL one, zero
210  parameter ( one = 1.0e+0, zero = 0.0e+0 )
211  INTEGER ntypes
212  parameter ( ntypes = 11 )
213  INTEGER ntests
214  parameter ( ntests = 8 )
215  INTEGER ntran
216  parameter ( ntran = 3 )
217 * ..
218 * .. Local Scalars ..
219  LOGICAL trfcon, zerot
220  CHARACTER dist, norm, trans, TYPE, xtype
221  CHARACTER*3 path
222  INTEGER i, im, imat, in, inb, info, ioff, irhs, itran,
223  $ izero, k, kl, ku, lda, lwork, m, mode, n, nb,
224  $ nerrs, nfail, nimat, nrhs, nrun, nt
225  REAL ainvnm, anorm, anormi, anormo, cndnum, dummy,
226  $ rcond, rcondc, rcondi, rcondo
227 * ..
228 * .. Local Arrays ..
229  CHARACTER transs( ntran )
230  INTEGER iseed( 4 ), iseedy( 4 )
231  REAL result( ntests )
232 * ..
233 * .. External Functions ..
234  REAL sget06, slange
235  EXTERNAL sget06, slange
236 * ..
237 * .. External Subroutines ..
238  EXTERNAL alaerh, alahd, alasum, serrge, sgecon, sgerfs,
241  $ slatms, xlaenv
242 * ..
243 * .. Intrinsic Functions ..
244  INTRINSIC max, min
245 * ..
246 * .. Scalars in Common ..
247  LOGICAL lerr, ok
248  CHARACTER*32 srnamt
249  INTEGER infot, nunit
250 * ..
251 * .. Common blocks ..
252  COMMON / infoc / infot, nunit, ok, lerr
253  COMMON / srnamc / srnamt
254 * ..
255 * .. Data statements ..
256  DATA iseedy / 1988, 1989, 1990, 1991 / ,
257  $ transs / 'N', 'T', 'C' /
258 * ..
259 * .. Executable Statements ..
260 *
261 * Initialize constants and the random number seed.
262 *
263  path( 1: 1 ) = 'Single precision'
264  path( 2: 3 ) = 'GE'
265  nrun = 0
266  nfail = 0
267  nerrs = 0
268  DO 10 i = 1, 4
269  iseed( i ) = iseedy( i )
270  10 CONTINUE
271 *
272 * Test the error exits
273 *
274  CALL xlaenv( 1, 1 )
275  IF( tsterr )
276  $ CALL serrge( path, nout )
277  infot = 0
278  CALL xlaenv( 2, 2 )
279 *
280 * Do for each value of M in MVAL
281 *
282  DO 120 im = 1, nm
283  m = mval( im )
284  lda = max( 1, m )
285 *
286 * Do for each value of N in NVAL
287 *
288  DO 110 in = 1, nn
289  n = nval( in )
290  xtype = 'N'
291  nimat = ntypes
292  IF( m.LE.0 .OR. n.LE.0 )
293  $ nimat = 1
294 *
295  DO 100 imat = 1, nimat
296 *
297 * Do the tests only if DOTYPE( IMAT ) is true.
298 *
299  IF( .NOT.dotype( imat ) )
300  $ GO TO 100
301 *
302 * Skip types 5, 6, or 7 if the matrix size is too small.
303 *
304  zerot = imat.GE.5 .AND. imat.LE.7
305  IF( zerot .AND. n.LT.imat-4 )
306  $ GO TO 100
307 *
308 * Set up parameters with SLATB4 and generate a test matrix
309 * with SLATMS.
310 *
311  CALL slatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
312  $ cndnum, dist )
313 *
314  srnamt = 'SLATMS'
315  CALL slatms( m, n, dist, iseed, TYPE, rwork, mode,
316  $ cndnum, anorm, kl, ku, 'No packing', a, lda,
317  $ work, info )
318 *
319 * Check error code from SLATMS.
320 *
321  IF( info.NE.0 ) THEN
322  CALL alaerh( path, 'SLATMS', info, 0, ' ', m, n, -1,
323  $ -1, -1, imat, nfail, nerrs, nout )
324  GO TO 100
325  END IF
326 *
327 * For types 5-7, zero one or more columns of the matrix to
328 * test that INFO is returned correctly.
329 *
330  IF( zerot ) THEN
331  IF( imat.EQ.5 ) THEN
332  izero = 1
333  ELSE IF( imat.EQ.6 ) THEN
334  izero = min( m, n )
335  ELSE
336  izero = min( m, n ) / 2 + 1
337  END IF
338  ioff = ( izero-1 )*lda
339  IF( imat.LT.7 ) THEN
340  DO 20 i = 1, m
341  a( ioff+i ) = zero
342  20 CONTINUE
343  ELSE
344  CALL slaset( 'Full', m, n-izero+1, zero, zero,
345  $ a( ioff+1 ), lda )
346  END IF
347  ELSE
348  izero = 0
349  END IF
350 *
351 * These lines, if used in place of the calls in the DO 60
352 * loop, cause the code to bomb on a Sun SPARCstation.
353 *
354 * ANORMO = SLANGE( 'O', M, N, A, LDA, RWORK )
355 * ANORMI = SLANGE( 'I', M, N, A, LDA, RWORK )
356 *
357 * Do for each blocksize in NBVAL
358 *
359  DO 90 inb = 1, nnb
360  nb = nbval( inb )
361  CALL xlaenv( 1, nb )
362 *
363 * Compute the LU factorization of the matrix.
364 *
365  CALL slacpy( 'Full', m, n, a, lda, afac, lda )
366  srnamt = 'SGETRF'
367  CALL sgetrf( m, n, afac, lda, iwork, info )
368 *
369 * Check error code from SGETRF.
370 *
371  IF( info.NE.izero )
372  $ CALL alaerh( path, 'SGETRF', info, izero, ' ', m,
373  $ n, -1, -1, nb, imat, nfail, nerrs,
374  $ nout )
375  trfcon = .false.
376 *
377 *+ TEST 1
378 * Reconstruct matrix from factors and compute residual.
379 *
380  CALL slacpy( 'Full', m, n, afac, lda, ainv, lda )
381  CALL sget01( m, n, a, lda, ainv, lda, iwork, rwork,
382  $ result( 1 ) )
383  nt = 1
384 *
385 *+ TEST 2
386 * Form the inverse if the factorization was successful
387 * and compute the residual.
388 *
389  IF( m.EQ.n .AND. info.EQ.0 ) THEN
390  CALL slacpy( 'Full', n, n, afac, lda, ainv, lda )
391  srnamt = 'SGETRI'
392  nrhs = nsval( 1 )
393  lwork = nmax*max( 3, nrhs )
394  CALL sgetri( n, ainv, lda, iwork, work, lwork,
395  $ info )
396 *
397 * Check error code from SGETRI.
398 *
399  IF( info.NE.0 )
400  $ CALL alaerh( path, 'SGETRI', info, 0, ' ', n, n,
401  $ -1, -1, nb, imat, nfail, nerrs,
402  $ nout )
403 *
404 * Compute the residual for the matrix times its
405 * inverse. Also compute the 1-norm condition number
406 * of A.
407 *
408  CALL sget03( n, a, lda, ainv, lda, work, lda,
409  $ rwork, rcondo, result( 2 ) )
410  anormo = slange( 'O', m, n, a, lda, rwork )
411 *
412 * Compute the infinity-norm condition number of A.
413 *
414  anormi = slange( 'I', m, n, a, lda, rwork )
415  ainvnm = slange( 'I', n, n, ainv, lda, rwork )
416  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
417  rcondi = one
418  ELSE
419  rcondi = ( one / anormi ) / ainvnm
420  END IF
421  nt = 2
422  ELSE
423 *
424 * Do only the condition estimate if INFO > 0.
425 *
426  trfcon = .true.
427  anormo = slange( 'O', m, n, a, lda, rwork )
428  anormi = slange( 'I', m, n, a, lda, rwork )
429  rcondo = zero
430  rcondi = zero
431  END IF
432 *
433 * Print information about the tests so far that did not
434 * pass the threshold.
435 *
436  DO 30 k = 1, nt
437  IF( result( k ).GE.thresh ) THEN
438  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
439  $ CALL alahd( nout, path )
440  WRITE( nout, fmt = 9999 )m, n, nb, imat, k,
441  $ result( k )
442  nfail = nfail + 1
443  END IF
444  30 CONTINUE
445  nrun = nrun + nt
446 *
447 * Skip the remaining tests if this is not the first
448 * block size or if M .ne. N. Skip the solve tests if
449 * the matrix is singular.
450 *
451  IF( inb.GT.1 .OR. m.NE.n )
452  $ GO TO 90
453  IF( trfcon )
454  $ GO TO 70
455 *
456  DO 60 irhs = 1, nns
457  nrhs = nsval( irhs )
458  xtype = 'N'
459 *
460  DO 50 itran = 1, ntran
461  trans = transs( itran )
462  IF( itran.EQ.1 ) THEN
463  rcondc = rcondo
464  ELSE
465  rcondc = rcondi
466  END IF
467 *
468 *+ TEST 3
469 * Solve and compute residual for A * X = B.
470 *
471  srnamt = 'SLARHS'
472  CALL slarhs( path, xtype, ' ', trans, n, n, kl,
473  $ ku, nrhs, a, lda, xact, lda, b,
474  $ lda, iseed, info )
475  xtype = 'C'
476 *
477  CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
478  srnamt = 'SGETRS'
479  CALL sgetrs( trans, n, nrhs, afac, lda, iwork,
480  $ x, lda, info )
481 *
482 * Check error code from SGETRS.
483 *
484  IF( info.NE.0 )
485  $ CALL alaerh( path, 'SGETRS', info, 0, trans,
486  $ n, n, -1, -1, nrhs, imat, nfail,
487  $ nerrs, nout )
488 *
489  CALL slacpy( 'Full', n, nrhs, b, lda, work,
490  $ lda )
491  CALL sget02( trans, n, n, nrhs, a, lda, x, lda,
492  $ work, lda, rwork, result( 3 ) )
493 *
494 *+ TEST 4
495 * Check solution from generated exact solution.
496 *
497  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
498  $ result( 4 ) )
499 *
500 *+ TESTS 5, 6, and 7
501 * Use iterative refinement to improve the
502 * solution.
503 *
504  srnamt = 'SGERFS'
505  CALL sgerfs( trans, n, nrhs, a, lda, afac, lda,
506  $ iwork, b, lda, x, lda, rwork,
507  $ rwork( nrhs+1 ), work,
508  $ iwork( n+1 ), info )
509 *
510 * Check error code from SGERFS.
511 *
512  IF( info.NE.0 )
513  $ CALL alaerh( path, 'SGERFS', info, 0, trans,
514  $ n, n, -1, -1, nrhs, imat, nfail,
515  $ nerrs, nout )
516 *
517  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
518  $ result( 5 ) )
519  CALL sget07( trans, n, nrhs, a, lda, b, lda, x,
520  $ lda, xact, lda, rwork, .true.,
521  $ rwork( nrhs+1 ), result( 6 ) )
522 *
523 * Print information about the tests that did not
524 * pass the threshold.
525 *
526  DO 40 k = 3, 7
527  IF( result( k ).GE.thresh ) THEN
528  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
529  $ CALL alahd( nout, path )
530  WRITE( nout, fmt = 9998 )trans, n, nrhs,
531  $ imat, k, result( k )
532  nfail = nfail + 1
533  END IF
534  40 CONTINUE
535  nrun = nrun + 5
536  50 CONTINUE
537  60 CONTINUE
538 *
539 *+ TEST 8
540 * Get an estimate of RCOND = 1/CNDNUM.
541 *
542  70 CONTINUE
543  DO 80 itran = 1, 2
544  IF( itran.EQ.1 ) THEN
545  anorm = anormo
546  rcondc = rcondo
547  norm = 'O'
548  ELSE
549  anorm = anormi
550  rcondc = rcondi
551  norm = 'I'
552  END IF
553  srnamt = 'SGECON'
554  CALL sgecon( norm, n, afac, lda, anorm, rcond,
555  $ work, iwork( n+1 ), info )
556 *
557 * Check error code from SGECON.
558 *
559  IF( info.NE.0 )
560  $ CALL alaerh( path, 'SGECON', info, 0, norm, n,
561  $ n, -1, -1, -1, imat, nfail, nerrs,
562  $ nout )
563 *
564 * This line is needed on a Sun SPARCstation.
565 *
566  dummy = rcond
567 *
568  result( 8 ) = sget06( rcond, rcondc )
569 *
570 * Print information about the tests that did not pass
571 * the threshold.
572 *
573  IF( result( 8 ).GE.thresh ) THEN
574  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
575  $ CALL alahd( nout, path )
576  WRITE( nout, fmt = 9997 )norm, n, imat, 8,
577  $ result( 8 )
578  nfail = nfail + 1
579  END IF
580  nrun = nrun + 1
581  80 CONTINUE
582  90 CONTINUE
583  100 CONTINUE
584  110 CONTINUE
585  120 CONTINUE
586 *
587 * Print a summary of the results.
588 *
589  CALL alasum( path, nout, nfail, nrun, nerrs )
590 *
591  9999 FORMAT( ' M = ', i5, ', N =', i5, ', NB =', i4, ', type ', i2,
592  $ ', test(', i2, ') =', g12.5 )
593  9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
594  $ i2, ', test(', i2, ') =', g12.5 )
595  9997 FORMAT( ' NORM =''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
596  $ ', test(', i2, ') =', g12.5 )
597  RETURN
598 *
599 * End of SCHKGE
600 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine sgetrs(TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
SGETRS
Definition: sgetrs.f:123
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine sget07(TRANS, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, CHKFERR, BERR, RESLTS)
SGET07
Definition: sget07.f:167
subroutine sget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
SGET02
Definition: sget02.f:135
subroutine slatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
SLATB4
Definition: slatb4.f:122
subroutine slarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
SLARHS
Definition: slarhs.f:206
subroutine sgecon(NORM, N, A, LDA, ANORM, RCOND, WORK, IWORK, INFO)
SGECON
Definition: sgecon.f:126
subroutine serrge(PATH, NUNIT)
SERRGE
Definition: serrge.f:57
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
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
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
subroutine sgetrf(M, N, A, LDA, IPIV, INFO)
SGETRF
Definition: sgetrf.f:110
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 sget03(N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
SGET03
Definition: sget03.f:111
subroutine sget01(M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK, RESID)
SGET01
Definition: sget01.f:109
subroutine sgerfs(TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
SGERFS
Definition: sgerfs.f:187
subroutine sgetri(N, A, LDA, IPIV, WORK, LWORK, INFO)
SGETRI
Definition: sgetri.f:116
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: