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