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

◆ zchksp()

subroutine zchksp ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
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 
)

ZCHKSP

Purpose:
 ZCHKSP tests ZSPTRF, -TRI, -TRS, -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]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+1)/2)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension
                      (NMAX*(NMAX+1)/2)
[out]AINV
          AINV is COMPLEX*16 array, dimension
                      (NMAX*(NMAX+1)/2)
[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 161 of file zchksp.f.

164*
165* -- LAPACK test routine --
166* -- LAPACK is a software package provided by Univ. of Tennessee, --
167* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
168*
169* .. Scalar Arguments ..
170 LOGICAL TSTERR
171 INTEGER NMAX, NN, NNS, NOUT
172 DOUBLE PRECISION THRESH
173* ..
174* .. Array Arguments ..
175 LOGICAL DOTYPE( * )
176 INTEGER IWORK( * ), NSVAL( * ), NVAL( * )
177 DOUBLE PRECISION RWORK( * )
178 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
179 $ WORK( * ), X( * ), XACT( * )
180* ..
181*
182* =====================================================================
183*
184* .. Parameters ..
185 DOUBLE PRECISION ZERO
186 parameter( zero = 0.0d+0 )
187 INTEGER NTYPES
188 parameter( ntypes = 11 )
189 INTEGER NTESTS
190 parameter( ntests = 8 )
191* ..
192* .. Local Scalars ..
193 LOGICAL TRFCON, ZEROT
194 CHARACTER DIST, PACKIT, TYPE, UPLO, XTYPE
195 CHARACTER*3 PATH
196 INTEGER I, I1, I2, IMAT, IN, INFO, IOFF, IRHS, IUPLO,
197 $ IZERO, J, K, KL, KU, LDA, MODE, N, NERRS,
198 $ NFAIL, NIMAT, NPP, NRHS, NRUN, NT
199 DOUBLE PRECISION ANORM, CNDNUM, RCOND, RCONDC
200* ..
201* .. Local Arrays ..
202 CHARACTER UPLOS( 2 )
203 INTEGER ISEED( 4 ), ISEEDY( 4 )
204 DOUBLE PRECISION RESULT( NTESTS )
205* ..
206* .. External Functions ..
207 LOGICAL LSAME
208 DOUBLE PRECISION DGET06, ZLANSP
209 EXTERNAL lsame, dget06, zlansp
210* ..
211* .. External Subroutines ..
212 EXTERNAL alaerh, alahd, alasum, zcopy, zerrsy, zget04,
215 $ zsptri, zsptrs
216* ..
217* .. Intrinsic Functions ..
218 INTRINSIC max, min
219* ..
220* .. Scalars in Common ..
221 LOGICAL LERR, OK
222 CHARACTER*32 SRNAMT
223 INTEGER INFOT, NUNIT
224* ..
225* .. Common blocks ..
226 COMMON / infoc / infot, nunit, ok, lerr
227 COMMON / srnamc / srnamt
228* ..
229* .. Data statements ..
230 DATA iseedy / 1988, 1989, 1990, 1991 /
231 DATA uplos / 'U', 'L' /
232* ..
233* .. Executable Statements ..
234*
235* Initialize constants and the random number seed.
236*
237 path( 1: 1 ) = 'Zomplex precision'
238 path( 2: 3 ) = 'SP'
239 nrun = 0
240 nfail = 0
241 nerrs = 0
242 DO 10 i = 1, 4
243 iseed( i ) = iseedy( i )
244 10 CONTINUE
245*
246* Test the error exits
247*
248 IF( tsterr )
249 $ CALL zerrsy( path, nout )
250 infot = 0
251*
252* Do for each value of N in NVAL
253*
254 DO 170 in = 1, nn
255 n = nval( in )
256 lda = max( n, 1 )
257 xtype = 'N'
258 nimat = ntypes
259 IF( n.LE.0 )
260 $ nimat = 1
261*
262 DO 160 imat = 1, nimat
263*
264* Do the tests only if DOTYPE( IMAT ) is true.
265*
266 IF( .NOT.dotype( imat ) )
267 $ GO TO 160
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 160
274*
275* Do first for UPLO = 'U', then for UPLO = 'L'
276*
277 DO 150 iuplo = 1, 2
278 uplo = uplos( iuplo )
279 IF( lsame( uplo, 'U' ) ) THEN
280 packit = 'C'
281 ELSE
282 packit = 'R'
283 END IF
284*
285 IF( imat.NE.ntypes ) THEN
286*
287* Set up parameters with ZLATB4 and generate a test
288* matrix with ZLATMS.
289*
290 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM,
291 $ MODE, CNDNUM, DIST )
292*
293 srnamt = 'ZLATMS'
294 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
295 $ CNDNUM, ANORM, KL, KU, PACKIT, A, LDA,
296 $ WORK, INFO )
297*
298* Check error code from ZLATMS.
299*
300 IF( info.NE.0 ) THEN
301 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
302 $ -1, -1, -1, imat, nfail, nerrs, nout )
303 GO TO 150
304 END IF
305*
306* For types 3-6, zero one or more rows and columns of
307* the matrix to test that INFO is returned correctly.
308*
309 IF( zerot ) THEN
310 IF( imat.EQ.3 ) THEN
311 izero = 1
312 ELSE IF( imat.EQ.4 ) THEN
313 izero = n
314 ELSE
315 izero = n / 2 + 1
316 END IF
317*
318 IF( imat.LT.6 ) THEN
319*
320* Set row and column IZERO to zero.
321*
322 IF( iuplo.EQ.1 ) THEN
323 ioff = ( izero-1 )*izero / 2
324 DO 20 i = 1, izero - 1
325 a( ioff+i ) = zero
326 20 CONTINUE
327 ioff = ioff + izero
328 DO 30 i = izero, n
329 a( ioff ) = zero
330 ioff = ioff + i
331 30 CONTINUE
332 ELSE
333 ioff = izero
334 DO 40 i = 1, izero - 1
335 a( ioff ) = zero
336 ioff = ioff + n - i
337 40 CONTINUE
338 ioff = ioff - izero
339 DO 50 i = izero, n
340 a( ioff+i ) = zero
341 50 CONTINUE
342 END IF
343 ELSE
344 IF( iuplo.EQ.1 ) THEN
345*
346* Set the first IZERO rows and columns to zero.
347*
348 ioff = 0
349 DO 70 j = 1, n
350 i2 = min( j, izero )
351 DO 60 i = 1, i2
352 a( ioff+i ) = zero
353 60 CONTINUE
354 ioff = ioff + j
355 70 CONTINUE
356 ELSE
357*
358* Set the last IZERO rows and columns to zero.
359*
360 ioff = 0
361 DO 90 j = 1, n
362 i1 = max( j, izero )
363 DO 80 i = i1, n
364 a( ioff+i ) = zero
365 80 CONTINUE
366 ioff = ioff + n - j
367 90 CONTINUE
368 END IF
369 END IF
370 ELSE
371 izero = 0
372 END IF
373 ELSE
374*
375* Use a special block diagonal matrix to test alternate
376* code for the 2 x 2 blocks.
377*
378 CALL zlatsp( uplo, n, a, iseed )
379 END IF
380*
381* Compute the L*D*L' or U*D*U' factorization of the matrix.
382*
383 npp = n*( n+1 ) / 2
384 CALL zcopy( npp, a, 1, afac, 1 )
385 srnamt = 'ZSPTRF'
386 CALL zsptrf( uplo, n, afac, iwork, info )
387*
388* Adjust the expected value of INFO to account for
389* pivoting.
390*
391 k = izero
392 IF( k.GT.0 ) THEN
393 100 CONTINUE
394 IF( iwork( k ).LT.0 ) THEN
395 IF( iwork( k ).NE.-k ) THEN
396 k = -iwork( k )
397 GO TO 100
398 END IF
399 ELSE IF( iwork( k ).NE.k ) THEN
400 k = iwork( k )
401 GO TO 100
402 END IF
403 END IF
404*
405* Check error code from ZSPTRF.
406*
407 IF( info.NE.k )
408 $ CALL alaerh( path, 'ZSPTRF', info, k, uplo, n, n, -1,
409 $ -1, -1, imat, nfail, nerrs, nout )
410 IF( info.NE.0 ) THEN
411 trfcon = .true.
412 ELSE
413 trfcon = .false.
414 END IF
415*
416*+ TEST 1
417* Reconstruct matrix from factors and compute residual.
418*
419 CALL zspt01( uplo, n, a, afac, iwork, ainv, lda, rwork,
420 $ result( 1 ) )
421 nt = 1
422*
423*+ TEST 2
424* Form the inverse and compute the residual.
425*
426 IF( .NOT.trfcon ) THEN
427 CALL zcopy( npp, afac, 1, ainv, 1 )
428 srnamt = 'ZSPTRI'
429 CALL zsptri( uplo, n, ainv, iwork, work, info )
430*
431* Check error code from ZSPTRI.
432*
433 IF( info.NE.0 )
434 $ CALL alaerh( path, 'ZSPTRI', info, 0, uplo, n, n,
435 $ -1, -1, -1, imat, nfail, nerrs, nout )
436*
437 CALL zspt03( uplo, n, a, ainv, work, lda, rwork,
438 $ rcondc, result( 2 ) )
439 nt = 2
440 END IF
441*
442* Print information about the tests that did not pass
443* the threshold.
444*
445 DO 110 k = 1, nt
446 IF( result( k ).GE.thresh ) THEN
447 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
448 $ CALL alahd( nout, path )
449 WRITE( nout, fmt = 9999 )uplo, n, imat, k,
450 $ result( k )
451 nfail = nfail + 1
452 END IF
453 110 CONTINUE
454 nrun = nrun + nt
455*
456* Do only the condition estimate if INFO is not 0.
457*
458 IF( trfcon ) THEN
459 rcondc = zero
460 GO TO 140
461 END IF
462*
463 DO 130 irhs = 1, nns
464 nrhs = nsval( irhs )
465*
466*+ TEST 3
467* Solve and compute residual for A * X = B.
468*
469 srnamt = 'ZLARHS'
470 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
471 $ nrhs, a, lda, xact, lda, b, lda, iseed,
472 $ info )
473 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
474*
475 srnamt = 'ZSPTRS'
476 CALL zsptrs( uplo, n, nrhs, afac, iwork, x, lda,
477 $ info )
478*
479* Check error code from ZSPTRS.
480*
481 IF( info.NE.0 )
482 $ CALL alaerh( path, 'ZSPTRS', info, 0, uplo, n, n,
483 $ -1, -1, nrhs, imat, nfail, nerrs,
484 $ nout )
485*
486 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
487 CALL zspt02( uplo, n, nrhs, a, x, lda, work, lda,
488 $ rwork, result( 3 ) )
489*
490*+ TEST 4
491* Check solution from generated exact solution.
492*
493 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
494 $ result( 4 ) )
495*
496*+ TESTS 5, 6, and 7
497* Use iterative refinement to improve the solution.
498*
499 srnamt = 'ZSPRFS'
500 CALL zsprfs( uplo, n, nrhs, a, afac, iwork, b, lda, x,
501 $ lda, rwork, rwork( nrhs+1 ), work,
502 $ rwork( 2*nrhs+1 ), info )
503*
504* Check error code from ZSPRFS.
505*
506 IF( info.NE.0 )
507 $ CALL alaerh( path, 'ZSPRFS', info, 0, uplo, n, n,
508 $ -1, -1, nrhs, imat, nfail, nerrs,
509 $ nout )
510*
511 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
512 $ result( 5 ) )
513 CALL zppt05( uplo, n, nrhs, a, b, lda, x, lda, xact,
514 $ lda, rwork, rwork( nrhs+1 ),
515 $ result( 6 ) )
516*
517* Print information about the tests that did not pass
518* the threshold.
519*
520 DO 120 k = 3, 7
521 IF( result( k ).GE.thresh ) THEN
522 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
523 $ CALL alahd( nout, path )
524 WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat,
525 $ k, result( k )
526 nfail = nfail + 1
527 END IF
528 120 CONTINUE
529 nrun = nrun + 5
530 130 CONTINUE
531*
532*+ TEST 8
533* Get an estimate of RCOND = 1/CNDNUM.
534*
535 140 CONTINUE
536 anorm = zlansp( '1', uplo, n, a, rwork )
537 srnamt = 'ZSPCON'
538 CALL zspcon( uplo, n, afac, iwork, anorm, rcond, work,
539 $ info )
540*
541* Check error code from ZSPCON.
542*
543 IF( info.NE.0 )
544 $ CALL alaerh( path, 'ZSPCON', info, 0, uplo, n, n, -1,
545 $ -1, -1, imat, nfail, nerrs, nout )
546*
547 result( 8 ) = dget06( rcond, rcondc )
548*
549* Print the test ratio if it is .GE. THRESH.
550*
551 IF( result( 8 ).GE.thresh ) THEN
552 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
553 $ CALL alahd( nout, path )
554 WRITE( nout, fmt = 9999 )uplo, n, imat, 8,
555 $ result( 8 )
556 nfail = nfail + 1
557 END IF
558 nrun = nrun + 1
559 150 CONTINUE
560 160 CONTINUE
561 170 CONTINUE
562*
563* Print a summary of the results.
564*
565 CALL alasum( path, nout, nfail, nrun, nerrs )
566*
567 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', type ', i2, ', test ',
568 $ i2, ', ratio =', g12.5 )
569 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
570 $ i2, ', test(', i2, ') =', g12.5 )
571 RETURN
572*
573* End of ZCHKSP
574*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
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 zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zspcon(uplo, n, ap, ipiv, anorm, rcond, work, info)
ZSPCON
Definition zspcon.f:118
subroutine zsprfs(uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZSPRFS
Definition zsprfs.f:180
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 zsptrs(uplo, n, nrhs, ap, ipiv, b, ldb, info)
ZSPTRS
Definition zsptrs.f:115
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 zlansp(norm, uplo, n, ap, work)
ZLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlansp.f:115
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
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 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
subroutine zspt03(uplo, n, a, ainv, work, ldw, rwork, rcond, resid)
ZSPT03
Definition zspt03.f:110
Here is the call graph for this function:
Here is the caller graph for this function: