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

◆ cchksy_aa()

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

CCHKSY_AA

Purpose:
 CCHKSY_AA tests CSYTRF_AA, -TRS_AA.
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 REAL
          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 REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]AINV
          AINV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is REAL array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is REAL array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is REAL array, dimension (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL array, dimension (max(NMAX,2*NSMAX))
[out]IWORK
          IWORK is INTEGER array, dimension (2*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 167 of file cchksy_aa.f.

170*
171* -- LAPACK test routine --
172* -- LAPACK is a software package provided by Univ. of Tennessee, --
173* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
174*
175 IMPLICIT NONE
176*
177* .. Scalar Arguments ..
178 LOGICAL TSTERR
179 INTEGER NN, NNB, NNS, NMAX, NOUT
180 REAL THRESH
181* ..
182* .. Array Arguments ..
183 LOGICAL DOTYPE( * )
184 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
185 REAL RWORK( * )
186 COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ),
187 $ WORK( * ), X( * ), XACT( * )
188* ..
189*
190* =====================================================================
191*
192* .. Parameters ..
193 REAL ZERO
194 parameter( zero = 0.0d+0 )
195 COMPLEX CZERO
196 parameter( czero = 0.0e+0 )
197 INTEGER NTYPES
198 parameter( ntypes = 10 )
199 INTEGER NTESTS
200 parameter( ntests = 9 )
201* ..
202* .. Local Scalars ..
203 LOGICAL ZEROT
204 CHARACTER DIST, TYPE, UPLO, XTYPE
205 CHARACTER*3 PATH, MATPATH
206 INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
207 $ IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
208 $ N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
209 REAL ANORM, CNDNUM
210* ..
211* .. Local Arrays ..
212 CHARACTER UPLOS( 2 )
213 INTEGER ISEED( 4 ), ISEEDY( 4 )
214 REAL RESULT( NTESTS )
215* ..
216* .. External Subroutines ..
217 EXTERNAL alaerh, alahd, alasum, cerrsy, clacpy, clarhs,
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC max, min
223* ..
224* .. Scalars in Common ..
225 LOGICAL LERR, OK
226 CHARACTER*32 SRNAMT
227 INTEGER INFOT, NUNIT
228* ..
229* .. Common blocks ..
230 COMMON / infoc / infot, nunit, ok, lerr
231 COMMON / srnamc / srnamt
232* ..
233* .. Data statements ..
234 DATA iseedy / 1988, 1989, 1990, 1991 /
235 DATA uplos / 'U', 'L' /
236* ..
237* .. Executable Statements ..
238*
239* Initialize constants and the random number seed.
240*
241* Test path
242*
243 path( 1: 1 ) = 'Complex precision'
244 path( 2: 3 ) = 'SA'
245*
246* Path to generate matrices
247*
248 matpath( 1: 1 ) = 'Complex precision'
249 matpath( 2: 3 ) = 'SY'
250 nrun = 0
251 nfail = 0
252 nerrs = 0
253 DO 10 i = 1, 4
254 iseed( i ) = iseedy( i )
255 10 CONTINUE
256*
257* Test the error exits
258*
259 IF( tsterr )
260 $ CALL cerrsy( path, nout )
261 infot = 0
262*
263* Set the minimum block size for which the block routine should
264* be used, which will be later returned by ILAENV
265*
266 CALL xlaenv( 2, 2 )
267*
268* Do for each value of N in NVAL
269*
270 DO 180 in = 1, nn
271 n = nval( in )
272 IF( n .GT. nmax ) THEN
273 nfail = nfail + 1
274 WRITE(nout, 9995) 'M ', n, nmax
275 GO TO 180
276 END IF
277 lda = max( n, 1 )
278 xtype = 'N'
279 nimat = ntypes
280 IF( n.LE.0 )
281 $ nimat = 1
282*
283 izero = 0
284*
285* Do for each value of matrix type IMAT
286*
287 DO 170 imat = 1, nimat
288*
289* Do the tests only if DOTYPE( IMAT ) is true.
290*
291 IF( .NOT.dotype( imat ) )
292 $ GO TO 170
293*
294* Skip types 3, 4, 5, or 6 if the matrix size is too small.
295*
296 zerot = imat.GE.3 .AND. imat.LE.6
297 IF( zerot .AND. n.LT.imat-2 )
298 $ GO TO 170
299*
300* Do first for UPLO = 'U', then for UPLO = 'L'
301*
302 DO 160 iuplo = 1, 2
303 uplo = uplos( iuplo )
304*
305* Begin generate the test matrix A.
306*
307*
308* Set up parameters with CLATB4 for the matrix generator
309* based on the type of matrix to be generated.
310*
311 CALL clatb4( matpath, imat, n, n, TYPE, KL, KU,
312 $ ANORM, MODE, CNDNUM, DIST )
313*
314* Generate a matrix with CLATMS.
315*
316 srnamt = 'CLATMS'
317 CALL clatms( n, n, dist, iseed, TYPE, RWORK, MODE,
318 $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
319 $ INFO )
320*
321* Check error code from CLATMS and handle error.
322*
323 IF( info.NE.0 ) THEN
324 CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n, -1,
325 $ -1, -1, imat, nfail, nerrs, nout )
326*
327* Skip all tests for this generated matrix
328*
329 GO TO 160
330 END IF
331*
332* For matrix types 3-6, zero one or more rows and
333* columns of the matrix to test that INFO is returned
334* correctly.
335*
336 IF( zerot ) THEN
337 IF( imat.EQ.3 ) THEN
338 izero = 1
339 ELSE IF( imat.EQ.4 ) THEN
340 izero = n
341 ELSE
342 izero = n / 2 + 1
343 END IF
344*
345 IF( imat.LT.6 ) THEN
346*
347* Set row and column IZERO to zero.
348*
349 IF( iuplo.EQ.1 ) THEN
350 ioff = ( izero-1 )*lda
351 DO 20 i = 1, izero - 1
352 a( ioff+i ) = czero
353 20 CONTINUE
354 ioff = ioff + izero
355 DO 30 i = izero, n
356 a( ioff ) = czero
357 ioff = ioff + lda
358 30 CONTINUE
359 ELSE
360 ioff = izero
361 DO 40 i = 1, izero - 1
362 a( ioff ) = czero
363 ioff = ioff + lda
364 40 CONTINUE
365 ioff = ioff - izero
366 DO 50 i = izero, n
367 a( ioff+i ) = czero
368 50 CONTINUE
369 END IF
370 ELSE
371 IF( iuplo.EQ.1 ) THEN
372*
373* Set the first IZERO rows and columns to zero.
374*
375 ioff = 0
376 DO 70 j = 1, n
377 i2 = min( j, izero )
378 DO 60 i = 1, i2
379 a( ioff+i ) = czero
380 60 CONTINUE
381 ioff = ioff + lda
382 70 CONTINUE
383 izero = 1
384 ELSE
385*
386* Set the last IZERO rows and columns to zero.
387*
388 ioff = 0
389 DO 90 j = 1, n
390 i1 = max( j, izero )
391 DO 80 i = i1, n
392 a( ioff+i ) = czero
393 80 CONTINUE
394 ioff = ioff + lda
395 90 CONTINUE
396 END IF
397 END IF
398 ELSE
399 izero = 0
400 END IF
401*
402* End generate the test matrix A.
403*
404* Do for each value of NB in NBVAL
405*
406 DO 150 inb = 1, nnb
407*
408* Set the optimal blocksize, which will be later
409* returned by ILAENV.
410*
411 nb = nbval( inb )
412 CALL xlaenv( 1, nb )
413*
414* Copy the test matrix A into matrix AFAC which
415* will be factorized in place. This is needed to
416* preserve the test matrix A for subsequent tests.
417*
418 CALL clacpy( uplo, n, n, a, lda, afac, lda )
419*
420* Compute the L*D*L**T or U*D*U**T factorization of the
421* matrix. IWORK stores details of the interchanges and
422* the block structure of D. AINV is a work array for
423* block factorization, LWORK is the length of AINV.
424*
425 srnamt = 'CSYTRF_AA'
426 lwork = max( 1, n*nb + n )
427 CALL csytrf_aa( uplo, n, afac, lda, iwork, ainv,
428 $ lwork, info )
429*
430* Adjust the expected value of INFO to account for
431* pivoting.
432*
433c IF( IZERO.GT.0 ) THEN
434c J = 1
435c K = IZERO
436c 100 CONTINUE
437c IF( J.EQ.K ) THEN
438c K = IWORK( J )
439c ELSE IF( IWORK( J ).EQ.K ) THEN
440c K = J
441c END IF
442c IF( J.LT.K ) THEN
443c J = J + 1
444c GO TO 100
445c END IF
446c ELSE
447 k = 0
448c END IF
449*
450* Check error code from CSYTRF and handle error.
451*
452 IF( info.NE.k ) THEN
453 CALL alaerh( path, 'CSYTRF_AA', info, k, uplo,
454 $ n, n, -1, -1, nb, imat, nfail, nerrs,
455 $ nout )
456 END IF
457*
458*+ TEST 1
459* Reconstruct matrix from factors and compute residual.
460*
461 CALL csyt01_aa( uplo, n, a, lda, afac, lda, iwork,
462 $ ainv, lda, rwork, result( 1 ) )
463 nt = 1
464*
465*
466* Print information about the tests that did not pass
467* the threshold.
468*
469 DO 110 k = 1, nt
470 IF( result( k ).GE.thresh ) THEN
471 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
472 $ CALL alahd( nout, path )
473 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
474 $ result( k )
475 nfail = nfail + 1
476 END IF
477 110 CONTINUE
478 nrun = nrun + nt
479*
480* Skip solver test if INFO is not 0.
481*
482 IF( info.NE.0 ) THEN
483 GO TO 140
484 END IF
485*
486* Do for each value of NRHS in NSVAL.
487*
488 DO 130 irhs = 1, nns
489 nrhs = nsval( irhs )
490*
491*+ TEST 2 (Using TRS)
492* Solve and compute residual for A * X = B.
493*
494* Choose a set of NRHS random solution vectors
495* stored in XACT and set up the right hand side B
496*
497 srnamt = 'CLARHS'
498 CALL clarhs( matpath, xtype, uplo, ' ', n, n,
499 $ kl, ku, nrhs, a, lda, xact, lda,
500 $ b, lda, iseed, info )
501 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
502*
503 srnamt = 'CSYTRS_AA'
504 lwork = max( 1, 3*n-2 )
505 CALL csytrs_aa( uplo, n, nrhs, afac, lda,
506 $ iwork, x, lda, work, lwork,
507 $ info )
508*
509* Check error code from CSYTRS and handle error.
510*
511 IF( info.NE.0 ) THEN
512 IF( izero.EQ.0 ) THEN
513 CALL alaerh( path, 'CSYTRS_AA', info, 0,
514 $ uplo, n, n, -1, -1, nrhs, imat,
515 $ nfail, nerrs, nout )
516 END IF
517 ELSE
518 CALL clacpy( 'Full', n, nrhs, b, lda, work, lda
519 $ )
520*
521* Compute the residual for the solution
522*
523 CALL csyt02( uplo, n, nrhs, a, lda, x, lda,
524 $ work, lda, rwork, result( 2 ) )
525*
526*
527* Print information about the tests that did not pass
528* the threshold.
529*
530 DO 120 k = 2, 2
531 IF( result( k ).GE.thresh ) THEN
532 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
533 $ CALL alahd( nout, path )
534 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
535 $ imat, k, result( k )
536 nfail = nfail + 1
537 END IF
538 120 CONTINUE
539 END IF
540 nrun = nrun + 1
541*
542* End do for each value of NRHS in NSVAL.
543*
544 130 CONTINUE
545 140 CONTINUE
546 150 CONTINUE
547 160 CONTINUE
548 170 CONTINUE
549 180 CONTINUE
550*
551* Print a summary of the results.
552*
553 CALL alasum( path, nout, nfail, nrun, nerrs )
554*
555 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
556 $ i2, ', test ', i2, ', ratio =', g12.5 )
557 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
558 $ i2, ', test(', i2, ') =', g12.5 )
559 9995 FORMAT( ' Invalid input value: ', a4, '=', i6, '; must be <=',
560 $ i6 )
561 RETURN
562*
563* End of CCHKSY_AA
564*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine clarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
CLARHS
Definition clarhs.f:208
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
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
subroutine cerrsy(path, nunit)
CERRSY
Definition cerrsy.f:55
subroutine clatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
CLATB4
Definition clatb4.f:121
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine csyt01_aa(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
CSYT01
Definition csyt01_aa.f:124
subroutine csyt02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
CSYT02
Definition csyt02.f:127
subroutine csytrf_aa(uplo, n, a, lda, ipiv, work, lwork, info)
CSYTRF_AA
Definition csytrf_aa.f:132
subroutine csytrs_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
CSYTRS_AA
Definition csytrs_aa.f:131
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
Here is the call graph for this function:
Here is the caller graph for this function: