1*> \brief \b ZDRVPO
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 ZDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12*                          A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
13*                          RWORK, NOUT )
14*
15*       .. Scalar Arguments ..
16*       LOGICAL            TSTERR
17*       INTEGER            NMAX, NN, NOUT, NRHS
18*       DOUBLE PRECISION   THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            NVAL( * )
23*       DOUBLE PRECISION   RWORK( * ), S( * )
24*       COMPLEX*16         A( * ), AFAC( * ), ASAV( * ), B( * ),
25*      $                   BSAV( * ), WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> ZDRVPO tests the driver routines ZPOSV and -SVX.
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] NN
49*> \verbatim
50*>          NN is INTEGER
51*>          The number of values of N contained in the vector NVAL.
52*> \endverbatim
53*>
54*> \param[in] NVAL
55*> \verbatim
56*>          NVAL is INTEGER array, dimension (NN)
57*>          The values of the matrix dimension N.
58*> \endverbatim
59*>
60*> \param[in] NRHS
61*> \verbatim
62*>          NRHS is INTEGER
63*>          The number of right hand side vectors to be generated for
64*>          each linear system.
65*> \endverbatim
66*>
67*> \param[in] THRESH
68*> \verbatim
69*>          THRESH is DOUBLE PRECISION
70*>          The threshold value for the test ratios.  A result is
71*>          included in the output file if RESULT >= THRESH.  To have
72*>          every test ratio printed, use THRESH = 0.
73*> \endverbatim
74*>
75*> \param[in] TSTERR
76*> \verbatim
77*>          TSTERR is LOGICAL
78*>          Flag that indicates whether error exits are to be tested.
79*> \endverbatim
80*>
81*> \param[in] NMAX
82*> \verbatim
83*>          NMAX is INTEGER
84*>          The maximum value permitted for N, used in dimensioning the
85*>          work arrays.
86*> \endverbatim
87*>
88*> \param[out] A
89*> \verbatim
90*>          A is COMPLEX*16 array, dimension (NMAX*NMAX)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*>          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
96*> \endverbatim
97*>
98*> \param[out] ASAV
99*> \verbatim
100*>          ASAV is COMPLEX*16 array, dimension (NMAX*NMAX)
101*> \endverbatim
102*>
103*> \param[out] B
104*> \verbatim
105*>          B is COMPLEX*16 array, dimension (NMAX*NRHS)
106*> \endverbatim
107*>
108*> \param[out] BSAV
109*> \verbatim
110*>          BSAV is COMPLEX*16 array, dimension (NMAX*NRHS)
111*> \endverbatim
112*>
113*> \param[out] X
114*> \verbatim
115*>          X is COMPLEX*16 array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] XACT
119*> \verbatim
120*>          XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
121*> \endverbatim
122*>
123*> \param[out] S
124*> \verbatim
125*>          S is DOUBLE PRECISION array, dimension (NMAX)
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*>          WORK is COMPLEX*16 array, dimension
131*>                      (NMAX*max(3,NRHS))
132*> \endverbatim
133*>
134*> \param[out] RWORK
135*> \verbatim
136*>          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
137*> \endverbatim
138*>
139*> \param[in] NOUT
140*> \verbatim
141*>          NOUT is INTEGER
142*>          The unit number for output.
143*> \endverbatim
144*
145*  Authors:
146*  ========
147*
148*> \author Univ. of Tennessee
149*> \author Univ. of California Berkeley
150*> \author Univ. of Colorado Denver
151*> \author NAG Ltd.
152*
153*> \date December 2016
154*
155*> \ingroup complex16_lin
156*
157*  =====================================================================
158      SUBROUTINE ZDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
159     $                   A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
160     $                   RWORK, NOUT )
161*
162*  -- LAPACK test routine (version 3.7.0) --
163*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
164*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
165*     December 2016
166*
167*     .. Scalar Arguments ..
168      LOGICAL            TSTERR
169      INTEGER            NMAX, NN, NOUT, NRHS
170      DOUBLE PRECISION   THRESH
171*     ..
172*     .. Array Arguments ..
173      LOGICAL            DOTYPE( * )
174      INTEGER            NVAL( * )
175      DOUBLE PRECISION   RWORK( * ), S( * )
176      COMPLEX*16         A( * ), AFAC( * ), ASAV( * ), B( * ),
177     $                   BSAV( * ), WORK( * ), X( * ), XACT( * )
178*     ..
179*
180*  =====================================================================
181*
182*     .. Parameters ..
183      DOUBLE PRECISION   ONE, ZERO
184      PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
185      INTEGER            NTYPES
186      PARAMETER          ( NTYPES = 9 )
187      INTEGER            NTESTS
188      PARAMETER          ( NTESTS = 6 )
189*     ..
190*     .. Local Scalars ..
191      LOGICAL            EQUIL, NOFACT, PREFAC, ZEROT
192      CHARACTER          DIST, EQUED, FACT, TYPE, UPLO, XTYPE
193      CHARACTER*3        PATH
194      INTEGER            I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
195     $                   IZERO, K, K1, KL, KU, LDA, MODE, N, NB, NBMIN,
196     $                   NERRS, NFACT, NFAIL, NIMAT, NRUN, NT
197      DOUBLE PRECISION   AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
198     $                   ROLDC, SCOND
199*     ..
200*     .. Local Arrays ..
201      CHARACTER          EQUEDS( 2 ), FACTS( 3 ), UPLOS( 2 )
202      INTEGER            ISEED( 4 ), ISEEDY( 4 )
203      DOUBLE PRECISION   RESULT( NTESTS )
204*     ..
205*     .. External Functions ..
206      LOGICAL            LSAME
207      DOUBLE PRECISION   DGET06, ZLANHE
208      EXTERNAL           LSAME, DGET06, ZLANHE
209*     ..
210*     .. External Subroutines ..
211      EXTERNAL           ALADHD, ALAERH, ALASVM, XLAENV, ZERRVX, ZGET04,
212     $                   ZLACPY, ZLAIPD, ZLAQHE, ZLARHS, ZLASET, ZLATB4,
213     $                   ZLATMS, ZPOEQU, ZPOSV, ZPOSVX, ZPOT01, ZPOT02,
214     $                   ZPOT05, ZPOTRF, ZPOTRI
215*     ..
216*     .. Scalars in Common ..
217      LOGICAL            LERR, OK
218      CHARACTER*32       SRNAMT
219      INTEGER            INFOT, NUNIT
220*     ..
221*     .. Common blocks ..
222      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
223      COMMON             / SRNAMC / SRNAMT
224*     ..
225*     .. Intrinsic Functions ..
226      INTRINSIC          DCMPLX, MAX
227*     ..
228*     .. Data statements ..
229      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
230      DATA               UPLOS / 'U', 'L' /
231      DATA               FACTS / 'F', 'N', 'E' /
232      DATA               EQUEDS / 'N', 'Y' /
233*     ..
234*     .. Executable Statements ..
235*
236*     Initialize constants and the random number seed.
237*
238      PATH( 1: 1 ) = 'Zomplex precision'
239      PATH( 2: 3 ) = 'PO'
240      NRUN = 0
241      NFAIL = 0
242      NERRS = 0
243      DO 10 I = 1, 4
244         ISEED( I ) = ISEEDY( I )
245   10 CONTINUE
246*
247*     Test the error exits
248*
249      IF( TSTERR )
250     $   CALL ZERRVX( PATH, NOUT )
251      INFOT = 0
252*
253*     Set the block size and minimum block size for testing.
254*
255      NB = 1
256      NBMIN = 2
257      CALL XLAENV( 1, NB )
258      CALL XLAENV( 2, NBMIN )
259*
260*     Do for each value of N in NVAL
261*
262      DO 130 IN = 1, NN
263         N = NVAL( IN )
264         LDA = MAX( N, 1 )
265         XTYPE = 'N'
266         NIMAT = NTYPES
267         IF( N.LE.0 )
268     $      NIMAT = 1
269*
270         DO 120 IMAT = 1, NIMAT
271*
272*           Do the tests only if DOTYPE( IMAT ) is true.
273*
274            IF( .NOT.DOTYPE( IMAT ) )
275     $         GO TO 120
276*
277*           Skip types 3, 4, or 5 if the matrix size is too small.
278*
279            ZEROT = IMAT.GE.3 .AND. IMAT.LE.5
280            IF( ZEROT .AND. N.LT.IMAT-2 )
281     $         GO TO 120
282*
283*           Do first for UPLO = 'U', then for UPLO = 'L'
284*
285            DO 110 IUPLO = 1, 2
286               UPLO = UPLOS( IUPLO )
287*
288*              Set up parameters with ZLATB4 and generate a test matrix
289*              with ZLATMS.
290*
291               CALL ZLATB4( PATH, IMAT, N, N, TYPE, KL, KU, ANORM, MODE,
292     $                      CNDNUM, DIST )
293*
294               SRNAMT = 'ZLATMS'
295               CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
296     $                      CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
297     $                      INFO )
298*
299*              Check error code from ZLATMS.
300*
301               IF( INFO.NE.0 ) THEN
302                  CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N, -1,
303     $                         -1, -1, IMAT, NFAIL, NERRS, NOUT )
304                  GO TO 110
305               END IF
306*
307*              For types 3-5, zero one row and column of the matrix to
308*              test that INFO is returned correctly.
309*
310               IF( ZEROT ) THEN
311                  IF( IMAT.EQ.3 ) THEN
312                     IZERO = 1
313                  ELSE IF( IMAT.EQ.4 ) THEN
314                     IZERO = N
315                  ELSE
316                     IZERO = N / 2 + 1
317                  END IF
318                  IOFF = ( IZERO-1 )*LDA
319*
320*                 Set row and column IZERO of A to 0.
321*
322                  IF( IUPLO.EQ.1 ) THEN
323                     DO 20 I = 1, IZERO - 1
324                        A( IOFF+I ) = ZERO
325   20                CONTINUE
326                     IOFF = IOFF + IZERO
327                     DO 30 I = IZERO, N
328                        A( IOFF ) = ZERO
329                        IOFF = IOFF + LDA
330   30                CONTINUE
331                  ELSE
332                     IOFF = IZERO
333                     DO 40 I = 1, IZERO - 1
334                        A( IOFF ) = ZERO
335                        IOFF = IOFF + LDA
336   40                CONTINUE
337                     IOFF = IOFF - IZERO
338                     DO 50 I = IZERO, N
339                        A( IOFF+I ) = ZERO
340   50                CONTINUE
341                  END IF
342               ELSE
343                  IZERO = 0
344               END IF
345*
346*              Set the imaginary part of the diagonals.
347*
348               CALL ZLAIPD( N, A, LDA+1, 0 )
349*
350*              Save a copy of the matrix A in ASAV.
351*
352               CALL ZLACPY( UPLO, N, N, A, LDA, ASAV, LDA )
353*
354               DO 100 IEQUED = 1, 2
355                  EQUED = EQUEDS( IEQUED )
356                  IF( IEQUED.EQ.1 ) THEN
357                     NFACT = 3
358                  ELSE
359                     NFACT = 1
360                  END IF
361*
362                  DO 90 IFACT = 1, NFACT
363                     FACT = FACTS( IFACT )
364                     PREFAC = LSAME( FACT, 'F' )
365                     NOFACT = LSAME( FACT, 'N' )
366                     EQUIL = LSAME( FACT, 'E' )
367*
368                     IF( ZEROT ) THEN
369                        IF( PREFAC )
370     $                     GO TO 90
371                        RCONDC = ZERO
372*
373                     ELSE IF( .NOT.LSAME( FACT, 'N' ) ) THEN
374*
375*                       Compute the condition number for comparison with
376*                       the value returned by ZPOSVX (FACT = 'N' reuses
377*                       the condition number from the previous iteration
378*                       with FACT = 'F').
379*
380                        CALL ZLACPY( UPLO, N, N, ASAV, LDA, AFAC, LDA )
381                        IF( EQUIL .OR. IEQUED.GT.1 ) THEN
382*
383*                          Compute row and column scale factors to
384*                          equilibrate the matrix A.
385*
386                           CALL ZPOEQU( N, AFAC, LDA, S, SCOND, AMAX,
387     $                                  INFO )
388                           IF( INFO.EQ.0 .AND. N.GT.0 ) THEN
389                              IF( IEQUED.GT.1 )
390     $                           SCOND = ZERO
391*
392*                             Equilibrate the matrix.
393*
394                              CALL ZLAQHE( UPLO, N, AFAC, LDA, S, SCOND,
395     $                                     AMAX, EQUED )
396                           END IF
397                        END IF
398*
399*                       Save the condition number of the
400*                       non-equilibrated system for use in ZGET04.
401*
402                        IF( EQUIL )
403     $                     ROLDC = RCONDC
404*
405*                       Compute the 1-norm of A.
406*
407                        ANORM = ZLANHE( '1', UPLO, N, AFAC, LDA, RWORK )
408*
409*                       Factor the matrix A.
410*
411                        CALL ZPOTRF( UPLO, N, AFAC, LDA, INFO )
412*
413*                       Form the inverse of A.
414*
415                        CALL ZLACPY( UPLO, N, N, AFAC, LDA, A, LDA )
416                        CALL ZPOTRI( UPLO, N, A, LDA, INFO )
417*
418*                       Compute the 1-norm condition number of A.
419*
420                        AINVNM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
421                        IF( ANORM.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
422                           RCONDC = ONE
423                        ELSE
424                           RCONDC = ( ONE / ANORM ) / AINVNM
425                        END IF
426                     END IF
427*
428*                    Restore the matrix A.
429*
430                     CALL ZLACPY( UPLO, N, N, ASAV, LDA, A, LDA )
431*
432*                    Form an exact solution and set the right hand side.
433*
434                     SRNAMT = 'ZLARHS'
435                     CALL ZLARHS( PATH, XTYPE, UPLO, ' ', N, N, KL, KU,
436     $                            NRHS, A, LDA, XACT, LDA, B, LDA,
437     $                            ISEED, INFO )
438                     XTYPE = 'C'
439                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, BSAV, LDA )
440*
441                     IF( NOFACT ) THEN
442*
443*                       --- Test ZPOSV  ---
444*
445*                       Compute the L*L' or U'*U factorization of the
446*                       matrix and solve the system.
447*
448                        CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
449                        CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
450*
451                        SRNAMT = 'ZPOSV '
452                        CALL ZPOSV( UPLO, N, NRHS, AFAC, LDA, X, LDA,
453     $                              INFO )
454*
455*                       Check error code from ZPOSV .
456*
457                        IF( INFO.NE.IZERO ) THEN
458                           CALL ALAERH( PATH, 'ZPOSV ', INFO, IZERO,
459     $                                  UPLO, N, N, -1, -1, NRHS, IMAT,
460     $                                  NFAIL, NERRS, NOUT )
461                           GO TO 70
462                        ELSE IF( INFO.NE.0 ) THEN
463                           GO TO 70
464                        END IF
465*
466*                       Reconstruct matrix from factors and compute
467*                       residual.
468*
469                        CALL ZPOT01( UPLO, N, A, LDA, AFAC, LDA, RWORK,
470     $                               RESULT( 1 ) )
471*
472*                       Compute residual of the computed solution.
473*
474                        CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK,
475     $                               LDA )
476                        CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA,
477     $                               WORK, LDA, RWORK, RESULT( 2 ) )
478*
479*                       Check solution from generated exact solution.
480*
481                        CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
482     $                               RESULT( 3 ) )
483                        NT = 3
484*
485*                       Print information about the tests that did not
486*                       pass the threshold.
487*
488                        DO 60 K = 1, NT
489                           IF( RESULT( K ).GE.THRESH ) THEN
490                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
491     $                           CALL ALADHD( NOUT, PATH )
492                              WRITE( NOUT, FMT = 9999 )'ZPOSV ', UPLO,
493     $                           N, IMAT, K, RESULT( K )
494                              NFAIL = NFAIL + 1
495                           END IF
496   60                   CONTINUE
497                        NRUN = NRUN + NT
498   70                   CONTINUE
499                     END IF
500*
501*                    --- Test ZPOSVX ---
502*
503                     IF( .NOT.PREFAC )
504     $                  CALL ZLASET( UPLO, N, N, DCMPLX( ZERO ),
505     $                               DCMPLX( ZERO ), AFAC, LDA )
506                     CALL ZLASET( 'Full', N, NRHS, DCMPLX( ZERO ),
507     $                            DCMPLX( ZERO ), X, LDA )
508                     IF( IEQUED.GT.1 .AND. N.GT.0 ) THEN
509*
510*                       Equilibrate the matrix if FACT='F' and
511*                       EQUED='Y'.
512*
513                        CALL ZLAQHE( UPLO, N, A, LDA, S, SCOND, AMAX,
514     $                               EQUED )
515                     END IF
516*
517*                    Solve the system and compute the condition number
518*                    and error bounds using ZPOSVX.
519*
520                     SRNAMT = 'ZPOSVX'
521                     CALL ZPOSVX( FACT, UPLO, N, NRHS, A, LDA, AFAC,
522     $                            LDA, EQUED, S, B, LDA, X, LDA, RCOND,
523     $                            RWORK, RWORK( NRHS+1 ), WORK,
524     $                            RWORK( 2*NRHS+1 ), INFO )
525*
526*                    Check the error code from ZPOSVX.
527*
528                     IF( INFO.NE.IZERO ) THEN
529                        CALL ALAERH( PATH, 'ZPOSVX', INFO, IZERO,
530     $                               FACT // UPLO, N, N, -1, -1, NRHS,
531     $                               IMAT, NFAIL, NERRS, NOUT )
532                        GO TO 90
533                     END IF
534*
535                     IF( INFO.EQ.0 ) THEN
536                        IF( .NOT.PREFAC ) THEN
537*
538*                          Reconstruct matrix from factors and compute
539*                          residual.
540*
541                           CALL ZPOT01( UPLO, N, A, LDA, AFAC, LDA,
542     $                                  RWORK( 2*NRHS+1 ), RESULT( 1 ) )
543                           K1 = 1
544                        ELSE
545                           K1 = 2
546                        END IF
547*
548*                       Compute residual of the computed solution.
549*
550                        CALL ZLACPY( 'Full', N, NRHS, BSAV, LDA, WORK,
551     $                               LDA )
552                        CALL ZPOT02( UPLO, N, NRHS, ASAV, LDA, X, LDA,
553     $                               WORK, LDA, RWORK( 2*NRHS+1 ),
554     $                               RESULT( 2 ) )
555*
556*                       Check solution from generated exact solution.
557*
558                        IF( NOFACT .OR. ( PREFAC .AND. LSAME( EQUED,
559     $                      'N' ) ) ) THEN
560                           CALL ZGET04( N, NRHS, X, LDA, XACT, LDA,
561     $                                  RCONDC, RESULT( 3 ) )
562                        ELSE
563                           CALL ZGET04( N, NRHS, X, LDA, XACT, LDA,
564     $                                  ROLDC, RESULT( 3 ) )
565                        END IF
566*
567*                       Check the error bounds from iterative
568*                       refinement.
569*
570                        CALL ZPOT05( UPLO, N, NRHS, ASAV, LDA, B, LDA,
571     $                               X, LDA, XACT, LDA, RWORK,
572     $                               RWORK( NRHS+1 ), RESULT( 4 ) )
573                     ELSE
574                        K1 = 6
575                     END IF
576*
577*                    Compare RCOND from ZPOSVX with the computed value
578*                    in RCONDC.
579*
580                     RESULT( 6 ) = DGET06( RCOND, RCONDC )
581*
582*                    Print information about the tests that did not pass
583*                    the threshold.
584*
585                     DO 80 K = K1, 6
586                        IF( RESULT( K ).GE.THRESH ) THEN
587                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
588     $                        CALL ALADHD( NOUT, PATH )
589                           IF( PREFAC ) THEN
590                              WRITE( NOUT, FMT = 9997 )'ZPOSVX', FACT,
591     $                           UPLO, N, EQUED, IMAT, K, RESULT( K )
592                           ELSE
593                              WRITE( NOUT, FMT = 9998 )'ZPOSVX', FACT,
594     $                           UPLO, N, IMAT, K, RESULT( K )
595                           END IF
596                           NFAIL = NFAIL + 1
597                        END IF
598   80                CONTINUE
599                     NRUN = NRUN + 7 - K1
600   90             CONTINUE
601  100          CONTINUE
602  110       CONTINUE
603  120    CONTINUE
604  130 CONTINUE
605*
606*     Print a summary of the results.
607*
608      CALL ALASVM( PATH, NOUT, NFAIL, NRUN, NERRS )
609*
610 9999 FORMAT( 1X, A, ', UPLO=''', A1, ''', N =', I5, ', type ', I1,
611     $      ', test(', I1, ')=', G12.5 )
612 9998 FORMAT( 1X, A, ', FACT=''', A1, ''', UPLO=''', A1, ''', N=', I5,
613     $      ', type ', I1, ', test(', I1, ')=', G12.5 )
614 9997 FORMAT( 1X, A, ', FACT=''', A1, ''', UPLO=''', A1, ''', N=', I5,
615     $      ', EQUED=''', A1, ''', type ', I1, ', test(', I1, ') =',
616     $      G12.5 )
617      RETURN
618*
619*     End of ZDRVPO
620*
621      END
622