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