LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zchksy()

subroutine zchksy ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
integer  nns,
integer, dimension( * )  nsval,
double precision  thresh,
logical  tsterr,
integer  nmax,
complex*16, dimension( * )  a,
complex*16, dimension( * )  afac,
complex*16, dimension( * )  ainv,
complex*16, dimension( * )  b,
complex*16, dimension( * )  x,
complex*16, dimension( * )  xact,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

ZCHKSY

Purpose:
 ZCHKSY tests ZSYTRF, -TRI2, -TRS, -TRS2,  -RFS, and -CON.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix dimension N.
[in]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the blocksize NB.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension (NMAX*max(2,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NSMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 168 of file zchksy.f.

171*
172* -- LAPACK test routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 LOGICAL TSTERR
178 INTEGER NMAX, NN, NNB, NNS, NOUT
179 DOUBLE PRECISION THRESH
180* ..
181* .. Array Arguments ..
182 LOGICAL DOTYPE( * )
183 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
184 DOUBLE PRECISION RWORK( * )
185 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
186 $ WORK( * ), X( * ), XACT( * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 DOUBLE PRECISION ZERO
193 parameter( zero = 0.0d+0 )
194 COMPLEX*16 CZERO
195 parameter( czero = ( 0.0d+0, 0.0d+0 ) )
196 INTEGER NTYPES
197 parameter( ntypes = 11 )
198 INTEGER NTESTS
199 parameter( ntests = 9 )
200* ..
201* .. Local Scalars ..
202 LOGICAL TRFCON, ZEROT
203 CHARACTER DIST, TYPE, UPLO, XTYPE
204 CHARACTER*3 PATH
205 INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
206 $ IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
207 $ N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
208 DOUBLE PRECISION ANORM, CNDNUM, RCOND, RCONDC
209* ..
210* .. Local Arrays ..
211 CHARACTER UPLOS( 2 )
212 INTEGER ISEED( 4 ), ISEEDY( 4 )
213 DOUBLE PRECISION RESULT( NTESTS )
214* ..
215* .. External Functions ..
216 DOUBLE PRECISION DGET06, ZLANSY
217 EXTERNAL dget06, zlansy
218* ..
219* .. External Subroutines ..
220 EXTERNAL alaerh, alahd, alasum, xlaenv, zerrsy, zget04,
224* ..
225* .. Intrinsic Functions ..
226 INTRINSIC max, min
227* ..
228* .. Scalars in Common ..
229 LOGICAL LERR, OK
230 CHARACTER*32 SRNAMT
231 INTEGER INFOT, NUNIT
232* ..
233* .. Common blocks ..
234 COMMON / infoc / infot, nunit, ok, lerr
235 COMMON / srnamc / srnamt
236* ..
237* .. Data statements ..
238 DATA iseedy / 1988, 1989, 1990, 1991 /
239 DATA uplos / 'U', 'L' /
240* ..
241* .. Executable Statements ..
242*
243* Initialize constants and the random number seed.
244*
245 path( 1: 1 ) = 'Zomplex precision'
246 path( 2: 3 ) = 'SY'
247 nrun = 0
248 nfail = 0
249 nerrs = 0
250 DO 10 i = 1, 4
251 iseed( i ) = iseedy( i )
252 10 CONTINUE
253*
254* Test the error exits
255*
256 IF( tsterr )
257 $ CALL zerrsy( path, nout )
258 infot = 0
259*
260* Set the minimum block size for which the block routine should
261* be used, which will be later returned by ILAENV
262*
263 CALL xlaenv( 2, 2 )
264*
265* Do for each value of N in NVAL
266*
267 DO 180 in = 1, nn
268 n = nval( in )
269 lda = max( n, 1 )
270 xtype = 'N'
271 nimat = ntypes
272 IF( n.LE.0 )
273 $ nimat = 1
274*
275 izero = 0
276*
277* Do for each value of matrix type IMAT
278*
279 DO 170 imat = 1, nimat
280*
281* Do the tests only if DOTYPE( IMAT ) is true.
282*
283 IF( .NOT.dotype( imat ) )
284 $ GO TO 170
285*
286* Skip types 3, 4, 5, or 6 if the matrix size is too small.
287*
288 zerot = imat.GE.3 .AND. imat.LE.6
289 IF( zerot .AND. n.LT.imat-2 )
290 $ GO TO 170
291*
292* Do first for UPLO = 'U', then for UPLO = 'L'
293*
294 DO 160 iuplo = 1, 2
295 uplo = uplos( iuplo )
296*
297* Begin generate test matrix A.
298*
299 IF( imat.NE.ntypes ) THEN
300*
301* Set up parameters with ZLATB4 for the matrix generator
302* based on the type of matrix to be generated.
303*
304 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM,
305 $ MODE, CNDNUM, DIST )
306*
307* Generate a matrix with ZLATMS.
308*
309 srnamt = 'ZLATMS'
310 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
311 $ CNDNUM, ANORM, KL, KU, 'N', A, LDA, WORK,
312 $ INFO )
313*
314* Check error code from ZLATMS and handle error.
315*
316 IF( info.NE.0 ) THEN
317 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
318 $ -1, -1, -1, imat, nfail, nerrs, nout )
319*
320* Skip all tests for this generated matrix
321*
322 GO TO 160
323 END IF
324*
325* For matrix types 3-6, zero one or more rows and
326* columns of the matrix to test that INFO is returned
327* correctly.
328*
329 IF( zerot ) THEN
330 IF( imat.EQ.3 ) THEN
331 izero = 1
332 ELSE IF( imat.EQ.4 ) THEN
333 izero = n
334 ELSE
335 izero = n / 2 + 1
336 END IF
337*
338 IF( imat.LT.6 ) THEN
339*
340* Set row and column IZERO to zero.
341*
342 IF( iuplo.EQ.1 ) THEN
343 ioff = ( izero-1 )*lda
344 DO 20 i = 1, izero - 1
345 a( ioff+i ) = czero
346 20 CONTINUE
347 ioff = ioff + izero
348 DO 30 i = izero, n
349 a( ioff ) = czero
350 ioff = ioff + lda
351 30 CONTINUE
352 ELSE
353 ioff = izero
354 DO 40 i = 1, izero - 1
355 a( ioff ) = czero
356 ioff = ioff + lda
357 40 CONTINUE
358 ioff = ioff - izero
359 DO 50 i = izero, n
360 a( ioff+i ) = czero
361 50 CONTINUE
362 END IF
363 ELSE
364 IF( iuplo.EQ.1 ) THEN
365*
366* Set the first IZERO rows to zero.
367*
368 ioff = 0
369 DO 70 j = 1, n
370 i2 = min( j, izero )
371 DO 60 i = 1, i2
372 a( ioff+i ) = czero
373 60 CONTINUE
374 ioff = ioff + lda
375 70 CONTINUE
376 ELSE
377*
378* Set the last IZERO rows to zero.
379*
380 ioff = 0
381 DO 90 j = 1, n
382 i1 = max( j, izero )
383 DO 80 i = i1, n
384 a( ioff+i ) = czero
385 80 CONTINUE
386 ioff = ioff + lda
387 90 CONTINUE
388 END IF
389 END IF
390 ELSE
391 izero = 0
392 END IF
393*
394 ELSE
395*
396* For matrix kind IMAT = 11, generate special block
397* diagonal matrix to test alternate code
398* for the 2 x 2 blocks.
399*
400 CALL zlatsy( uplo, n, a, lda, iseed )
401*
402 END IF
403*
404* End generate test matrix A.
405*
406*
407* Do for each value of NB in NBVAL
408*
409 DO 150 inb = 1, nnb
410*
411* Set the optimal blocksize, which will be later
412* returned by ILAENV.
413*
414 nb = nbval( inb )
415 CALL xlaenv( 1, nb )
416*
417* Copy the test matrix A into matrix AFAC which
418* will be factorized in place. This is needed to
419* preserve the test matrix A for subsequent tests.
420*
421 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
422*
423* Compute the L*D*L**T or U*D*U**T factorization of the
424* matrix. IWORK stores details of the interchanges and
425* the block structure of D. AINV is a work array for
426* block factorization, LWORK is the length of AINV.
427*
428 lwork = max( 2, nb )*lda
429 srnamt = 'ZSYTRF'
430 CALL zsytrf( uplo, n, afac, lda, iwork, ainv, lwork,
431 $ info )
432*
433* Adjust the expected value of INFO to account for
434* pivoting.
435*
436 k = izero
437 IF( k.GT.0 ) THEN
438 100 CONTINUE
439 IF( iwork( k ).LT.0 ) THEN
440 IF( iwork( k ).NE.-k ) THEN
441 k = -iwork( k )
442 GO TO 100
443 END IF
444 ELSE IF( iwork( k ).NE.k ) THEN
445 k = iwork( k )
446 GO TO 100
447 END IF
448 END IF
449*
450* Check error code from ZSYTRF and handle error.
451*
452 IF( info.NE.k )
453 $ CALL alaerh( path, 'ZSYTRF', info, k, uplo, n, n,
454 $ -1, -1, nb, imat, nfail, nerrs, nout )
455*
456* Set the condition estimate flag if the INFO is not 0.
457*
458 IF( info.NE.0 ) THEN
459 trfcon = .true.
460 ELSE
461 trfcon = .false.
462 END IF
463*
464*+ TEST 1
465* Reconstruct matrix from factors and compute residual.
466*
467 CALL zsyt01( uplo, n, a, lda, afac, lda, iwork, ainv,
468 $ lda, rwork, result( 1 ) )
469 nt = 1
470*
471*+ TEST 2
472* Form the inverse and compute the residual,
473* if the factorization was competed without INFO > 0
474* (i.e. there is no zero rows and columns).
475* Do it only for the first block size.
476*
477 IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
478 CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
479 srnamt = 'ZSYTRI2'
480 lwork = (n+nb+1)*(nb+3)
481 CALL zsytri2( uplo, n, ainv, lda, iwork, work,
482 $ lwork, info )
483*
484* Check error code from ZSYTRI2 and handle error.
485*
486 IF( info.NE.0 )
487 $ CALL alaerh( path, 'ZSYTRI2', info, 0, uplo, n,
488 $ n, -1, -1, -1, imat, nfail, nerrs,
489 $ nout )
490*
491* Compute the residual for a symmetric matrix times
492* its inverse.
493*
494 CALL zsyt03( uplo, n, a, lda, ainv, lda, work, lda,
495 $ rwork, rcondc, result( 2 ) )
496 nt = 2
497 END IF
498*
499* Print information about the tests that did not pass
500* the threshold.
501*
502 DO 110 k = 1, nt
503 IF( result( k ).GE.thresh ) THEN
504 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
505 $ CALL alahd( nout, path )
506 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
507 $ result( k )
508 nfail = nfail + 1
509 END IF
510 110 CONTINUE
511 nrun = nrun + nt
512*
513* Skip the other tests if this is not the first block
514* size.
515*
516 IF( inb.GT.1 )
517 $ GO TO 150
518*
519* Do only the condition estimate if INFO is not 0.
520*
521 IF( trfcon ) THEN
522 rcondc = zero
523 GO TO 140
524 END IF
525*
526* Do for each value of NRHS in NSVAL.
527*
528 DO 130 irhs = 1, nns
529 nrhs = nsval( irhs )
530*
531*+ TEST 3 (Using TRS)
532* Solve and compute residual for A * X = B.
533*
534* Choose a set of NRHS random solution vectors
535* stored in XACT and set up the right hand side B
536*
537 srnamt = 'ZLARHS'
538 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
539 $ nrhs, a, lda, xact, lda, b, lda,
540 $ iseed, info )
541 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
542*
543 srnamt = 'ZSYTRS'
544 CALL zsytrs( uplo, n, nrhs, afac, lda, iwork, x,
545 $ lda, info )
546*
547* Check error code from ZSYTRS and handle error.
548*
549 IF( info.NE.0 )
550 $ CALL alaerh( path, 'ZSYTRS', info, 0, uplo, n,
551 $ n, -1, -1, nrhs, imat, nfail,
552 $ nerrs, nout )
553*
554 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
555*
556* Compute the residual for the solution
557*
558 CALL zsyt02( uplo, n, nrhs, a, lda, x, lda, work,
559 $ lda, rwork, result( 3 ) )
560*
561*+ TEST 4 (Using TRS2)
562* Solve and compute residual for A * X = B.
563*
564* Choose a set of NRHS random solution vectors
565* stored in XACT and set up the right hand side B
566*
567 srnamt = 'ZLARHS'
568 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
569 $ nrhs, a, lda, xact, lda, b, lda,
570 $ iseed, info )
571 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
572*
573 srnamt = 'ZSYTRS2'
574 CALL zsytrs2( uplo, n, nrhs, afac, lda, iwork, x,
575 $ lda, work, info )
576*
577* Check error code from ZSYTRS2 and handle error.
578*
579 IF( info.NE.0 )
580 $ CALL alaerh( path, 'ZSYTRS', info, 0, uplo, n,
581 $ n, -1, -1, nrhs, imat, nfail,
582 $ nerrs, nout )
583*
584 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
585*
586* Compute the residual for the solution
587*
588 CALL zsyt02( uplo, n, nrhs, a, lda, x, lda, work,
589 $ lda, rwork, result( 4 ) )
590*
591*
592*+ TEST 5
593* Check solution from generated exact solution.
594*
595 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
596 $ result( 5 ) )
597*
598*+ TESTS 6, 7, and 8
599* Use iterative refinement to improve the solution.
600*
601 srnamt = 'ZSYRFS'
602 CALL zsyrfs( uplo, n, nrhs, a, lda, afac, lda,
603 $ iwork, b, lda, x, lda, rwork,
604 $ rwork( nrhs+1 ), work,
605 $ rwork( 2*nrhs+1 ), info )
606*
607* Check error code from ZSYRFS and handle error.
608*
609 IF( info.NE.0 )
610 $ CALL alaerh( path, 'ZSYRFS', info, 0, uplo, n,
611 $ n, -1, -1, nrhs, imat, nfail,
612 $ nerrs, nout )
613*
614 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
615 $ result( 6 ) )
616 CALL zpot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
617 $ xact, lda, rwork, rwork( nrhs+1 ),
618 $ result( 7 ) )
619*
620* Print information about the tests that did not pass
621* the threshold.
622*
623 DO 120 k = 3, 8
624 IF( result( k ).GE.thresh ) THEN
625 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
626 $ CALL alahd( nout, path )
627 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
628 $ imat, k, result( k )
629 nfail = nfail + 1
630 END IF
631 120 CONTINUE
632 nrun = nrun + 6
633*
634* End do for each value of NRHS in NSVAL.
635*
636 130 CONTINUE
637*
638*+ TEST 9
639* Get an estimate of RCOND = 1/CNDNUM.
640*
641 140 CONTINUE
642 anorm = zlansy( '1', uplo, n, a, lda, rwork )
643 srnamt = 'ZSYCON'
644 CALL zsycon( uplo, n, afac, lda, iwork, anorm, rcond,
645 $ work, info )
646*
647* Check error code from ZSYCON and handle error.
648*
649 IF( info.NE.0 )
650 $ CALL alaerh( path, 'ZSYCON', info, 0, uplo, n, n,
651 $ -1, -1, -1, imat, nfail, nerrs, nout )
652*
653* Compute the test ratio to compare values of RCOND
654*
655 result( 9 ) = dget06( rcond, rcondc )
656*
657* Print information about the tests that did not pass
658* the threshold.
659*
660 IF( result( 9 ).GE.thresh ) THEN
661 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
662 $ CALL alahd( nout, path )
663 WRITE( nout, fmt = 9997 )uplo, n, imat, 9,
664 $ result( 9 )
665 nfail = nfail + 1
666 END IF
667 nrun = nrun + 1
668 150 CONTINUE
669 160 CONTINUE
670 170 CONTINUE
671 180 CONTINUE
672*
673* Print a summary of the results.
674*
675 CALL alasum( path, nout, nfail, nrun, nerrs )
676*
677 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
678 $ i2, ', test ', i2, ', ratio =', g12.5 )
679 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
680 $ i2, ', test(', i2, ') =', g12.5 )
681 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
682 $ ', test(', i2, ') =', g12.5 )
683 RETURN
684*
685* End of ZCHKSY
686*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.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 alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
subroutine zsycon(uplo, n, a, lda, ipiv, anorm, rcond, work, info)
ZSYCON
Definition zsycon.f:125
subroutine zsyrfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZSYRFS
Definition zsyrfs.f:192
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 zsytrs2(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, info)
ZSYTRS2
Definition zsytrs2.f:132
subroutine zsytrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
ZSYTRS
Definition zsytrs.f:120
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
double precision function zlansy(norm, uplo, n, a, lda, work)
ZLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlansy.f:123
subroutine zerrsy(path, nunit)
ZERRSY
Definition zerrsy.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
subroutine zsyt03(uplo, n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
ZSYT03
Definition zsyt03.f:126
Here is the call graph for this function:
Here is the caller graph for this function: