LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dchksy_aa.f
Go to the documentation of this file.
1*> \brief \b DCHKSY_AA
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 DCHKSY_AA( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
12* THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
13* XACT, WORK, RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NNB, NNS, NOUT
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
23* DOUBLE PRECISION A( * ), AFAC( * ), AINV( * ), B( * ),
24* $ RWORK( * ), WORK( * ), X( * ), XACT( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> DCHKSY_AA tests DSYTRF_AA, -TRS_AA.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] DOTYPE
40*> \verbatim
41*> DOTYPE is LOGICAL array, dimension (NTYPES)
42*> The matrix types to be used for testing. Matrices of type j
43*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
44*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
45*> \endverbatim
46*>
47*> \param[in] NN
48*> \verbatim
49*> NN is INTEGER
50*> The number of values of N contained in the vector NVAL.
51*> \endverbatim
52*>
53*> \param[in] NVAL
54*> \verbatim
55*> NVAL is INTEGER array, dimension (NN)
56*> The values of the matrix dimension N.
57*> \endverbatim
58*>
59*> \param[in] NNB
60*> \verbatim
61*> NNB is INTEGER
62*> The number of values of NB contained in the vector NBVAL.
63*> \endverbatim
64*>
65*> \param[in] NBVAL
66*> \verbatim
67*> NBVAL is INTEGER array, dimension (NNB)
68*> The values of the blocksize NB.
69*> \endverbatim
70*>
71*> \param[in] NNS
72*> \verbatim
73*> NNS is INTEGER
74*> The number of values of NRHS contained in the vector NSVAL.
75*> \endverbatim
76*>
77*> \param[in] NSVAL
78*> \verbatim
79*> NSVAL is INTEGER array, dimension (NNS)
80*> The values of the number of right hand sides NRHS.
81*> \endverbatim
82*>
83*> \param[in] THRESH
84*> \verbatim
85*> THRESH is DOUBLE PRECISION
86*> The threshold value for the test ratios. A result is
87*> included in the output file if RESULT >= THRESH. To have
88*> every test ratio printed, use THRESH = 0.
89*> \endverbatim
90*>
91*> \param[in] TSTERR
92*> \verbatim
93*> TSTERR is LOGICAL
94*> Flag that indicates whether error exits are to be tested.
95*> \endverbatim
96*>
97*> \param[in] NMAX
98*> \verbatim
99*> NMAX is INTEGER
100*> The maximum value permitted for N, used in dimensioning the
101*> work arrays.
102*> \endverbatim
103*>
104*> \param[out] A
105*> \verbatim
106*> A is DOUBLE PRECISION array, dimension (NMAX*NMAX)
107*> \endverbatim
108*>
109*> \param[out] AFAC
110*> \verbatim
111*> AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
112*> \endverbatim
113*>
114*> \param[out] AINV
115*> \verbatim
116*> AINV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
117*> \endverbatim
118*>
119*> \param[out] B
120*> \verbatim
121*> B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
122*> where NSMAX is the largest entry in NSVAL.
123*> \endverbatim
124*>
125*> \param[out] X
126*> \verbatim
127*> X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
128*> \endverbatim
129*>
130*> \param[out] XACT
131*> \verbatim
132*> XACT is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
133*> \endverbatim
134*>
135*> \param[out] WORK
136*> \verbatim
137*> WORK is DOUBLE PRECISION array, dimension (NMAX*max(3,NSMAX))
138*> \endverbatim
139*>
140*> \param[out] RWORK
141*> \verbatim
142*> RWORK is DOUBLE PRECISION array, dimension (max(NMAX,2*NSMAX))
143*> \endverbatim
144*>
145*> \param[out] IWORK
146*> \verbatim
147*> IWORK is INTEGER array, dimension (2*NMAX)
148*> \endverbatim
149*>
150*> \param[in] NOUT
151*> \verbatim
152*> NOUT is INTEGER
153*> The unit number for output.
154*> \endverbatim
155*
156* Authors:
157* ========
158*
159*> \author Univ. of Tennessee
160*> \author Univ. of California Berkeley
161*> \author Univ. of Colorado Denver
162*> \author NAG Ltd.
163*
164*> \ingroup double_lin
165*
166* =====================================================================
167 SUBROUTINE dchksy_aa( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
168 $ THRESH, TSTERR, NMAX, A, AFAC, AINV, B,
169 $ X, XACT, WORK, RWORK, IWORK, NOUT )
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 DOUBLE PRECISION THRESH
181* ..
182* .. Array Arguments ..
183 LOGICAL DOTYPE( * )
184 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
185 DOUBLE PRECISION A( * ), AFAC( * ), AINV( * ), B( * ),
186 $ rwork( * ), work( * ), x( * ), xact( * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 DOUBLE PRECISION ZERO, ONE
193 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
194 INTEGER NTYPES
195 parameter( ntypes = 10 )
196 INTEGER NTESTS
197 parameter( ntests = 9 )
198* ..
199* .. Local Scalars ..
200 LOGICAL ZEROT
201 CHARACTER DIST, TYPE, UPLO, XTYPE
202 CHARACTER*3 PATH, MATPATH
203 INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
204 $ iuplo, izero, j, k, kl, ku, lda, lwork, mode,
205 $ n, nb, nerrs, nfail, nimat, nrhs, nrun, nt
206 DOUBLE PRECISION ANORM, CNDNUM
207* ..
208* .. Local Arrays ..
209 CHARACTER UPLOS( 2 )
210 INTEGER ISEED( 4 ), ISEEDY( 4 )
211 DOUBLE PRECISION RESULT( NTESTS )
212* ..
213* .. External Subroutines ..
214 EXTERNAL alaerh, alahd, alasum, derrsy, dlacpy, dlarhs,
217* ..
218* .. Intrinsic Functions ..
219 INTRINSIC max, min
220* ..
221* .. Scalars in Common ..
222 LOGICAL LERR, OK
223 CHARACTER*32 SRNAMT
224 INTEGER INFOT, NUNIT
225* ..
226* .. Common blocks ..
227 COMMON / infoc / infot, nunit, ok, lerr
228 COMMON / srnamc / srnamt
229* ..
230* .. Data statements ..
231 DATA iseedy / 1988, 1989, 1990, 1991 /
232 DATA uplos / 'U', 'L' /
233* ..
234* .. Executable Statements ..
235*
236* Initialize constants and the random number seed.
237*
238* Test path
239*
240 path( 1: 1 ) = 'Double precision'
241 path( 2: 3 ) = 'SA'
242*
243* Path to generate matrices
244*
245 matpath( 1: 1 ) = 'Double precision'
246 matpath( 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 derrsy( 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 IF( n .GT. nmax ) THEN
270 nfail = nfail + 1
271 WRITE(nout, 9995) 'M ', n, nmax
272 GO TO 180
273 END IF
274 lda = max( n, 1 )
275 xtype = 'N'
276 nimat = ntypes
277 IF( n.LE.0 )
278 $ nimat = 1
279*
280 izero = 0
281*
282* Do for each value of matrix type IMAT
283*
284 DO 170 imat = 1, nimat
285*
286* Do the tests only if DOTYPE( IMAT ) is true.
287*
288 IF( .NOT.dotype( imat ) )
289 $ GO TO 170
290*
291* Skip types 3, 4, 5, or 6 if the matrix size is too small.
292*
293 zerot = imat.GE.3 .AND. imat.LE.6
294 IF( zerot .AND. n.LT.imat-2 )
295 $ GO TO 170
296*
297* Do first for UPLO = 'U', then for UPLO = 'L'
298*
299 DO 160 iuplo = 1, 2
300 uplo = uplos( iuplo )
301*
302* Begin generate the test matrix A.
303*
304*
305* Set up parameters with DLATB4 for the matrix generator
306* based on the type of matrix to be generated.
307*
308 CALL dlatb4( matpath, imat, n, n, TYPE, kl, ku,
309 $ anorm, mode, cndnum, dist )
310*
311* Generate a matrix with DLATMS.
312*
313 srnamt = 'DLATMS'
314 CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
315 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
316 $ info )
317*
318* Check error code from DLATMS and handle error.
319*
320 IF( info.NE.0 ) THEN
321 CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
322 $ -1, -1, imat, nfail, nerrs, nout )
323*
324* Skip all tests for this generated matrix
325*
326 GO TO 160
327 END IF
328*
329* For matrix types 3-6, zero one or more rows and
330* columns of the matrix to test that INFO is returned
331* correctly.
332*
333 IF( zerot ) THEN
334 IF( imat.EQ.3 ) THEN
335 izero = 1
336 ELSE IF( imat.EQ.4 ) THEN
337 izero = n
338 ELSE
339 izero = n / 2 + 1
340 END IF
341*
342 IF( imat.LT.6 ) THEN
343*
344* Set row and column IZERO to zero.
345*
346 IF( iuplo.EQ.1 ) THEN
347 ioff = ( izero-1 )*lda
348 DO 20 i = 1, izero - 1
349 a( ioff+i ) = zero
350 20 CONTINUE
351 ioff = ioff + izero
352 DO 30 i = izero, n
353 a( ioff ) = zero
354 ioff = ioff + lda
355 30 CONTINUE
356 ELSE
357 ioff = izero
358 DO 40 i = 1, izero - 1
359 a( ioff ) = zero
360 ioff = ioff + lda
361 40 CONTINUE
362 ioff = ioff - izero
363 DO 50 i = izero, n
364 a( ioff+i ) = zero
365 50 CONTINUE
366 END IF
367 ELSE
368 IF( iuplo.EQ.1 ) THEN
369*
370* Set the first IZERO rows and columns to zero.
371*
372 ioff = 0
373 DO 70 j = 1, n
374 i2 = min( j, izero )
375 DO 60 i = 1, i2
376 a( ioff+i ) = zero
377 60 CONTINUE
378 ioff = ioff + lda
379 70 CONTINUE
380 izero = 1
381 ELSE
382*
383* Set the last IZERO rows and columns to zero.
384*
385 ioff = 0
386 DO 90 j = 1, n
387 i1 = max( j, izero )
388 DO 80 i = i1, n
389 a( ioff+i ) = zero
390 80 CONTINUE
391 ioff = ioff + lda
392 90 CONTINUE
393 END IF
394 END IF
395 ELSE
396 izero = 0
397 END IF
398*
399* End generate the test matrix A.
400*
401* Do for each value of NB in NBVAL
402*
403 DO 150 inb = 1, nnb
404*
405* Set the optimal blocksize, which will be later
406* returned by ILAENV.
407*
408 nb = nbval( inb )
409 CALL xlaenv( 1, nb )
410*
411* Copy the test matrix A into matrix AFAC which
412* will be factorized in place. This is needed to
413* preserve the test matrix A for subsequent tests.
414*
415 CALL dlacpy( uplo, n, n, a, lda, afac, lda )
416*
417* Compute the L*D*L**T or U*D*U**T factorization of the
418* matrix. IWORK stores details of the interchanges and
419* the block structure of D. AINV is a work array for
420* block factorization, LWORK is the length of AINV.
421*
422 srnamt = 'DSYTRF_AA'
423 lwork = max( 1, n*nb + n )
424 CALL dsytrf_aa( uplo, n, afac, lda, iwork, ainv,
425 $ lwork, info )
426*
427* Adjust the expected value of INFO to account for
428* pivoting.
429*
430c IF( IZERO.GT.0 ) THEN
431c J = 1
432c K = IZERO
433c 100 CONTINUE
434c IF( J.EQ.K ) THEN
435c K = IWORK( J )
436c ELSE IF( IWORK( J ).EQ.K ) THEN
437c K = J
438c END IF
439c IF( J.LT.K ) THEN
440c J = J + 1
441c GO TO 100
442c END IF
443c ELSE
444 k = 0
445c END IF
446*
447* Check error code from DSYTRF and handle error.
448*
449 IF( info.NE.k ) THEN
450 CALL alaerh( path, 'DSYTRF_AA', info, k, uplo,
451 $ n, n, -1, -1, nb, imat, nfail, nerrs,
452 $ nout )
453 END IF
454*
455*+ TEST 1
456* Reconstruct matrix from factors and compute residual.
457*
458 CALL dsyt01_aa( uplo, n, a, lda, afac, lda, iwork,
459 $ ainv, lda, rwork, result( 1 ) )
460 nt = 1
461*
462*
463* Print information about the tests that did not pass
464* the threshold.
465*
466 DO 110 k = 1, nt
467 IF( result( k ).GE.thresh ) THEN
468 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
469 $ CALL alahd( nout, path )
470 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
471 $ result( k )
472 nfail = nfail + 1
473 END IF
474 110 CONTINUE
475 nrun = nrun + nt
476*
477* Skip solver test if INFO is not 0.
478*
479 IF( info.NE.0 ) THEN
480 GO TO 140
481 END IF
482*
483* Do for each value of NRHS in NSVAL.
484*
485 DO 130 irhs = 1, nns
486 nrhs = nsval( irhs )
487*
488*+ TEST 2 (Using TRS)
489* Solve and compute residual for A * X = B.
490*
491* Choose a set of NRHS random solution vectors
492* stored in XACT and set up the right hand side B
493*
494 srnamt = 'DLARHS'
495 CALL dlarhs( matpath, xtype, uplo, ' ', n, n,
496 $ kl, ku, nrhs, a, lda, xact, lda,
497 $ b, lda, iseed, info )
498 CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
499*
500 srnamt = 'DSYTRS_AA'
501 lwork = max( 1, 3*n-2 )
502 CALL dsytrs_aa( uplo, n, nrhs, afac, lda,
503 $ iwork, x, lda, work, lwork,
504 $ info )
505*
506* Check error code from DSYTRS and handle error.
507*
508 IF( info.NE.0 ) THEN
509 IF( izero.EQ.0 ) THEN
510 CALL alaerh( path, 'DSYTRS_AA', info, 0,
511 $ uplo, n, n, -1, -1, nrhs, imat,
512 $ nfail, nerrs, nout )
513 END IF
514 ELSE
515 CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda
516 $ )
517*
518* Compute the residual for the solution
519*
520 CALL dpot02( uplo, n, nrhs, a, lda, x, lda,
521 $ work, lda, rwork, result( 2 ) )
522*
523*
524* Print information about the tests that did not pass
525* the threshold.
526*
527 DO 120 k = 2, 2
528 IF( result( k ).GE.thresh ) THEN
529 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
530 $ CALL alahd( nout, path )
531 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
532 $ imat, k, result( k )
533 nfail = nfail + 1
534 END IF
535 120 CONTINUE
536 END IF
537 nrun = nrun + 1
538*
539* End do for each value of NRHS in NSVAL.
540*
541 130 CONTINUE
542 140 CONTINUE
543 150 CONTINUE
544 160 CONTINUE
545 170 CONTINUE
546 180 CONTINUE
547*
548* Print a summary of the results.
549*
550 CALL alasum( path, nout, nfail, nrun, nerrs )
551*
552 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
553 $ i2, ', test ', i2, ', ratio =', g12.5 )
554 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
555 $ i2, ', test(', i2, ') =', g12.5 )
556 9995 FORMAT( ' Invalid input value: ', a4, '=', i6, '; must be <=',
557 $ i6 )
558 RETURN
559*
560* End of DCHKSY_AA
561*
562 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine dlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
DLARHS
Definition dlarhs.f:205
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 dchksy_aa(dotype, nn, nval, nnb, nbval, nns, nsval, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
DCHKSY_AA
Definition dchksy_aa.f:170
subroutine derrsy(path, nunit)
DERRSY
Definition derrsy.f:55
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
Definition dlatb4.f:120
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DPOT02
Definition dpot02.f:127
subroutine dsyt01_aa(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
DSYT01
Definition dsyt01_aa.f:124
subroutine dsytrf_aa(uplo, n, a, lda, ipiv, work, lwork, info)
DSYTRF_AA
Definition dsytrf_aa.f:132
subroutine dsytrs_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
DSYTRS_AA
Definition dsytrs_aa.f:131
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103