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

◆ zchkpp()

subroutine zchkpp ( 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  nout 
)

ZCHKPP

Purpose:
 ZCHKPP tests ZPPTRF, -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(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(NMAX,2*NSMAX))
[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 156 of file zchkpp.f.

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, NNS, NOUT
167 DOUBLE PRECISION THRESH
168* ..
169* .. Array Arguments ..
170 LOGICAL DOTYPE( * )
171 INTEGER NSVAL( * ), NVAL( * )
172 DOUBLE PRECISION RWORK( * )
173 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
174 $ WORK( * ), X( * ), XACT( * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 DOUBLE PRECISION ZERO
181 parameter( zero = 0.0d+0 )
182 INTEGER NTYPES
183 parameter( ntypes = 9 )
184 INTEGER NTESTS
185 parameter( ntests = 8 )
186* ..
187* .. Local Scalars ..
188 LOGICAL ZEROT
189 CHARACTER DIST, PACKIT, TYPE, UPLO, XTYPE
190 CHARACTER*3 PATH
191 INTEGER I, IMAT, IN, INFO, IOFF, IRHS, IUPLO, IZERO, K,
192 $ KL, KU, LDA, MODE, N, NERRS, NFAIL, NIMAT, NPP,
193 $ NRHS, NRUN
194 DOUBLE PRECISION ANORM, CNDNUM, RCOND, RCONDC
195* ..
196* .. Local Arrays ..
197 CHARACTER PACKS( 2 ), UPLOS( 2 )
198 INTEGER ISEED( 4 ), ISEEDY( 4 )
199 DOUBLE PRECISION RESULT( NTESTS )
200* ..
201* .. External Functions ..
202 DOUBLE PRECISION DGET06, ZLANHP
203 EXTERNAL dget06, zlanhp
204* ..
205* .. External Subroutines ..
206 EXTERNAL alaerh, alahd, alasum, zcopy, zerrpo, zget04,
209 $ zpptri, zpptrs
210* ..
211* .. Scalars in Common ..
212 LOGICAL LERR, OK
213 CHARACTER*32 SRNAMT
214 INTEGER INFOT, NUNIT
215* ..
216* .. Common blocks ..
217 COMMON / infoc / infot, nunit, ok, lerr
218 COMMON / srnamc / srnamt
219* ..
220* .. Intrinsic Functions ..
221 INTRINSIC max
222* ..
223* .. Data statements ..
224 DATA iseedy / 1988, 1989, 1990, 1991 /
225 DATA uplos / 'U', 'L' / , packs / 'C', 'R' /
226* ..
227* .. Executable Statements ..
228*
229* Initialize constants and the random number seed.
230*
231 path( 1: 1 ) = 'Zomplex precision'
232 path( 2: 3 ) = 'PP'
233 nrun = 0
234 nfail = 0
235 nerrs = 0
236 DO 10 i = 1, 4
237 iseed( i ) = iseedy( i )
238 10 CONTINUE
239*
240* Test the error exits
241*
242 IF( tsterr )
243 $ CALL zerrpo( path, nout )
244 infot = 0
245*
246* Do for each value of N in NVAL
247*
248 DO 110 in = 1, nn
249 n = nval( in )
250 lda = max( n, 1 )
251 xtype = 'N'
252 nimat = ntypes
253 IF( n.LE.0 )
254 $ nimat = 1
255*
256 DO 100 imat = 1, nimat
257*
258* Do the tests only if DOTYPE( IMAT ) is true.
259*
260 IF( .NOT.dotype( imat ) )
261 $ GO TO 100
262*
263* Skip types 3, 4, or 5 if the matrix size is too small.
264*
265 zerot = imat.GE.3 .AND. imat.LE.5
266 IF( zerot .AND. n.LT.imat-2 )
267 $ GO TO 100
268*
269* Do first for UPLO = 'U', then for UPLO = 'L'
270*
271 DO 90 iuplo = 1, 2
272 uplo = uplos( iuplo )
273 packit = packs( iuplo )
274*
275* Set up parameters with ZLATB4 and generate a test matrix
276* with ZLATMS.
277*
278 CALL zlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
279 $ CNDNUM, DIST )
280*
281 srnamt = 'ZLATMS'
282 CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
283 $ CNDNUM, ANORM, KL, KU, PACKIT, A, LDA, WORK,
284 $ INFO )
285*
286* Check error code from ZLATMS.
287*
288 IF( info.NE.0 ) THEN
289 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
290 $ -1, -1, imat, nfail, nerrs, nout )
291 GO TO 90
292 END IF
293*
294* For types 3-5, zero one row and column of the matrix to
295* test that INFO is returned correctly.
296*
297 IF( zerot ) THEN
298 IF( imat.EQ.3 ) THEN
299 izero = 1
300 ELSE IF( imat.EQ.4 ) THEN
301 izero = n
302 ELSE
303 izero = n / 2 + 1
304 END IF
305*
306* Set row and column IZERO of A to 0.
307*
308 IF( iuplo.EQ.1 ) THEN
309 ioff = ( izero-1 )*izero / 2
310 DO 20 i = 1, izero - 1
311 a( ioff+i ) = zero
312 20 CONTINUE
313 ioff = ioff + izero
314 DO 30 i = izero, n
315 a( ioff ) = zero
316 ioff = ioff + i
317 30 CONTINUE
318 ELSE
319 ioff = izero
320 DO 40 i = 1, izero - 1
321 a( ioff ) = zero
322 ioff = ioff + n - i
323 40 CONTINUE
324 ioff = ioff - izero
325 DO 50 i = izero, n
326 a( ioff+i ) = zero
327 50 CONTINUE
328 END IF
329 ELSE
330 izero = 0
331 END IF
332*
333* Set the imaginary part of the diagonals.
334*
335 IF( iuplo.EQ.1 ) THEN
336 CALL zlaipd( n, a, 2, 1 )
337 ELSE
338 CALL zlaipd( n, a, n, -1 )
339 END IF
340*
341* Compute the L*L' or U'*U factorization of the matrix.
342*
343 npp = n*( n+1 ) / 2
344 CALL zcopy( npp, a, 1, afac, 1 )
345 srnamt = 'ZPPTRF'
346 CALL zpptrf( uplo, n, afac, info )
347*
348* Check error code from ZPPTRF.
349*
350 IF( info.NE.izero ) THEN
351 CALL alaerh( path, 'ZPPTRF', info, izero, uplo, n, n,
352 $ -1, -1, -1, imat, nfail, nerrs, nout )
353 GO TO 90
354 END IF
355*
356* Skip the tests if INFO is not 0.
357*
358 IF( info.NE.0 )
359 $ GO TO 90
360*
361*+ TEST 1
362* Reconstruct matrix from factors and compute residual.
363*
364 CALL zcopy( npp, afac, 1, ainv, 1 )
365 CALL zppt01( uplo, n, a, ainv, rwork, result( 1 ) )
366*
367*+ TEST 2
368* Form the inverse and compute the residual.
369*
370 CALL zcopy( npp, afac, 1, ainv, 1 )
371 srnamt = 'ZPPTRI'
372 CALL zpptri( uplo, n, ainv, info )
373*
374* Check error code from ZPPTRI.
375*
376 IF( info.NE.0 )
377 $ CALL alaerh( path, 'ZPPTRI', info, 0, uplo, n, n, -1,
378 $ -1, -1, imat, nfail, nerrs, nout )
379*
380 CALL zppt03( uplo, n, a, ainv, work, lda, rwork, rcondc,
381 $ result( 2 ) )
382*
383* Print information about the tests that did not pass
384* the threshold.
385*
386 DO 60 k = 1, 2
387 IF( result( k ).GE.thresh ) THEN
388 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
389 $ CALL alahd( nout, path )
390 WRITE( nout, fmt = 9999 )uplo, n, imat, k,
391 $ result( k )
392 nfail = nfail + 1
393 END IF
394 60 CONTINUE
395 nrun = nrun + 2
396*
397 DO 80 irhs = 1, nns
398 nrhs = nsval( irhs )
399*
400*+ TEST 3
401* Solve and compute residual for A * X = B.
402*
403 srnamt = 'ZLARHS'
404 CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
405 $ nrhs, a, lda, xact, lda, b, lda, iseed,
406 $ info )
407 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
408*
409 srnamt = 'ZPPTRS'
410 CALL zpptrs( uplo, n, nrhs, afac, x, lda, info )
411*
412* Check error code from ZPPTRS.
413*
414 IF( info.NE.0 )
415 $ CALL alaerh( path, 'ZPPTRS', info, 0, uplo, n, n,
416 $ -1, -1, nrhs, imat, nfail, nerrs,
417 $ nout )
418*
419 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
420 CALL zppt02( uplo, n, nrhs, a, x, lda, work, lda,
421 $ rwork, result( 3 ) )
422*
423*+ TEST 4
424* Check solution from generated exact solution.
425*
426 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
427 $ result( 4 ) )
428*
429*+ TESTS 5, 6, and 7
430* Use iterative refinement to improve the solution.
431*
432 srnamt = 'ZPPRFS'
433 CALL zpprfs( uplo, n, nrhs, a, afac, b, lda, x, lda,
434 $ rwork, rwork( nrhs+1 ), work,
435 $ rwork( 2*nrhs+1 ), info )
436*
437* Check error code from ZPPRFS.
438*
439 IF( info.NE.0 )
440 $ CALL alaerh( path, 'ZPPRFS', info, 0, uplo, n, n,
441 $ -1, -1, nrhs, imat, nfail, nerrs,
442 $ nout )
443*
444 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
445 $ result( 5 ) )
446 CALL zppt05( uplo, n, nrhs, a, b, lda, x, lda, xact,
447 $ lda, rwork, rwork( nrhs+1 ),
448 $ result( 6 ) )
449*
450* Print information about the tests that did not pass
451* the threshold.
452*
453 DO 70 k = 3, 7
454 IF( result( k ).GE.thresh ) THEN
455 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
456 $ CALL alahd( nout, path )
457 WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat,
458 $ k, result( k )
459 nfail = nfail + 1
460 END IF
461 70 CONTINUE
462 nrun = nrun + 5
463 80 CONTINUE
464*
465*+ TEST 8
466* Get an estimate of RCOND = 1/CNDNUM.
467*
468 anorm = zlanhp( '1', uplo, n, a, rwork )
469 srnamt = 'ZPPCON'
470 CALL zppcon( uplo, n, afac, anorm, rcond, work, rwork,
471 $ info )
472*
473* Check error code from ZPPCON.
474*
475 IF( info.NE.0 )
476 $ CALL alaerh( path, 'ZPPCON', info, 0, uplo, n, n, -1,
477 $ -1, -1, imat, nfail, nerrs, nout )
478*
479 result( 8 ) = dget06( rcond, rcondc )
480*
481* Print the test ratio if greater than or equal to THRESH.
482*
483 IF( result( 8 ).GE.thresh ) THEN
484 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
485 $ CALL alahd( nout, path )
486 WRITE( nout, fmt = 9999 )uplo, n, imat, 8,
487 $ result( 8 )
488 nfail = nfail + 1
489 END IF
490 nrun = nrun + 1
491*
492 90 CONTINUE
493 100 CONTINUE
494 110 CONTINUE
495*
496* Print a summary of the results.
497*
498 CALL alasum( path, nout, nfail, nrun, nerrs )
499*
500 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', type ', i2, ', test ',
501 $ i2, ', ratio =', g12.5 )
502 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
503 $ i2, ', test(', i2, ') =', g12.5 )
504 RETURN
505*
506* End of ZCHKPP
507*
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 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 zlanhp(norm, uplo, n, ap, work)
ZLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhp.f:117
subroutine zppcon(uplo, n, ap, anorm, rcond, work, rwork, info)
ZPPCON
Definition zppcon.f:118
subroutine zpprfs(uplo, n, nrhs, ap, afp, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPPRFS
Definition zpprfs.f:171
subroutine zpptrf(uplo, n, ap, info)
ZPPTRF
Definition zpptrf.f:119
subroutine zpptri(uplo, n, ap, info)
ZPPTRI
Definition zpptri.f:93
subroutine zpptrs(uplo, n, nrhs, ap, b, ldb, info)
ZPPTRS
Definition zpptrs.f:108
subroutine zerrpo(path, nunit)
ZERRPO
Definition zerrpo.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlaipd(n, a, inda, vinda)
ZLAIPD
Definition zlaipd.f:83
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 zppt01(uplo, n, a, afac, rwork, resid)
ZPPT01
Definition zppt01.f:95
subroutine zppt02(uplo, n, nrhs, a, x, ldx, b, ldb, rwork, resid)
ZPPT02
Definition zppt02.f:123
subroutine zppt03(uplo, n, a, ainv, work, ldwork, rwork, rcond, resid)
ZPPT03
Definition zppt03.f:110
subroutine zppt05(uplo, n, nrhs, ap, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZPPT05
Definition zppt05.f:157
Here is the call graph for this function:
Here is the caller graph for this function: