1*> \brief \b CCHKGE
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 CCHKGE( 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               RWORK( * )
25*       COMPLEX            A( * ), AFAC( * ), AINV( * ), B( * ),
26*      $                   WORK( * ), X( * ), XACT( * )
27*       ..
28*
29*
30*> \par Purpose:
31*  =============
32*>
33*> \verbatim
34*>
35*> CCHKGE tests CGETRF, -TRI, -TRS, -RFS, and -CON.
36*> \endverbatim
37*
38*  Arguments:
39*  ==========
40*
41*> \param[in] DOTYPE
42*> \verbatim
43*>          DOTYPE is LOGICAL array, dimension (NTYPES)
44*>          The matrix types to be used for testing.  Matrices of type j
45*>          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
46*>          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
47*> \endverbatim
48*>
49*> \param[in] NM
50*> \verbatim
51*>          NM is INTEGER
52*>          The number of values of M contained in the vector MVAL.
53*> \endverbatim
54*>
55*> \param[in] MVAL
56*> \verbatim
57*>          MVAL is INTEGER array, dimension (NM)
58*>          The values of the matrix row dimension M.
59*> \endverbatim
60*>
61*> \param[in] NN
62*> \verbatim
63*>          NN is INTEGER
64*>          The number of values of N contained in the vector NVAL.
65*> \endverbatim
66*>
67*> \param[in] NVAL
68*> \verbatim
69*>          NVAL is INTEGER array, dimension (NN)
70*>          The values of the matrix column dimension N.
71*> \endverbatim
72*>
73*> \param[in] NNB
74*> \verbatim
75*>          NNB is INTEGER
76*>          The number of values of NB contained in the vector NBVAL.
77*> \endverbatim
78*>
79*> \param[in] NBVAL
80*> \verbatim
81*>          NBVAL is INTEGER array, dimension (NBVAL)
82*>          The values of the blocksize NB.
83*> \endverbatim
84*>
85*> \param[in] NNS
86*> \verbatim
87*>          NNS is INTEGER
88*>          The number of values of NRHS contained in the vector NSVAL.
89*> \endverbatim
90*>
91*> \param[in] NSVAL
92*> \verbatim
93*>          NSVAL is INTEGER array, dimension (NNS)
94*>          The values of the number of right hand sides NRHS.
95*> \endverbatim
96*>
97*> \param[in] THRESH
98*> \verbatim
99*>          THRESH is REAL
100*>          The threshold value for the test ratios.  A result is
101*>          included in the output file if RESULT >= THRESH.  To have
102*>          every test ratio printed, use THRESH = 0.
103*> \endverbatim
104*>
105*> \param[in] TSTERR
106*> \verbatim
107*>          TSTERR is LOGICAL
108*>          Flag that indicates whether error exits are to be tested.
109*> \endverbatim
110*>
111*> \param[in] NMAX
112*> \verbatim
113*>          NMAX is INTEGER
114*>          The maximum value permitted for M or N, used in dimensioning
115*>          the work arrays.
116*> \endverbatim
117*>
118*> \param[out] A
119*> \verbatim
120*>          A is COMPLEX array, dimension (NMAX*NMAX)
121*> \endverbatim
122*>
123*> \param[out] AFAC
124*> \verbatim
125*>          AFAC is COMPLEX array, dimension (NMAX*NMAX)
126*> \endverbatim
127*>
128*> \param[out] AINV
129*> \verbatim
130*>          AINV is COMPLEX array, dimension (NMAX*NMAX)
131*> \endverbatim
132*>
133*> \param[out] B
134*> \verbatim
135*>          B is COMPLEX array, dimension (NMAX*NSMAX)
136*>          where NSMAX is the largest entry in NSVAL.
137*> \endverbatim
138*>
139*> \param[out] X
140*> \verbatim
141*>          X is COMPLEX array, dimension (NMAX*NSMAX)
142*> \endverbatim
143*>
144*> \param[out] XACT
145*> \verbatim
146*>          XACT is COMPLEX array, dimension (NMAX*NSMAX)
147*> \endverbatim
148*>
149*> \param[out] WORK
150*> \verbatim
151*>          WORK is COMPLEX array, dimension
152*>                      (NMAX*max(3,NSMAX))
153*> \endverbatim
154*>
155*> \param[out] RWORK
156*> \verbatim
157*>          RWORK is REAL array, dimension
158*>                      (max(2*NMAX,2*NSMAX+NWORK))
159*> \endverbatim
160*>
161*> \param[out] IWORK
162*> \verbatim
163*>          IWORK is INTEGER array, dimension (NMAX)
164*> \endverbatim
165*>
166*> \param[in] NOUT
167*> \verbatim
168*>          NOUT is INTEGER
169*>          The unit number for output.
170*> \endverbatim
171*
172*  Authors:
173*  ========
174*
175*> \author Univ. of Tennessee
176*> \author Univ. of California Berkeley
177*> \author Univ. of Colorado Denver
178*> \author NAG Ltd.
179*
180*> \date November 2011
181*
182*> \ingroup complex_lin
183*
184*  =====================================================================
185      SUBROUTINE CCHKGE( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
186     $                   NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B,
187     $                   X, XACT, WORK, RWORK, IWORK, NOUT )
188*
189*  -- LAPACK test routine (version 3.4.0) --
190*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
191*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
192*     November 2011
193*
194*     .. Scalar Arguments ..
195      LOGICAL            TSTERR
196      INTEGER            NM, NMAX, NN, NNB, NNS, NOUT
197      REAL               THRESH
198*     ..
199*     .. Array Arguments ..
200      LOGICAL            DOTYPE( * )
201      INTEGER            IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
202     $                   NVAL( * )
203      REAL               RWORK( * )
204      COMPLEX            A( * ), AFAC( * ), AINV( * ), B( * ),
205     $                   WORK( * ), X( * ), XACT( * )
206*     ..
207*
208*  =====================================================================
209*
210*     .. Parameters ..
211      REAL               ONE, ZERO
212      PARAMETER          ( ONE = 1.0E+0, ZERO = 0.0E+0 )
213      INTEGER            NTYPES
214      PARAMETER          ( NTYPES = 11 )
215      INTEGER            NTESTS
216      PARAMETER          ( NTESTS = 8 )
217      INTEGER            NTRAN
218      PARAMETER          ( NTRAN = 3 )
219*     ..
220*     .. Local Scalars ..
221      LOGICAL            TRFCON, ZEROT
222      CHARACTER          DIST, NORM, TRANS, TYPE, XTYPE
223      CHARACTER*3        PATH
224      INTEGER            I, IM, IMAT, IN, INB, INFO, IOFF, IRHS, ITRAN,
225     $                   IZERO, K, KL, KU, LDA, LWORK, M, MODE, N, NB,
226     $                   NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
227      REAL               AINVNM, ANORM, ANORMI, ANORMO, CNDNUM, DUMMY,
228     $                   RCOND, RCONDC, RCONDI, RCONDO
229*     ..
230*     .. Local Arrays ..
231      CHARACTER          TRANSS( NTRAN )
232      INTEGER            ISEED( 4 ), ISEEDY( 4 )
233      REAL               RESULT( NTESTS )
234*     ..
235*     .. External Functions ..
236      REAL               CLANGE, SGET06
237      EXTERNAL           CLANGE, SGET06
238*     ..
239*     .. External Subroutines ..
240      EXTERNAL           ALAERH, ALAHD, ALASUM, CERRGE, CGECON, CGERFS,
241     $                   CGET01, CGET02, CGET03, CGET04, CGET07, CGETRF,
242     $                   CGETRI, CGETRS, CLACPY, CLARHS, CLASET, CLATB4,
243     $                   CLATMS, XLAENV
244*     ..
245*     .. Intrinsic Functions ..
246      INTRINSIC          CMPLX, MAX, MIN
247*     ..
248*     .. Scalars in Common ..
249      LOGICAL            LERR, OK
250      CHARACTER*32       SRNAMT
251      INTEGER            INFOT, NUNIT
252*     ..
253*     .. Common blocks ..
254      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
255      COMMON             / SRNAMC / SRNAMT
256*     ..
257*     .. Data statements ..
258      DATA               ISEEDY / 1988, 1989, 1990, 1991 / ,
259     $                   TRANSS / 'N', 'T', 'C' /
260*     ..
261*     .. Executable Statements ..
262*
263*     Initialize constants and the random number seed.
264*
265      PATH( 1: 1 ) = 'Complex precision'
266      PATH( 2: 3 ) = 'GE'
267      NRUN = 0
268      NFAIL = 0
269      NERRS = 0
270      DO 10 I = 1, 4
271         ISEED( I ) = ISEEDY( I )
272   10 CONTINUE
273*
274*     Test the error exits
275*
276      CALL XLAENV( 1, 1 )
277      IF( TSTERR )
278     $   CALL CERRGE( PATH, NOUT )
279      INFOT = 0
280      CALL XLAENV( 2, 2 )
281*
282*     Do for each value of M in MVAL
283*
284      DO 120 IM = 1, NM
285         M = MVAL( IM )
286         LDA = MAX( 1, M )
287*
288*        Do for each value of N in NVAL
289*
290         DO 110 IN = 1, NN
291            N = NVAL( IN )
292            XTYPE = 'N'
293            NIMAT = NTYPES
294            IF( M.LE.0 .OR. N.LE.0 )
295     $         NIMAT = 1
296*
297            DO 100 IMAT = 1, NIMAT
298*
299*              Do the tests only if DOTYPE( IMAT ) is true.
300*
301               IF( .NOT.DOTYPE( IMAT ) )
302     $            GO TO 100
303*
304*              Skip types 5, 6, or 7 if the matrix size is too small.
305*
306               ZEROT = IMAT.GE.5 .AND. IMAT.LE.7
307               IF( ZEROT .AND. N.LT.IMAT-4 )
308     $            GO TO 100
309*
310*              Set up parameters with CLATB4 and generate a test matrix
311*              with CLATMS.
312*
313               CALL CLATB4( PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE,
314     $                      CNDNUM, DIST )
315*
316               SRNAMT = 'CLATMS'
317               CALL CLATMS( M, N, DIST, ISEED, TYPE, RWORK, MODE,
318     $                      CNDNUM, ANORM, KL, KU, 'No packing', A, LDA,
319     $                      WORK, INFO )
320*
321*              Check error code from CLATMS.
322*
323               IF( INFO.NE.0 ) THEN
324                  CALL ALAERH( PATH, 'CLATMS', INFO, 0, ' ', M, N, -1,
325     $                         -1, -1, IMAT, NFAIL, NERRS, NOUT )
326                  GO TO 100
327               END IF
328*
329*              For types 5-7, zero one or more columns of the matrix to
330*              test that INFO is returned correctly.
331*
332               IF( ZEROT ) THEN
333                  IF( IMAT.EQ.5 ) THEN
334                     IZERO = 1
335                  ELSE IF( IMAT.EQ.6 ) THEN
336                     IZERO = MIN( M, N )
337                  ELSE
338                     IZERO = MIN( M, N ) / 2 + 1
339                  END IF
340                  IOFF = ( IZERO-1 )*LDA
341                  IF( IMAT.LT.7 ) THEN
342                     DO 20 I = 1, M
343                        A( IOFF+I ) = ZERO
344   20                CONTINUE
345                  ELSE
346                     CALL CLASET( 'Full', M, N-IZERO+1, CMPLX( ZERO ),
347     $                            CMPLX( ZERO ), A( IOFF+1 ), LDA )
348                  END IF
349               ELSE
350                  IZERO = 0
351               END IF
352*
353*              These lines, if used in place of the calls in the DO 60
354*              loop, cause the code to bomb on a Sun SPARCstation.
355*
356*               ANORMO = CLANGE( 'O', M, N, A, LDA, RWORK )
357*               ANORMI = CLANGE( 'I', M, N, A, LDA, RWORK )
358*
359*              Do for each blocksize in NBVAL
360*
361               DO 90 INB = 1, NNB
362                  NB = NBVAL( INB )
363                  CALL XLAENV( 1, NB )
364*
365*                 Compute the LU factorization of the matrix.
366*
367                  CALL CLACPY( 'Full', M, N, A, LDA, AFAC, LDA )
368                  SRNAMT = 'CGETRF'
369                  CALL CGETRF( M, N, AFAC, LDA, IWORK, INFO )
370*
371*                 Check error code from CGETRF.
372*
373                  IF( INFO.NE.IZERO )
374     $               CALL ALAERH( PATH, 'CGETRF', INFO, IZERO, ' ', M,
375     $                            N, -1, -1, NB, IMAT, NFAIL, NERRS,
376     $                            NOUT )
377                  TRFCON = .FALSE.
378*
379*+    TEST 1
380*                 Reconstruct matrix from factors and compute residual.
381*
382                  CALL CLACPY( 'Full', M, N, AFAC, LDA, AINV, LDA )
383                  CALL CGET01( M, N, A, LDA, AINV, LDA, IWORK, RWORK,
384     $                         RESULT( 1 ) )
385                  NT = 1
386*
387*+    TEST 2
388*                 Form the inverse if the factorization was successful
389*                 and compute the residual.
390*
391                  IF( M.EQ.N .AND. INFO.EQ.0 ) THEN
392                     CALL CLACPY( 'Full', N, N, AFAC, LDA, AINV, LDA )
393                     SRNAMT = 'CGETRI'
394                     NRHS = NSVAL( 1 )
395                     LWORK = NMAX*MAX( 3, NRHS )
396                     CALL CGETRI( N, AINV, LDA, IWORK, WORK, LWORK,
397     $                            INFO )
398*
399*                    Check error code from CGETRI.
400*
401                     IF( INFO.NE.0 )
402     $                  CALL ALAERH( PATH, 'CGETRI', INFO, 0, ' ', N, N,
403     $                               -1, -1, NB, IMAT, NFAIL, NERRS,
404     $                               NOUT )
405*
406*                    Compute the residual for the matrix times its
407*                    inverse.  Also compute the 1-norm condition number
408*                    of A.
409*
410                     CALL CGET03( N, A, LDA, AINV, LDA, WORK, LDA,
411     $                            RWORK, RCONDO, RESULT( 2 ) )
412                     ANORMO = CLANGE( 'O', M, N, A, LDA, RWORK )
413*
414*                    Compute the infinity-norm condition number of A.
415*
416                     ANORMI = CLANGE( 'I', M, N, A, LDA, RWORK )
417                     AINVNM = CLANGE( 'I', N, N, AINV, LDA, RWORK )
418                     IF( ANORMI.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
419                        RCONDI = ONE
420                     ELSE
421                        RCONDI = ( ONE / ANORMI ) / AINVNM
422                     END IF
423                     NT = 2
424                  ELSE
425*
426*                    Do only the condition estimate if INFO > 0.
427*
428                     TRFCON = .TRUE.
429                     ANORMO = CLANGE( 'O', M, N, A, LDA, RWORK )
430                     ANORMI = CLANGE( 'I', M, N, A, LDA, RWORK )
431                     RCONDO = ZERO
432                     RCONDI = ZERO
433                  END IF
434*
435*                 Print information about the tests so far that did not
436*                 pass the threshold.
437*
438                  DO 30 K = 1, NT
439                     IF( RESULT( K ).GE.THRESH ) THEN
440                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
441     $                     CALL ALAHD( NOUT, PATH )
442                        WRITE( NOUT, FMT = 9999 )M, N, NB, IMAT, K,
443     $                     RESULT( K )
444                        NFAIL = NFAIL + 1
445                     END IF
446   30             CONTINUE
447                  NRUN = NRUN + NT
448*
449*                 Skip the remaining tests if this is not the first
450*                 block size or if M .ne. N.  Skip the solve tests if
451*                 the matrix is singular.
452*
453                  IF( INB.GT.1 .OR. M.NE.N )
454     $               GO TO 90
455                  IF( TRFCON )
456     $               GO TO 70
457*
458                  DO 60 IRHS = 1, NNS
459                     NRHS = NSVAL( IRHS )
460                     XTYPE = 'N'
461*
462                     DO 50 ITRAN = 1, NTRAN
463                        TRANS = TRANSS( ITRAN )
464                        IF( ITRAN.EQ.1 ) THEN
465                           RCONDC = RCONDO
466                        ELSE
467                           RCONDC = RCONDI
468                        END IF
469*
470*+    TEST 3
471*                       Solve and compute residual for A * X = B.
472*
473                        SRNAMT = 'CLARHS'
474                        CALL CLARHS( PATH, XTYPE, ' ', TRANS, N, N, KL,
475     $                               KU, NRHS, A, LDA, XACT, LDA, B,
476     $                               LDA, ISEED, INFO )
477                        XTYPE = 'C'
478*
479                        CALL CLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
480                        SRNAMT = 'CGETRS'
481                        CALL CGETRS( TRANS, N, NRHS, AFAC, LDA, IWORK,
482     $                               X, LDA, INFO )
483*
484*                       Check error code from CGETRS.
485*
486                        IF( INFO.NE.0 )
487     $                     CALL ALAERH( PATH, 'CGETRS', INFO, 0, TRANS,
488     $                                  N, N, -1, -1, NRHS, IMAT, NFAIL,
489     $                                  NERRS, NOUT )
490*
491                        CALL CLACPY( 'Full', N, NRHS, B, LDA, WORK,
492     $                               LDA )
493                        CALL CGET02( TRANS, N, N, NRHS, A, LDA, X, LDA,
494     $                               WORK, LDA, RWORK, RESULT( 3 ) )
495*
496*+    TEST 4
497*                       Check solution from generated exact solution.
498*
499                        CALL CGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
500     $                               RESULT( 4 ) )
501*
502*+    TESTS 5, 6, and 7
503*                       Use iterative refinement to improve the
504*                       solution.
505*
506                        SRNAMT = 'CGERFS'
507                        CALL CGERFS( TRANS, N, NRHS, A, LDA, AFAC, LDA,
508     $                               IWORK, B, LDA, X, LDA, RWORK,
509     $                               RWORK( NRHS+1 ), WORK,
510     $                               RWORK( 2*NRHS+1 ), INFO )
511*
512*                       Check error code from CGERFS.
513*
514                        IF( INFO.NE.0 )
515     $                     CALL ALAERH( PATH, 'CGERFS', INFO, 0, TRANS,
516     $                                  N, N, -1, -1, NRHS, IMAT, NFAIL,
517     $                                  NERRS, NOUT )
518*
519                        CALL CGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
520     $                               RESULT( 5 ) )
521                        CALL CGET07( TRANS, N, NRHS, A, LDA, B, LDA, X,
522     $                               LDA, XACT, LDA, RWORK, .TRUE.,
523     $                               RWORK( NRHS+1 ), RESULT( 6 ) )
524*
525*                       Print information about the tests that did not
526*                       pass the threshold.
527*
528                        DO 40 K = 3, 7
529                           IF( RESULT( K ).GE.THRESH ) THEN
530                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
531     $                           CALL ALAHD( NOUT, PATH )
532                              WRITE( NOUT, FMT = 9998 )TRANS, N, NRHS,
533     $                           IMAT, K, RESULT( K )
534                              NFAIL = NFAIL + 1
535                           END IF
536   40                   CONTINUE
537                        NRUN = NRUN + 5
538   50                CONTINUE
539   60             CONTINUE
540*
541*+    TEST 8
542*                    Get an estimate of RCOND = 1/CNDNUM.
543*
544   70             CONTINUE
545                  DO 80 ITRAN = 1, 2
546                     IF( ITRAN.EQ.1 ) THEN
547                        ANORM = ANORMO
548                        RCONDC = RCONDO
549                        NORM = 'O'
550                     ELSE
551                        ANORM = ANORMI
552                        RCONDC = RCONDI
553                        NORM = 'I'
554                     END IF
555                     SRNAMT = 'CGECON'
556                     CALL CGECON( NORM, N, AFAC, LDA, ANORM, RCOND,
557     $                            WORK, RWORK, INFO )
558*
559*                       Check error code from CGECON.
560*
561                     IF( INFO.NE.0 )
562     $                  CALL ALAERH( PATH, 'CGECON', INFO, 0, NORM, N,
563     $                               N, -1, -1, -1, IMAT, NFAIL, NERRS,
564     $                               NOUT )
565*
566*                       This line is needed on a Sun SPARCstation.
567*
568                     DUMMY = RCOND
569*
570                     RESULT( 8 ) = SGET06( RCOND, RCONDC )
571*
572*                    Print information about the tests that did not pass
573*                    the threshold.
574*
575                     IF( RESULT( 8 ).GE.THRESH ) THEN
576                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
577     $                     CALL ALAHD( NOUT, PATH )
578                        WRITE( NOUT, FMT = 9997 )NORM, N, IMAT, 8,
579     $                     RESULT( 8 )
580                        NFAIL = NFAIL + 1
581                     END IF
582                     NRUN = NRUN + 1
583   80             CONTINUE
584   90          CONTINUE
585  100       CONTINUE
586*
587  110    CONTINUE
588  120 CONTINUE
589*
590*     Print a summary of the results.
591*
592      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
593*
594 9999 FORMAT( ' M = ', I5, ', N =', I5, ', NB =', I4, ', type ', I2,
595     $      ', test(', I2, ') =', G12.5 )
596 9998 FORMAT( ' TRANS=''', A1, ''', N =', I5, ', NRHS=', I3, ', type ',
597     $      I2, ', test(', I2, ') =', G12.5 )
598 9997 FORMAT( ' NORM =''', A1, ''', N =', I5, ',', 10X, ' type ', I2,
599     $      ', test(', I2, ') =', G12.5 )
600      RETURN
601*
602*     End of CCHKGE
603*
604      END
605