LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cdrvpp.f
Go to the documentation of this file.
1*> \brief \b CDRVPP
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 CDRVPP( 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* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER NVAL( * )
23* REAL RWORK( * ), S( * )
24* COMPLEX A( * ), AFAC( * ), ASAV( * ), B( * ),
25* $ BSAV( * ), WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> CDRVPP tests the driver routines CPPSV 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 REAL
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 array, dimension (NMAX*(NMAX+1)/2)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*> AFAC is COMPLEX array, dimension (NMAX*(NMAX+1)/2)
96*> \endverbatim
97*>
98*> \param[out] ASAV
99*> \verbatim
100*> ASAV is COMPLEX array, dimension (NMAX*(NMAX+1)/2)
101*> \endverbatim
102*>
103*> \param[out] B
104*> \verbatim
105*> B is COMPLEX array, dimension (NMAX*NRHS)
106*> \endverbatim
107*>
108*> \param[out] BSAV
109*> \verbatim
110*> BSAV is COMPLEX array, dimension (NMAX*NRHS)
111*> \endverbatim
112*>
113*> \param[out] X
114*> \verbatim
115*> X is COMPLEX array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] XACT
119*> \verbatim
120*> XACT is COMPLEX array, dimension (NMAX*NRHS)
121*> \endverbatim
122*>
123*> \param[out] S
124*> \verbatim
125*> S is REAL array, dimension (NMAX)
126*> \endverbatim
127*>
128*> \param[out] WORK
129*> \verbatim
130*> WORK is COMPLEX array, dimension
131*> (NMAX*max(3,NRHS))
132*> \endverbatim
133*>
134*> \param[out] RWORK
135*> \verbatim
136*> RWORK is REAL 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*> \ingroup complex_lin
154*
155* =====================================================================
156 SUBROUTINE cdrvpp( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
157 $ A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
158 $ RWORK, NOUT )
159*
160* -- LAPACK test routine --
161* -- LAPACK is a software package provided by Univ. of Tennessee, --
162* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163*
164* .. Scalar Arguments ..
165 LOGICAL TSTERR
166 INTEGER NMAX, NN, NOUT, NRHS
167 REAL THRESH
168* ..
169* .. Array Arguments ..
170 LOGICAL DOTYPE( * )
171 INTEGER NVAL( * )
172 REAL RWORK( * ), S( * )
173 COMPLEX A( * ), AFAC( * ), ASAV( * ), B( * ),
174 $ bsav( * ), work( * ), x( * ), xact( * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 REAL ONE, ZERO
181 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
182 INTEGER NTYPES
183 parameter( ntypes = 9 )
184 INTEGER NTESTS
185 parameter( ntests = 6 )
186* ..
187* .. Local Scalars ..
188 LOGICAL EQUIL, NOFACT, PREFAC, ZEROT
189 CHARACTER DIST, EQUED, FACT, PACKIT, TYPE, UPLO, XTYPE
190 CHARACTER*3 PATH
191 INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
192 $ izero, k, k1, kl, ku, lda, mode, n, nerrs,
193 $ nfact, nfail, nimat, npp, nrun, nt
194 REAL AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
195 $ ROLDC, SCOND
196* ..
197* .. Local Arrays ..
198 CHARACTER EQUEDS( 2 ), FACTS( 3 ), PACKS( 2 ), UPLOS( 2 )
199 INTEGER ISEED( 4 ), ISEEDY( 4 )
200 REAL RESULT( NTESTS )
201* ..
202* .. External Functions ..
203 LOGICAL LSAME
204 REAL CLANHP, SGET06
205 EXTERNAL lsame, clanhp, sget06
206* ..
207* .. External Subroutines ..
208 EXTERNAL aladhd, alaerh, alasvm, ccopy, cerrvx, cget04,
212* ..
213* .. Scalars in Common ..
214 LOGICAL LERR, OK
215 CHARACTER*32 SRNAMT
216 INTEGER INFOT, NUNIT
217* ..
218* .. Common blocks ..
219 COMMON / infoc / infot, nunit, ok, lerr
220 COMMON / srnamc / srnamt
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC cmplx, max
224* ..
225* .. Data statements ..
226 DATA iseedy / 1988, 1989, 1990, 1991 /
227 DATA uplos / 'U', 'L' / , facts / 'F', 'N', 'E' / ,
228 $ packs / 'C', 'R' / , equeds / 'N', 'Y' /
229* ..
230* .. Executable Statements ..
231*
232* Initialize constants and the random number seed.
233*
234 path( 1: 1 ) = 'Complex precision'
235 path( 2: 3 ) = 'PP'
236 nrun = 0
237 nfail = 0
238 nerrs = 0
239 DO 10 i = 1, 4
240 iseed( i ) = iseedy( i )
241 10 CONTINUE
242*
243* Test the error exits
244*
245 IF( tsterr )
246 $ CALL cerrvx( path, nout )
247 infot = 0
248*
249* Do for each value of N in NVAL
250*
251 DO 140 in = 1, nn
252 n = nval( in )
253 lda = max( n, 1 )
254 npp = n*( n+1 ) / 2
255 xtype = 'N'
256 nimat = ntypes
257 IF( n.LE.0 )
258 $ nimat = 1
259*
260 DO 130 imat = 1, nimat
261*
262* Do the tests only if DOTYPE( IMAT ) is true.
263*
264 IF( .NOT.dotype( imat ) )
265 $ GO TO 130
266*
267* Skip types 3, 4, or 5 if the matrix size is too small.
268*
269 zerot = imat.GE.3 .AND. imat.LE.5
270 IF( zerot .AND. n.LT.imat-2 )
271 $ GO TO 130
272*
273* Do first for UPLO = 'U', then for UPLO = 'L'
274*
275 DO 120 iuplo = 1, 2
276 uplo = uplos( iuplo )
277 packit = packs( iuplo )
278*
279* Set up parameters with CLATB4 and generate a test matrix
280* with CLATMS.
281*
282 CALL clatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
283 $ cndnum, dist )
284 rcondc = one / cndnum
285*
286 srnamt = 'CLATMS'
287 CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
288 $ cndnum, anorm, kl, ku, packit, a, lda, work,
289 $ info )
290*
291* Check error code from CLATMS.
292*
293 IF( info.NE.0 ) THEN
294 CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n, -1,
295 $ -1, -1, imat, nfail, nerrs, nout )
296 GO TO 120
297 END IF
298*
299* For types 3-5, zero one row and column of the matrix to
300* test that INFO is returned correctly.
301*
302 IF( zerot ) THEN
303 IF( imat.EQ.3 ) THEN
304 izero = 1
305 ELSE IF( imat.EQ.4 ) THEN
306 izero = n
307 ELSE
308 izero = n / 2 + 1
309 END IF
310*
311* Set row and column IZERO of A to 0.
312*
313 IF( iuplo.EQ.1 ) THEN
314 ioff = ( izero-1 )*izero / 2
315 DO 20 i = 1, izero - 1
316 a( ioff+i ) = zero
317 20 CONTINUE
318 ioff = ioff + izero
319 DO 30 i = izero, n
320 a( ioff ) = zero
321 ioff = ioff + i
322 30 CONTINUE
323 ELSE
324 ioff = izero
325 DO 40 i = 1, izero - 1
326 a( ioff ) = zero
327 ioff = ioff + n - i
328 40 CONTINUE
329 ioff = ioff - izero
330 DO 50 i = izero, n
331 a( ioff+i ) = zero
332 50 CONTINUE
333 END IF
334 ELSE
335 izero = 0
336 END IF
337*
338* Set the imaginary part of the diagonals.
339*
340 IF( iuplo.EQ.1 ) THEN
341 CALL claipd( n, a, 2, 1 )
342 ELSE
343 CALL claipd( n, a, n, -1 )
344 END IF
345*
346* Save a copy of the matrix A in ASAV.
347*
348 CALL ccopy( npp, a, 1, asav, 1 )
349*
350 DO 110 iequed = 1, 2
351 equed = equeds( iequed )
352 IF( iequed.EQ.1 ) THEN
353 nfact = 3
354 ELSE
355 nfact = 1
356 END IF
357*
358 DO 100 ifact = 1, nfact
359 fact = facts( ifact )
360 prefac = lsame( fact, 'F' )
361 nofact = lsame( fact, 'N' )
362 equil = lsame( fact, 'E' )
363*
364 IF( zerot ) THEN
365 IF( prefac )
366 $ GO TO 100
367 rcondc = zero
368*
369 ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
370*
371* Compute the condition number for comparison with
372* the value returned by CPPSVX (FACT = 'N' reuses
373* the condition number from the previous iteration
374* with FACT = 'F').
375*
376 CALL ccopy( npp, asav, 1, afac, 1 )
377 IF( equil .OR. iequed.GT.1 ) THEN
378*
379* Compute row and column scale factors to
380* equilibrate the matrix A.
381*
382 CALL cppequ( uplo, n, afac, s, scond, amax,
383 $ info )
384 IF( info.EQ.0 .AND. n.GT.0 ) THEN
385 IF( iequed.GT.1 )
386 $ scond = zero
387*
388* Equilibrate the matrix.
389*
390 CALL claqhp( uplo, n, afac, s, scond,
391 $ amax, equed )
392 END IF
393 END IF
394*
395* Save the condition number of the
396* non-equilibrated system for use in CGET04.
397*
398 IF( equil )
399 $ roldc = rcondc
400*
401* Compute the 1-norm of A.
402*
403 anorm = clanhp( '1', uplo, n, afac, rwork )
404*
405* Factor the matrix A.
406*
407 CALL cpptrf( uplo, n, afac, info )
408*
409* Form the inverse of A.
410*
411 CALL ccopy( npp, afac, 1, a, 1 )
412 CALL cpptri( uplo, n, a, info )
413*
414* Compute the 1-norm condition number of A.
415*
416 ainvnm = clanhp( '1', uplo, n, a, rwork )
417 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
418 rcondc = one
419 ELSE
420 rcondc = ( one / anorm ) / ainvnm
421 END IF
422 END IF
423*
424* Restore the matrix A.
425*
426 CALL ccopy( npp, asav, 1, a, 1 )
427*
428* Form an exact solution and set the right hand side.
429*
430 srnamt = 'CLARHS'
431 CALL clarhs( path, xtype, uplo, ' ', n, n, kl, ku,
432 $ nrhs, a, lda, xact, lda, b, lda,
433 $ iseed, info )
434 xtype = 'C'
435 CALL clacpy( 'Full', n, nrhs, b, lda, bsav, lda )
436*
437 IF( nofact ) THEN
438*
439* --- Test CPPSV ---
440*
441* Compute the L*L' or U'*U factorization of the
442* matrix and solve the system.
443*
444 CALL ccopy( npp, a, 1, afac, 1 )
445 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
446*
447 srnamt = 'CPPSV '
448 CALL cppsv( uplo, n, nrhs, afac, x, lda, info )
449*
450* Check error code from CPPSV .
451*
452 IF( info.NE.izero ) THEN
453 CALL alaerh( path, 'CPPSV ', info, izero,
454 $ uplo, n, n, -1, -1, nrhs, imat,
455 $ nfail, nerrs, nout )
456 GO TO 70
457 ELSE IF( info.NE.0 ) THEN
458 GO TO 70
459 END IF
460*
461* Reconstruct matrix from factors and compute
462* residual.
463*
464 CALL cppt01( uplo, n, a, afac, rwork,
465 $ result( 1 ) )
466*
467* Compute residual of the computed solution.
468*
469 CALL clacpy( 'Full', n, nrhs, b, lda, work,
470 $ lda )
471 CALL cppt02( uplo, n, nrhs, a, x, lda, work,
472 $ lda, rwork, result( 2 ) )
473*
474* Check solution from generated exact solution.
475*
476 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
477 $ result( 3 ) )
478 nt = 3
479*
480* Print information about the tests that did not
481* pass the threshold.
482*
483 DO 60 k = 1, nt
484 IF( result( k ).GE.thresh ) THEN
485 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
486 $ CALL aladhd( nout, path )
487 WRITE( nout, fmt = 9999 )'CPPSV ', uplo,
488 $ n, imat, k, result( k )
489 nfail = nfail + 1
490 END IF
491 60 CONTINUE
492 nrun = nrun + nt
493 70 CONTINUE
494 END IF
495*
496* --- Test CPPSVX ---
497*
498 IF( .NOT.prefac .AND. npp.GT.0 )
499 $ CALL claset( 'Full', npp, 1, cmplx( zero ),
500 $ cmplx( zero ), afac, npp )
501 CALL claset( 'Full', n, nrhs, cmplx( zero ),
502 $ cmplx( zero ), x, lda )
503 IF( iequed.GT.1 .AND. n.GT.0 ) THEN
504*
505* Equilibrate the matrix if FACT='F' and
506* EQUED='Y'.
507*
508 CALL claqhp( uplo, n, a, s, scond, amax, equed )
509 END IF
510*
511* Solve the system and compute the condition number
512* and error bounds using CPPSVX.
513*
514 srnamt = 'CPPSVX'
515 CALL cppsvx( fact, uplo, n, nrhs, a, afac, equed,
516 $ s, b, lda, x, lda, rcond, rwork,
517 $ rwork( nrhs+1 ), work,
518 $ rwork( 2*nrhs+1 ), info )
519*
520* Check the error code from CPPSVX.
521*
522 IF( info.NE.izero ) THEN
523 CALL alaerh( path, 'CPPSVX', info, izero,
524 $ fact // uplo, n, n, -1, -1, nrhs,
525 $ imat, nfail, nerrs, nout )
526 GO TO 90
527 END IF
528*
529 IF( info.EQ.0 ) THEN
530 IF( .NOT.prefac ) THEN
531*
532* Reconstruct matrix from factors and compute
533* residual.
534*
535 CALL cppt01( uplo, n, a, afac,
536 $ rwork( 2*nrhs+1 ), result( 1 ) )
537 k1 = 1
538 ELSE
539 k1 = 2
540 END IF
541*
542* Compute residual of the computed solution.
543*
544 CALL clacpy( 'Full', n, nrhs, bsav, lda, work,
545 $ lda )
546 CALL cppt02( uplo, n, nrhs, asav, x, lda, work,
547 $ lda, rwork( 2*nrhs+1 ),
548 $ result( 2 ) )
549*
550* Check solution from generated exact solution.
551*
552 IF( nofact .OR. ( prefac .AND. lsame( equed,
553 $ 'N' ) ) ) THEN
554 CALL cget04( n, nrhs, x, lda, xact, lda,
555 $ rcondc, result( 3 ) )
556 ELSE
557 CALL cget04( n, nrhs, x, lda, xact, lda,
558 $ roldc, result( 3 ) )
559 END IF
560*
561* Check the error bounds from iterative
562* refinement.
563*
564 CALL cppt05( uplo, n, nrhs, asav, b, lda, x,
565 $ lda, xact, lda, rwork,
566 $ rwork( nrhs+1 ), result( 4 ) )
567 ELSE
568 k1 = 6
569 END IF
570*
571* Compare RCOND from CPPSVX with the computed value
572* in RCONDC.
573*
574 result( 6 ) = sget06( rcond, rcondc )
575*
576* Print information about the tests that did not pass
577* the threshold.
578*
579 DO 80 k = k1, 6
580 IF( result( k ).GE.thresh ) THEN
581 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
582 $ CALL aladhd( nout, path )
583 IF( prefac ) THEN
584 WRITE( nout, fmt = 9997 )'CPPSVX', fact,
585 $ uplo, n, equed, imat, k, result( k )
586 ELSE
587 WRITE( nout, fmt = 9998 )'CPPSVX', fact,
588 $ uplo, n, imat, k, result( k )
589 END IF
590 nfail = nfail + 1
591 END IF
592 80 CONTINUE
593 nrun = nrun + 7 - k1
594 90 CONTINUE
595 100 CONTINUE
596 110 CONTINUE
597 120 CONTINUE
598 130 CONTINUE
599 140 CONTINUE
600*
601* Print a summary of the results.
602*
603 CALL alasvm( path, nout, nfail, nrun, nerrs )
604*
605 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
606 $ ', test(', i1, ')=', g12.5 )
607 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
608 $ ', type ', i1, ', test(', i1, ')=', g12.5 )
609 9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
610 $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ')=',
611 $ g12.5 )
612 RETURN
613*
614* End of CDRVPP
615*
616 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine clarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
CLARHS
Definition clarhs.f:208
subroutine aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine cdrvpp(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, asav, b, bsav, x, xact, s, work, rwork, nout)
CDRVPP
Definition cdrvpp.f:159
subroutine cerrvx(path, nunit)
CERRVX
Definition cerrvx.f:55
subroutine cget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
CGET04
Definition cget04.f:102
subroutine claipd(n, a, inda, vinda)
CLAIPD
Definition claipd.f:83
subroutine clatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
CLATB4
Definition clatb4.f:121
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine cppt01(uplo, n, a, afac, rwork, resid)
CPPT01
Definition cppt01.f:95
subroutine cppt02(uplo, n, nrhs, a, x, ldx, b, ldb, rwork, resid)
CPPT02
Definition cppt02.f:123
subroutine cppt05(uplo, n, nrhs, ap, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
CPPT05
Definition cppt05.f:157
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
subroutine claqhp(uplo, n, ap, s, scond, amax, equed)
CLAQHP scales a Hermitian matrix stored in packed form.
Definition claqhp.f:126
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:106
subroutine cppequ(uplo, n, ap, s, scond, amax, info)
CPPEQU
Definition cppequ.f:117
subroutine cppsv(uplo, n, nrhs, ap, b, ldb, info)
CPPSV computes the solution to system of linear equations A * X = B for OTHER matrices
Definition cppsv.f:144
subroutine cppsvx(fact, uplo, n, nrhs, ap, afp, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
CPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices
Definition cppsvx.f:311
subroutine cpptrf(uplo, n, ap, info)
CPPTRF
Definition cpptrf.f:119
subroutine cpptri(uplo, n, ap, info)
CPPTRI
Definition cpptri.f:93