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