1*> \brief \b SCHKGE
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*  Definition:
9*  ===========
10*
11*       SUBROUTINE SCHKGE( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
12*                          NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B,
13*                          X, XACT, WORK, RWORK, IWORK, NOUT )
14*
15*       .. Scalar Arguments ..
16*       LOGICAL            TSTERR
17*       INTEGER            NM, NMAX, NN, NNB, NNS, NOUT
18*       REAL               THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
23*      $                   NVAL( * )
24*       REAL               A( * ), AFAC( * ), AINV( * ), B( * ),
25*      $                   RWORK( * ), WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> SCHKGE tests SGETRF, -TRI, -TRS, -RFS, and -CON.
35*> \endverbatim
36*
37*  Arguments:
38*  ==========
39*
40*> \param[in] DOTYPE
41*> \verbatim
42*>          DOTYPE is LOGICAL array, dimension (NTYPES)
43*>          The matrix types to be used for testing.  Matrices of type j
44*>          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
45*>          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
46*> \endverbatim
47*>
48*> \param[in] NM
49*> \verbatim
50*>          NM is INTEGER
51*>          The number of values of M contained in the vector MVAL.
52*> \endverbatim
53*>
54*> \param[in] MVAL
55*> \verbatim
56*>          MVAL is INTEGER array, dimension (NM)
57*>          The values of the matrix row dimension M.
58*> \endverbatim
59*>
60*> \param[in] NN
61*> \verbatim
62*>          NN is INTEGER
63*>          The number of values of N contained in the vector NVAL.
64*> \endverbatim
65*>
66*> \param[in] NVAL
67*> \verbatim
68*>          NVAL is INTEGER array, dimension (NN)
69*>          The values of the matrix column dimension N.
70*> \endverbatim
71*>
72*> \param[in] NNB
73*> \verbatim
74*>          NNB is INTEGER
75*>          The number of values of NB contained in the vector NBVAL.
76*> \endverbatim
77*>
78*> \param[in] NBVAL
79*> \verbatim
80*>          NBVAL is INTEGER array, dimension (NBVAL)
81*>          The values of the blocksize NB.
82*> \endverbatim
83*>
84*> \param[in] NNS
85*> \verbatim
86*>          NNS is INTEGER
87*>          The number of values of NRHS contained in the vector NSVAL.
88*> \endverbatim
89*>
90*> \param[in] NSVAL
91*> \verbatim
92*>          NSVAL is INTEGER array, dimension (NNS)
93*>          The values of the number of right hand sides NRHS.
94*> \endverbatim
95*>
96*> \param[in] THRESH
97*> \verbatim
98*>          THRESH is REAL
99*>          The threshold value for the test ratios.  A result is
100*>          included in the output file if RESULT >= THRESH.  To have
101*>          every test ratio printed, use THRESH = 0.
102*> \endverbatim
103*>
104*> \param[in] TSTERR
105*> \verbatim
106*>          TSTERR is LOGICAL
107*>          Flag that indicates whether error exits are to be tested.
108*> \endverbatim
109*>
110*> \param[in] NMAX
111*> \verbatim
112*>          NMAX is INTEGER
113*>          The maximum value permitted for M or N, used in dimensioning
114*>          the work arrays.
115*> \endverbatim
116*>
117*> \param[out] A
118*> \verbatim
119*>          A is REAL array, dimension (NMAX*NMAX)
120*> \endverbatim
121*>
122*> \param[out] AFAC
123*> \verbatim
124*>          AFAC is REAL array, dimension (NMAX*NMAX)
125*> \endverbatim
126*>
127*> \param[out] AINV
128*> \verbatim
129*>          AINV is REAL array, dimension (NMAX*NMAX)
130*> \endverbatim
131*>
132*> \param[out] B
133*> \verbatim
134*>          B is REAL array, dimension (NMAX*NSMAX)
135*>          where NSMAX is the largest entry in NSVAL.
136*> \endverbatim
137*>
138*> \param[out] X
139*> \verbatim
140*>          X is REAL array, dimension (NMAX*NSMAX)
141*> \endverbatim
142*>
143*> \param[out] XACT
144*> \verbatim
145*>          XACT is REAL array, dimension (NMAX*NSMAX)
146*> \endverbatim
147*>
148*> \param[out] WORK
149*> \verbatim
150*>          WORK is REAL array, dimension
151*>                      (NMAX*max(3,NSMAX))
152*> \endverbatim
153*>
154*> \param[out] RWORK
155*> \verbatim
156*>          RWORK is REAL array, dimension
157*>                      (max(2*NMAX,2*NSMAX+NWORK))
158*> \endverbatim
159*>
160*> \param[out] IWORK
161*> \verbatim
162*>          IWORK is INTEGER array, dimension (2*NMAX)
163*> \endverbatim
164*>
165*> \param[in] NOUT
166*> \verbatim
167*>          NOUT is INTEGER
168*>          The unit number for output.
169*> \endverbatim
170*
171*  Authors:
172*  ========
173*
174*> \author Univ. of Tennessee
175*> \author Univ. of California Berkeley
176*> \author Univ. of Colorado Denver
177*> \author NAG Ltd.
178*
179*> \date April 2012
180*
181*> \ingroup single_lin
182*
183*  =====================================================================
184      SUBROUTINE SCHKGE( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
185     $                   NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B,
186     $                   X, XACT, WORK, RWORK, IWORK, NOUT )
187*
188*  -- LAPACK test routine (version 3.7.0) --
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,
239     $                   SGET01, SGET02, SGET03, SGET04, SGET07, SGETRF,
240     $                   SGETRI, SGETRS, SLACPY, SLARHS, SLASET, SLATB4,
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*
601      END
602