LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvsp.f
Go to the documentation of this file.
1*> \brief \b ZDRVSP
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 ZDRVSP( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12* A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK,
13* 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 IWORK( * ), NVAL( * )
23* DOUBLE PRECISION RWORK( * )
24* COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZDRVSP tests the driver routines ZSPSV 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
91*> (NMAX*(NMAX+1)/2)
92*> \endverbatim
93*>
94*> \param[out] AFAC
95*> \verbatim
96*> AFAC is COMPLEX*16 array, dimension
97*> (NMAX*(NMAX+1)/2)
98*> \endverbatim
99*>
100*> \param[out] AINV
101*> \verbatim
102*> AINV is COMPLEX*16 array, dimension
103*> (NMAX*(NMAX+1)/2)
104*> \endverbatim
105*>
106*> \param[out] B
107*> \verbatim
108*> B is COMPLEX*16 array, dimension (NMAX*NRHS)
109*> \endverbatim
110*>
111*> \param[out] X
112*> \verbatim
113*> X is COMPLEX*16 array, dimension (NMAX*NRHS)
114*> \endverbatim
115*>
116*> \param[out] XACT
117*> \verbatim
118*> XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
119*> \endverbatim
120*>
121*> \param[out] WORK
122*> \verbatim
123*> WORK is COMPLEX*16 array, dimension
124*> (NMAX*max(2,NRHS))
125*> \endverbatim
126*>
127*> \param[out] RWORK
128*> \verbatim
129*> RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
130*> \endverbatim
131*>
132*> \param[out] IWORK
133*> \verbatim
134*> IWORK is INTEGER array, dimension (NMAX)
135*> \endverbatim
136*>
137*> \param[in] NOUT
138*> \verbatim
139*> NOUT is INTEGER
140*> The unit number for output.
141*> \endverbatim
142*
143* Authors:
144* ========
145*
146*> \author Univ. of Tennessee
147*> \author Univ. of California Berkeley
148*> \author Univ. of Colorado Denver
149*> \author NAG Ltd.
150*
151*> \ingroup complex16_lin
152*
153* =====================================================================
154 SUBROUTINE zdrvsp( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
155 $ A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK,
156 $ NOUT )
157*
158* -- LAPACK test routine --
159* -- LAPACK is a software package provided by Univ. of Tennessee, --
160* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
161*
162* .. Scalar Arguments ..
163 LOGICAL TSTERR
164 INTEGER NMAX, NN, NOUT, NRHS
165 DOUBLE PRECISION THRESH
166* ..
167* .. Array Arguments ..
168 LOGICAL DOTYPE( * )
169 INTEGER IWORK( * ), NVAL( * )
170 DOUBLE PRECISION RWORK( * )
171 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
172 $ work( * ), x( * ), xact( * )
173* ..
174*
175* =====================================================================
176*
177* .. Parameters ..
178 DOUBLE PRECISION ONE, ZERO
179 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
180 INTEGER NTYPES, NTESTS
181 parameter( ntypes = 11, ntests = 6 )
182 INTEGER NFACT
183 parameter( nfact = 2 )
184* ..
185* .. Local Scalars ..
186 LOGICAL ZEROT
187 CHARACTER DIST, FACT, PACKIT, TYPE, UPLO, XTYPE
188 CHARACTER*3 PATH
189 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
190 $ izero, j, k, k1, kl, ku, lda, mode, n, nb,
191 $ nbmin, nerrs, nfail, nimat, npp, nrun, nt
192 DOUBLE PRECISION AINVNM, ANORM, CNDNUM, RCOND, RCONDC
193* ..
194* .. Local Arrays ..
195 CHARACTER FACTS( NFACT )
196 INTEGER ISEED( 4 ), ISEEDY( 4 )
197 DOUBLE PRECISION RESULT( NTESTS )
198* ..
199* .. External Functions ..
200 DOUBLE PRECISION DGET06, ZLANSP
201 EXTERNAL DGET06, ZLANSP
202* ..
203* .. External Subroutines ..
204 EXTERNAL aladhd, alaerh, alasvm, xlaenv, zcopy, zerrvx,
207 $ zsptrf, zsptri
208* ..
209* .. Scalars in Common ..
210 LOGICAL LERR, OK
211 CHARACTER*32 SRNAMT
212 INTEGER INFOT, NUNIT
213* ..
214* .. Common blocks ..
215 COMMON / infoc / infot, nunit, ok, lerr
216 COMMON / srnamc / srnamt
217* ..
218* .. Intrinsic Functions ..
219 INTRINSIC dcmplx, max, min
220* ..
221* .. Data statements ..
222 DATA iseedy / 1988, 1989, 1990, 1991 /
223 DATA facts / 'F', 'N' /
224* ..
225* .. Executable Statements ..
226*
227* Initialize constants and the random number seed.
228*
229 path( 1: 1 ) = 'Zomplex precision'
230 path( 2: 3 ) = 'SP'
231 nrun = 0
232 nfail = 0
233 nerrs = 0
234 DO 10 i = 1, 4
235 iseed( i ) = iseedy( i )
236 10 CONTINUE
237*
238* Test the error exits
239*
240 IF( tsterr )
241 $ CALL zerrvx( path, nout )
242 infot = 0
243*
244* Set the block size and minimum block size for testing.
245*
246 nb = 1
247 nbmin = 2
248 CALL xlaenv( 1, nb )
249 CALL xlaenv( 2, nbmin )
250*
251* Do for each value of N in NVAL
252*
253 DO 180 in = 1, nn
254 n = nval( in )
255 lda = max( n, 1 )
256 npp = n*( n+1 ) / 2
257 xtype = 'N'
258 nimat = ntypes
259 IF( n.LE.0 )
260 $ nimat = 1
261*
262 DO 170 imat = 1, nimat
263*
264* Do the tests only if DOTYPE( IMAT ) is true.
265*
266 IF( .NOT.dotype( imat ) )
267 $ GO TO 170
268*
269* Skip types 3, 4, 5, or 6 if the matrix size is too small.
270*
271 zerot = imat.GE.3 .AND. imat.LE.6
272 IF( zerot .AND. n.LT.imat-2 )
273 $ GO TO 170
274*
275* Do first for UPLO = 'U', then for UPLO = 'L'
276*
277 DO 160 iuplo = 1, 2
278 IF( iuplo.EQ.1 ) THEN
279 uplo = 'U'
280 packit = 'C'
281 ELSE
282 uplo = 'L'
283 packit = 'R'
284 END IF
285*
286 IF( imat.NE.ntypes ) THEN
287*
288* Set up parameters with ZLATB4 and generate a test
289* matrix with ZLATMS.
290*
291 CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm,
292 $ mode, cndnum, dist )
293*
294 srnamt = 'ZLATMS'
295 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
296 $ cndnum, anorm, kl, ku, packit, a, lda,
297 $ work, 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,
303 $ -1, -1, -1, imat, nfail, nerrs, nout )
304 GO TO 160
305 END IF
306*
307* For types 3-6, zero one or more rows and columns of
308* the matrix to 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*
319 IF( imat.LT.6 ) THEN
320*
321* Set row and column IZERO to zero.
322*
323 IF( iuplo.EQ.1 ) THEN
324 ioff = ( izero-1 )*izero / 2
325 DO 20 i = 1, izero - 1
326 a( ioff+i ) = zero
327 20 CONTINUE
328 ioff = ioff + izero
329 DO 30 i = izero, n
330 a( ioff ) = zero
331 ioff = ioff + i
332 30 CONTINUE
333 ELSE
334 ioff = izero
335 DO 40 i = 1, izero - 1
336 a( ioff ) = zero
337 ioff = ioff + n - i
338 40 CONTINUE
339 ioff = ioff - izero
340 DO 50 i = izero, n
341 a( ioff+i ) = zero
342 50 CONTINUE
343 END IF
344 ELSE
345 IF( iuplo.EQ.1 ) THEN
346*
347* Set the first IZERO rows and columns to zero.
348*
349 ioff = 0
350 DO 70 j = 1, n
351 i2 = min( j, izero )
352 DO 60 i = 1, i2
353 a( ioff+i ) = zero
354 60 CONTINUE
355 ioff = ioff + j
356 70 CONTINUE
357 ELSE
358*
359* Set the last IZERO rows and columns to zero.
360*
361 ioff = 0
362 DO 90 j = 1, n
363 i1 = max( j, izero )
364 DO 80 i = i1, n
365 a( ioff+i ) = zero
366 80 CONTINUE
367 ioff = ioff + n - j
368 90 CONTINUE
369 END IF
370 END IF
371 ELSE
372 izero = 0
373 END IF
374 ELSE
375*
376* Use a special block diagonal matrix to test alternate
377* code for the 2-by-2 blocks.
378*
379 CALL zlatsp( uplo, n, a, iseed )
380 END IF
381*
382 DO 150 ifact = 1, nfact
383*
384* Do first for FACT = 'F', then for other values.
385*
386 fact = facts( ifact )
387*
388* Compute the condition number for comparison with
389* the value returned by ZSPSVX.
390*
391 IF( zerot ) THEN
392 IF( ifact.EQ.1 )
393 $ GO TO 150
394 rcondc = zero
395*
396 ELSE IF( ifact.EQ.1 ) THEN
397*
398* Compute the 1-norm of A.
399*
400 anorm = zlansp( '1', uplo, n, a, rwork )
401*
402* Factor the matrix A.
403*
404 CALL zcopy( npp, a, 1, afac, 1 )
405 CALL zsptrf( uplo, n, afac, iwork, info )
406*
407* Compute inv(A) and take its norm.
408*
409 CALL zcopy( npp, afac, 1, ainv, 1 )
410 CALL zsptri( uplo, n, ainv, iwork, work, info )
411 ainvnm = zlansp( '1', uplo, n, ainv, rwork )
412*
413* Compute the 1-norm condition number of A.
414*
415 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
416 rcondc = one
417 ELSE
418 rcondc = ( one / anorm ) / ainvnm
419 END IF
420 END IF
421*
422* Form an exact solution and set the right hand side.
423*
424 srnamt = 'ZLARHS'
425 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
426 $ nrhs, a, lda, xact, lda, b, lda, iseed,
427 $ info )
428 xtype = 'C'
429*
430* --- Test ZSPSV ---
431*
432 IF( ifact.EQ.2 ) THEN
433 CALL zcopy( npp, a, 1, afac, 1 )
434 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
435*
436* Factor the matrix and solve the system using ZSPSV.
437*
438 srnamt = 'ZSPSV '
439 CALL zspsv( uplo, n, nrhs, afac, iwork, x, lda,
440 $ info )
441*
442* Adjust the expected value of INFO to account for
443* pivoting.
444*
445 k = izero
446 IF( k.GT.0 ) THEN
447 100 CONTINUE
448 IF( iwork( k ).LT.0 ) THEN
449 IF( iwork( k ).NE.-k ) THEN
450 k = -iwork( k )
451 GO TO 100
452 END IF
453 ELSE IF( iwork( k ).NE.k ) THEN
454 k = iwork( k )
455 GO TO 100
456 END IF
457 END IF
458*
459* Check error code from ZSPSV .
460*
461 IF( info.NE.k ) THEN
462 CALL alaerh( path, 'ZSPSV ', info, k, uplo, n,
463 $ n, -1, -1, nrhs, imat, nfail,
464 $ nerrs, nout )
465 GO TO 120
466 ELSE IF( info.NE.0 ) THEN
467 GO TO 120
468 END IF
469*
470* Reconstruct matrix from factors and compute
471* residual.
472*
473 CALL zspt01( uplo, n, a, afac, iwork, ainv, lda,
474 $ rwork, result( 1 ) )
475*
476* Compute residual of the computed solution.
477*
478 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
479 CALL zspt02( uplo, n, nrhs, a, x, lda, work, lda,
480 $ rwork, result( 2 ) )
481*
482* Check solution from generated exact solution.
483*
484 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
485 $ result( 3 ) )
486 nt = 3
487*
488* Print information about the tests that did not pass
489* the threshold.
490*
491 DO 110 k = 1, nt
492 IF( result( k ).GE.thresh ) THEN
493 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
494 $ CALL aladhd( nout, path )
495 WRITE( nout, fmt = 9999 )'ZSPSV ', uplo, n,
496 $ imat, k, result( k )
497 nfail = nfail + 1
498 END IF
499 110 CONTINUE
500 nrun = nrun + nt
501 120 CONTINUE
502 END IF
503*
504* --- Test ZSPSVX ---
505*
506 IF( ifact.EQ.2 .AND. npp.GT.0 )
507 $ CALL zlaset( 'Full', npp, 1, dcmplx( zero ),
508 $ dcmplx( zero ), afac, npp )
509 CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
510 $ dcmplx( zero ), x, lda )
511*
512* Solve the system and compute the condition number and
513* error bounds using ZSPSVX.
514*
515 srnamt = 'ZSPSVX'
516 CALL zspsvx( fact, uplo, n, nrhs, a, afac, iwork, b,
517 $ lda, x, lda, rcond, rwork,
518 $ rwork( nrhs+1 ), work, rwork( 2*nrhs+1 ),
519 $ info )
520*
521* Adjust the expected value of INFO to account for
522* pivoting.
523*
524 k = izero
525 IF( k.GT.0 ) THEN
526 130 CONTINUE
527 IF( iwork( k ).LT.0 ) THEN
528 IF( iwork( k ).NE.-k ) THEN
529 k = -iwork( k )
530 GO TO 130
531 END IF
532 ELSE IF( iwork( k ).NE.k ) THEN
533 k = iwork( k )
534 GO TO 130
535 END IF
536 END IF
537*
538* Check the error code from ZSPSVX.
539*
540 IF( info.NE.k ) THEN
541 CALL alaerh( path, 'ZSPSVX', info, k, fact // uplo,
542 $ n, n, -1, -1, nrhs, imat, nfail,
543 $ nerrs, nout )
544 GO TO 150
545 END IF
546*
547 IF( info.EQ.0 ) THEN
548 IF( ifact.GE.2 ) THEN
549*
550* Reconstruct matrix from factors and compute
551* residual.
552*
553 CALL zspt01( uplo, n, a, afac, iwork, ainv, lda,
554 $ rwork( 2*nrhs+1 ), result( 1 ) )
555 k1 = 1
556 ELSE
557 k1 = 2
558 END IF
559*
560* Compute residual of the computed solution.
561*
562 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
563 CALL zspt02( uplo, n, nrhs, a, x, lda, work, lda,
564 $ rwork( 2*nrhs+1 ), result( 2 ) )
565*
566* Check solution from generated exact solution.
567*
568 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
569 $ result( 3 ) )
570*
571* Check the error bounds from iterative refinement.
572*
573 CALL zppt05( uplo, n, nrhs, a, b, lda, x, lda,
574 $ xact, lda, rwork, rwork( nrhs+1 ),
575 $ result( 4 ) )
576 ELSE
577 k1 = 6
578 END IF
579*
580* Compare RCOND from ZSPSVX with the computed value
581* in RCONDC.
582*
583 result( 6 ) = dget06( rcond, rcondc )
584*
585* Print information about the tests that did not pass
586* the threshold.
587*
588 DO 140 k = k1, 6
589 IF( result( k ).GE.thresh ) THEN
590 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
591 $ CALL aladhd( nout, path )
592 WRITE( nout, fmt = 9998 )'ZSPSVX', fact, uplo,
593 $ n, imat, k, result( k )
594 nfail = nfail + 1
595 END IF
596 140 CONTINUE
597 nrun = nrun + 7 - k1
598*
599 150 CONTINUE
600*
601 160 CONTINUE
602 170 CONTINUE
603 180 CONTINUE
604*
605* Print a summary of the results.
606*
607 CALL alasvm( path, nout, nfail, nrun, nerrs )
608*
609 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
610 $ ', test ', i2, ', ratio =', g12.5 )
611 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N =', i5,
612 $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
613 RETURN
614*
615* End of ZDRVSP
616*
617 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.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 zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zspsv(uplo, n, nrhs, ap, ipiv, b, ldb, info)
ZSPSV computes the solution to system of linear equations A * X = B for OTHER matrices
Definition zspsv.f:162
subroutine zspsvx(fact, uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZSPSVX computes the solution to system of linear equations A * X = B for OTHER matrices
Definition zspsvx.f:277
subroutine zsptrf(uplo, n, ap, ipiv, info)
ZSPTRF
Definition zsptrf.f:158
subroutine zsptri(uplo, n, ap, ipiv, work, info)
ZSPTRI
Definition zsptri.f:109
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zdrvsp(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
ZDRVSP
Definition zdrvsp.f:157
subroutine zerrvx(path, nunit)
ZERRVX
Definition zerrvx.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332
subroutine zlatsp(uplo, n, x, iseed)
ZLATSP
Definition zlatsp.f:84
subroutine zppt05(uplo, n, nrhs, ap, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPPT05
Definition zppt05.f:157
subroutine zspt01(uplo, n, a, afac, ipiv, c, ldc, rwork, resid)
ZSPT01
Definition zspt01.f:112
subroutine zspt02(uplo, n, nrhs, a, x, ldx, b, ldb, rwork, resid)
ZSPT02
Definition zspt02.f:123