SCALAPACK 2.2.2
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
pdblas2tst.f
Go to the documentation of this file.
1 BLOCK DATA
2 INTEGER NSUBS
3 parameter(nsubs = 7)
4 CHARACTER*7 SNAMES( NSUBS )
5 COMMON /snamec/snames
6 DATA snames/'PDGEMV ', 'PDSYMV ', 'PDTRMV ',
7 $ 'PDTRSV ', 'PDGER ', 'PDSYR ',
8 $ 'PDSYR2 '/
9 END BLOCK DATA
10
11 PROGRAM pdbla2tst
12*
13* -- PBLAS testing driver (version 2.0.2) --
14* Univ. of Tennessee, Univ. of California Berkeley, Univ. of Colorado Denver
15* May 1 2012
16*
17* Purpose
18* =======
19*
20* PDBLA2TST is the main testing program for the PBLAS Level 2 routines.
21*
22* The program must be driven by a short data file. An annotated exam-
23* ple of a data file can be obtained by deleting the first 3 characters
24* from the following 60 lines:
25* 'Level 2 PBLAS, Testing input file'
26* 'Intel iPSC/860 hypercube, gamma model.'
27* 'PDBLAS2TST.SUMM' output file name (if any)
28* 6 device out
29* F logical flag, T to stop on failures
30* F logical flag, T to test error exits
31* 0 verbosity, 0 for pass/fail, 1-3 for matrix dump on errors
32* 10 the leading dimension gap
33* 16.0 threshold value of test ratio
34* 10 value of the logical computational blocksize NB
35* 1 number of process grids (ordered pairs of P & Q)
36* 2 2 1 4 2 3 8 values of P
37* 2 2 4 1 3 2 1 values of Q
38* 1.0D0 value of ALPHA
39* 1.0D0 value of BETA
40* 2 number of tests problems
41* 'U' 'L' values of UPLO
42* 'N' 'T' values of TRANS
43* 'N' 'U' values of DIAG
44* 3 4 values of M
45* 3 4 values of N
46* 6 10 values of M_A
47* 6 10 values of N_A
48* 2 5 values of IMB_A
49* 2 5 values of INB_A
50* 2 5 values of MB_A
51* 2 5 values of NB_A
52* 0 1 values of RSRC_A
53* 0 0 values of CSRC_A
54* 1 1 values of IA
55* 1 1 values of JA
56* 6 10 values of M_X
57* 6 10 values of N_X
58* 2 5 values of IMB_X
59* 2 5 values of INB_X
60* 2 5 values of MB_X
61* 2 5 values of NB_X
62* 0 1 values of RSRC_X
63* 0 0 values of CSRC_X
64* 1 1 values of IX
65* 1 1 values of JX
66* 1 1 values of INCX
67* 6 10 values of M_Y
68* 6 10 values of N_Y
69* 2 5 values of IMB_Y
70* 2 5 values of INB_Y
71* 2 5 values of MB_Y
72* 2 5 values of NB_Y
73* 0 1 values of RSRC_Y
74* 0 0 values of CSRC_Y
75* 1 1 values of IY
76* 1 1 values of JY
77* 6 1 values of INCY
78* PDGEMV T put F for no test in the same column
79* PDSYMV T put F for no test in the same column
80* PDTRMV T put F for no test in the same column
81* PDTRSV T put F for no test in the same column
82* PDGER T put F for no test in the same column
83* PDSYR T put F for no test in the same column
84* PDSYR2 T put F for no test in the same column
85*
86* Internal Parameters
87* ===================
88*
89* TOTMEM INTEGER
90* TOTMEM is a machine-specific parameter indicating the maxi-
91* mum amount of available memory per process in bytes. The
92* user should customize TOTMEM to his platform. Remember to
93* leave room in memory for the operating system, the BLACS
94* buffer, etc. For example, on a system with 8 MB of memory
95* per process (e.g., one processor on an Intel iPSC/860), the
96* parameters we use are TOTMEM=6200000 (leaving 1.8 MB for OS,
97* code, BLACS buffer, etc). However, for PVM, we usually set
98* TOTMEM = 2000000. Some experimenting with the maximum value
99* of TOTMEM may be required. By default, TOTMEM is 2000000.
100*
101* DBLESZ INTEGER
102* DBLESZ indicates the length in bytes on the given platform
103* for a double precision real. By default, DBLESZ is set to
104* eight.
105*
106* MEM DOUBLE PRECISION array
107* MEM is an array of dimension TOTMEM / DBLESZ.
108* All arrays used by SCALAPACK routines are allocated from this
109* array MEM and referenced by pointers. The integer IPA, for
110* example, is a pointer to the starting element of MEM for the
111* matrix A.
112*
113* -- Written on April 1, 1998 by
114* Antoine Petitet, University of Tennessee, Knoxville 37996, USA.
115*
116* =====================================================================
117*
118* .. Parameters ..
119 INTEGER maxtests, maxgrids, gapmul, dblesz, totmem,
120 $ memsiz, nsubs
121 DOUBLE PRECISION one, padval, zero, rogue
122 parameter( maxtests = 20, maxgrids = 20, gapmul = 10,
123 $ dblesz = 8, totmem = 2000000,
124 $ memsiz = totmem / dblesz, zero = 0.0d+0,
125 $ one = 1.0d+0, padval = -9923.0d+0,
126 $ nsubs = 7, rogue = -1.0d+10 )
127 INTEGER block_cyclic_2d_inb, csrc_, ctxt_, dlen_,
128 $ dtype_, imb_, inb_, lld_, mb_, m_, nb_, n_,
129 $ rsrc_
130 parameter( block_cyclic_2d_inb = 2, dlen_ = 11,
131 $ dtype_ = 1, ctxt_ = 2, m_ = 3, n_ = 4,
132 $ imb_ = 5, inb_ = 6, mb_ = 7, nb_ = 8,
133 $ rsrc_ = 9, csrc_ = 10, lld_ = 11 )
134* ..
135* .. Local Scalars ..
136 LOGICAL errflg, sof, tee
137 CHARACTER*1 aform, diag, diagdo, trans, uplo
138 INTEGER csrca, csrcx, csrcy, i, ia, iam, iaseed, ictxt,
139 $ igap, imba, imbx, imby, imida, imidx, imidy,
140 $ inba, inbx, inby, incx, incy, ipa, ipg, ipmata,
141 $ ipmatx, ipmaty, iposta, ipostx, iposty, iprea,
142 $ iprex, iprey, ipx, ipy, iverb, ix, ixseed, iy,
143 $ iyseed, j, ja, jx, jy, k, lda, ldx, ldy, m, ma,
144 $ mba, mbx, mby, memreqd, mpa, mpx, mpy, mx, my,
145 $ mycol, myrow, n, na, nba, nbx, nby, ncola,
146 $ ngrids, nlx, nly, nout, npcol, nprocs, nprow,
147 $ nqa, nqx, nqy, nrowa, ntests, nx, ny, offd,
148 $ rsrca, rsrcx, rsrcy, tskip, tstcnt
149 REAL thresh
150 DOUBLE PRECISION alpha, beta, scale
151* ..
152* .. Local Arrays ..
153 LOGICAL ltest( nsubs ), ycheck( nsubs )
154 CHARACTER*1 diagval( maxtests ), tranval( maxtests ),
155 $ uploval( maxtests )
156 CHARACTER*80 outfile
157 INTEGER cscaval( maxtests ), cscxval( maxtests ),
158 $ cscyval( maxtests ), desca( dlen_ ),
159 $ descar( dlen_ ), descx( dlen_ ),
160 $ descxr( dlen_ ), descy( dlen_ ),
161 $ descyr( dlen_ ), iaval( maxtests ), ierr( 6 ),
162 $ imbaval( maxtests ), imbxval( maxtests ),
163 $ imbyval( maxtests ), inbaval( maxtests ),
164 $ inbxval( maxtests ), inbyval( maxtests ),
165 $ incxval( maxtests ), incyval( maxtests ),
166 $ ixval( maxtests ), iyval( maxtests ),
167 $ javal( maxtests ), jxval( maxtests ),
168 $ jyval( maxtests )
169 INTEGER kfail( nsubs ), kpass( nsubs ), kskip( nsubs ),
170 $ ktests( nsubs ), maval( maxtests ),
171 $ mbaval( maxtests ), mbxval( maxtests ),
172 $ mbyval( maxtests ), mval( maxtests ),
173 $ mxval( maxtests ), myval( maxtests ),
174 $ naval( maxtests ), nbaval( maxtests ),
175 $ nbxval( maxtests ), nbyval( maxtests ),
176 $ nval( maxtests ), nxval( maxtests ),
177 $ nyval( maxtests ), pval( maxtests ),
178 $ qval( maxtests ), rscaval( maxtests ),
179 $ rscxval( maxtests ), rscyval( maxtests )
180 DOUBLE PRECISION mem( memsiz )
181* ..
182* .. External Subroutines ..
183 EXTERNAL blacs_exit, blacs_get, blacs_gridexit,
184 $ blacs_gridinfo, blacs_gridinit, blacs_pinfo,
188 $ pdchkarg2, pdchkvout, pdgemv, pdger, pdlagen,
189 $ pdlascal, pdlaset, pdmprnt, pdsymv, pdsyr,
190 $ pdsyr2, pdtrmv, pdtrsv, pdvprnt, pmdescchk,
192* ..
193* .. External Functions ..
194 LOGICAL lsame
195 EXTERNAL lsame
196* ..
197* .. Intrinsic Functions ..
198 INTRINSIC abs, dble, max, mod
199* ..
200* .. Common Blocks ..
201 CHARACTER*7 snames( nsubs )
202 LOGICAL abrtflg
203 INTEGER info, nblog
204 COMMON /snamec/snames
205 COMMON /infoc/info, nblog
206 COMMON /pberrorc/nout, abrtflg
207* ..
208* .. Data Statements ..
209 DATA ycheck/.true., .true., .false., .false.,
210 $ .true., .false., .true./
211* ..
212* .. Executable Statements ..
213*
214* Initialization
215*
216* Set flag so that the PBLAS error handler won't abort on errors, so
217* that the tester will detect unsupported operations.
218*
219 abrtflg = .false.
220*
221* So far no error, will become true as soon as one error is found.
222*
223 errflg = .false.
224*
225* Test counters
226*
227 tskip = 0
228 tstcnt = 0
229*
230* Seeds for random matrix generations.
231*
232 iaseed = 100
233 ixseed = 200
234 iyseed = 300
235*
236* So far no tests have been performed.
237*
238 DO 10 i = 1, nsubs
239 kpass( i ) = 0
240 kskip( i ) = 0
241 kfail( i ) = 0
242 ktests( i ) = 0
243 10 CONTINUE
244*
245* Get starting information
246*
247 CALL blacs_pinfo( iam, nprocs )
248 CALL pdbla2tstinfo( outfile, nout, ntests, diagval, tranval,
249 $ uploval, mval, nval, maval, naval, imbaval,
250 $ mbaval, inbaval, nbaval, rscaval, cscaval,
251 $ iaval, javal, mxval, nxval, imbxval, mbxval,
252 $ inbxval, nbxval, rscxval, cscxval, ixval,
253 $ jxval, incxval, myval, nyval, imbyval,
254 $ mbyval, inbyval, nbyval, rscyval, cscyval,
255 $ iyval, jyval, incyval, maxtests, ngrids,
256 $ pval, maxgrids, qval, maxgrids, nblog, ltest,
257 $ sof, tee, iam, igap, iverb, nprocs, thresh,
258 $ alpha, beta, mem )
259*
260 IF( iam.EQ.0 ) THEN
261 WRITE( nout, fmt = 9975 )
262 WRITE( nout, fmt = * )
263 END IF
264*
265* If TEE is set then Test Error Exits of routines.
266*
267 IF( tee )
268 $ CALL pdblas2tstchke( ltest, nout, nprocs )
269*
270* Loop over different process grids
271*
272 DO 60 i = 1, ngrids
273*
274 nprow = pval( i )
275 npcol = qval( i )
276*
277* Make sure grid information is correct
278*
279 ierr( 1 ) = 0
280 IF( nprow.LT.1 ) THEN
281 IF( iam.EQ.0 )
282 $ WRITE( nout, fmt = 9999 ) 'GRID SIZE', 'NPROW', nprow
283 ierr( 1 ) = 1
284 ELSE IF( npcol.LT.1 ) THEN
285 IF( iam.EQ.0 )
286 $ WRITE( nout, fmt = 9999 ) 'GRID SIZE', 'NPCOL', npcol
287 ierr( 1 ) = 1
288 ELSE IF( nprow*npcol.GT.nprocs ) THEN
289 IF( iam.EQ.0 )
290 $ WRITE( nout, fmt = 9998 ) nprow*npcol, nprocs
291 ierr( 1 ) = 1
292 END IF
293*
294 IF( ierr( 1 ).GT.0 ) THEN
295 IF( iam.EQ.0 )
296 $ WRITE( nout, fmt = 9997 ) 'GRID'
297 tskip = tskip + 1
298 GO TO 60
299 END IF
300*
301* Define process grid
302*
303 CALL blacs_get( -1, 0, ictxt )
304 CALL blacs_gridinit( ictxt, 'Row-major', nprow, npcol )
305 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
306*
307* Go to bottom of process grid loop if this case doesn't use my
308* process
309*
310 IF( myrow.GE.nprow .OR. mycol.GE.npcol )
311 $ GO TO 60
312*
313* Loop over number of tests
314*
315 DO 50 j = 1, ntests
316*
317* Get the test parameters
318*
319 diag = diagval( j )
320 trans = tranval( j )
321 uplo = uploval( j )
322*
323 m = mval( j )
324 n = nval( j )
325*
326 ma = maval( j )
327 na = naval( j )
328 imba = imbaval( j )
329 inba = inbaval( j )
330 mba = mbaval( j )
331 nba = nbaval( j )
332 rsrca = rscaval( j )
333 csrca = cscaval( j )
334 ia = iaval( j )
335 ja = javal( j )
336*
337 mx = mxval( j )
338 nx = nxval( j )
339 imbx = imbxval( j )
340 inbx = inbxval( j )
341 mbx = mbxval( j )
342 nbx = nbxval( j )
343 rsrcx = rscxval( j )
344 csrcx = cscxval( j )
345 ix = ixval( j )
346 jx = jxval( j )
347 incx = incxval( j )
348*
349 my = myval( j )
350 ny = nyval( j )
351 imby = imbyval( j )
352 inby = inbyval( j )
353 mby = mbyval( j )
354 nby = nbyval( j )
355 rsrcy = rscyval( j )
356 csrcy = cscyval( j )
357 iy = iyval( j )
358 jy = jyval( j )
359 incy = incyval( j )
360*
361 IF( iam.EQ.0 ) THEN
362 tstcnt = tstcnt + 1
363 WRITE( nout, fmt = * )
364 WRITE( nout, fmt = 9996 ) tstcnt, nprow, npcol
365 WRITE( nout, fmt = * )
366*
367 WRITE( nout, fmt = 9995 )
368 WRITE( nout, fmt = 9994 )
369 WRITE( nout, fmt = 9995 )
370 WRITE( nout, fmt = 9993 ) m, n, uplo, trans, diag
371*
372 WRITE( nout, fmt = 9995 )
373 WRITE( nout, fmt = 9992 )
374 WRITE( nout, fmt = 9995 )
375 WRITE( nout, fmt = 9991 ) ia, ja, ma, na, imba, inba,
376 $ mba, nba, rsrca, csrca
377*
378 WRITE( nout, fmt = 9995 )
379 WRITE( nout, fmt = 9990 )
380 WRITE( nout, fmt = 9995 )
381 WRITE( nout, fmt = 9989 ) ix, jx, mx, nx, imbx, inbx,
382 $ mbx, nbx, rsrcx, csrcx, incx
383*
384 WRITE( nout, fmt = 9995 )
385 WRITE( nout, fmt = 9988 )
386 WRITE( nout, fmt = 9995 )
387 WRITE( nout, fmt = 9989 ) iy, jy, my, ny, imby, inby,
388 $ mby, nby, rsrcy, csrcy, incy
389*
390 WRITE( nout, fmt = 9995 )
391*
392 END IF
393*
394* Check the validity of the input test parameters
395*
396 IF( .NOT.lsame( uplo, 'U' ).AND.
397 $ .NOT.lsame( uplo, 'L' ) ) THEN
398 IF( iam.EQ.0 )
399 $ WRITE( nout, fmt = 9997 ) 'UPLO'
400 tskip = tskip + 1
401 GO TO 40
402 END IF
403*
404 IF( .NOT.lsame( trans, 'N' ).AND.
405 $ .NOT.lsame( trans, 'T' ).AND.
406 $ .NOT.lsame( trans, 'C' ) ) THEN
407 IF( iam.EQ.0 )
408 $ WRITE( nout, fmt = 9997 ) 'TRANS'
409 tskip = tskip + 1
410 GO TO 40
411 END IF
412*
413 IF( .NOT.lsame( diag , 'U' ).AND.
414 $ .NOT.lsame( diag , 'N' ) )THEN
415 IF( iam.EQ.0 )
416 $ WRITE( nout, fmt = 9997 ) trans
417 WRITE( nout, fmt = 9997 ) 'DIAG'
418 tskip = tskip + 1
419 GO TO 40
420 END IF
421*
422* Check and initialize the matrix descriptors
423*
424 CALL pmdescchk( ictxt, nout, 'A', desca,
425 $ block_cyclic_2d_inb, ma, na, imba, inba,
426 $ mba, nba, rsrca, csrca, mpa, nqa, iprea,
427 $ imida, iposta, igap, gapmul, ierr( 1 ) )
428 CALL pvdescchk( ictxt, nout, 'X', descx,
429 $ block_cyclic_2d_inb, mx, nx, imbx, inbx,
430 $ mbx, nbx, rsrcx, csrcx, incx, mpx, nqx,
431 $ iprex, imidx, ipostx, igap, gapmul,
432 $ ierr( 2 ) )
433 CALL pvdescchk( ictxt, nout, 'Y', descy,
434 $ block_cyclic_2d_inb, my, ny, imby, inby,
435 $ mby, nby, rsrcy, csrcy, incy, mpy, nqy,
436 $ iprey, imidy, iposty, igap, gapmul,
437 $ ierr( 3 ) )
438*
439 IF( ierr( 1 ).GT.0 .OR. ierr( 2 ).GT.0 .OR.
440 $ ierr( 3 ).GT.0 ) THEN
441 tskip = tskip + 1
442 GO TO 40
443 END IF
444*
445 lda = max( 1, ma )
446 ldx = max( 1, mx )
447 ldy = max( 1, my )
448*
449* Assign pointers into MEM for matrices corresponding to
450* the distributed matrices A, X and Y.
451*
452 ipa = iprea + 1
453 ipx = ipa + desca( lld_ )*nqa + iposta + iprex
454 ipy = ipx + descx( lld_ )*nqx + ipostx + iprey
455 ipmata = ipy + descy( lld_ )*nqy + iposty
456 ipmatx = ipmata + ma*na
457 ipmaty = ipmatx + mx*nx
458 ipg = ipmaty + max( mx*nx, my*ny )
459*
460* Check if sufficient memory.
461* Requirement = mem for local part of parallel matrices +
462* mem for whole matrices for comp. check +
463* mem for recving comp. check error vals.
464*
465 memreqd = ipg + max( m, n ) - 1 +
466 $ max( max( imba, mba ),
467 $ max( max( imbx, mbx ),
468 $ max( imby, mby ) ) )
469 ierr( 1 ) = 0
470 IF( memreqd.GT.memsiz ) THEN
471 IF( iam.EQ.0 )
472 $ WRITE( nout, fmt = 9986 ) memreqd*dblesz
473 ierr( 1 ) = 1
474 END IF
475*
476* Check all processes for an error
477*
478 CALL igsum2d( ictxt, 'All', ' ', 1, 1, ierr, 1, -1, 0 )
479*
480 IF( ierr( 1 ).GT.0 ) THEN
481 IF( iam.EQ.0 )
482 $ WRITE( nout, fmt = 9987 )
483 tskip = tskip + 1
484 GO TO 40
485 END IF
486*
487* Loop over all PBLAS 2 routines
488*
489 DO 30 k = 1, nsubs
490*
491* Continue only if this subroutine has to be tested.
492*
493 IF( .NOT.ltest( k ) )
494 $ GO TO 30
495*
496 IF( iam.EQ.0 ) THEN
497 WRITE( nout, fmt = * )
498 WRITE( nout, fmt = 9985 ) snames( k )
499 END IF
500*
501* Define the size of the operands
502*
503 IF( k.EQ.1 ) THEN
504 nrowa = m
505 ncola = n
506 IF( lsame( trans, 'N' ) ) THEN
507 nlx = n
508 nly = m
509 ELSE
510 nlx = m
511 nly = n
512 END IF
513 ELSE IF( k.EQ.5 ) THEN
514 nrowa = m
515 ncola = n
516 nlx = m
517 nly = n
518 ELSE
519 nrowa = n
520 ncola = n
521 nlx = n
522 nly = n
523 END IF
524*
525* Check the validity of the operand sizes
526*
527 CALL pmdimchk( ictxt, nout, nrowa, ncola, 'A', ia, ja,
528 $ desca, ierr( 1 ) )
529 CALL pvdimchk( ictxt, nout, nlx, 'X', ix, jx, descx,
530 $ incx, ierr( 2 ) )
531 CALL pvdimchk( ictxt, nout, nly, 'Y', iy, jy, descy,
532 $ incy, ierr( 3 ) )
533*
534 IF( ierr( 1 ).NE.0 .OR. ierr( 2 ).NE.0 .OR.
535 $ ierr( 3 ).NE.0 ) THEN
536 kskip( k ) = kskip( k ) + 1
537 GO TO 30
538 END IF
539*
540* Generate distributed matrices A, X and Y
541*
542 IF( k.EQ.2 .OR. k.EQ.6 .OR. k.EQ.7 ) THEN
543 aform = 'S'
544 diagdo = 'N'
545 offd = ia - ja
546 ELSE IF( ( k.EQ.4 ).AND.( lsame( diag, 'N' ) ) ) THEN
547 aform = 'N'
548 diagdo = 'D'
549 offd = ia - ja
550 ELSE
551 aform = 'N'
552 diagdo = 'N'
553 offd = 0
554 END IF
555*
556 CALL pdlagen( .false., aform, diagdo, offd, ma, na,
557 $ 1, 1, desca, iaseed, mem( ipa ),
558 $ desca( lld_ ) )
559 CALL pdlagen( .false., 'None', 'No diag', 0, mx, nx, 1,
560 $ 1, descx, ixseed, mem( ipx ),
561 $ descx( lld_ ) )
562 IF( ycheck( k ) )
563 $ CALL pdlagen( .false., 'None', 'No diag', 0, my, ny,
564 $ 1, 1, descy, iyseed, mem( ipy ),
565 $ descy( lld_ ) )
566*
567* Generate entire matrices on each process.
568*
569 CALL pb_descset2( descar, ma, na, imba, inba, mba, nba,
570 $ -1, -1, ictxt, max( 1, ma ) )
571 CALL pdlagen( .false., aform, diagdo, offd, ma, na,
572 $ 1, 1, descar, iaseed, mem( ipmata ),
573 $ descar( lld_ ) )
574 CALL pb_descset2( descxr, mx, nx, imbx, inbx, mbx, nbx,
575 $ -1, -1, ictxt, max( 1, mx ) )
576 CALL pdlagen( .false., 'None', 'No diag', 0, mx, nx, 1,
577 $ 1, descxr, ixseed, mem( ipmatx ),
578 $ descxr( lld_ ) )
579 IF( ycheck( k ) ) THEN
580*
581 CALL pb_descset2( descyr, my, ny, imby, inby, mby,
582 $ nby, -1, -1, ictxt, max( 1, my ) )
583 CALL pdlagen( .false., 'None', 'No diag', 0, my, ny,
584 $ 1, 1, descyr, iyseed, mem( ipmaty ),
585 $ descyr( lld_ ) )
586*
587 ELSE
588*
589* If Y is not needed, generate a copy of X instead
590*
591 CALL pb_descset2( descyr, mx, nx, imbx, inbx, mbx,
592 $ nbx, -1, -1, ictxt, max( 1, mx ) )
593 CALL pdlagen( .false., 'None', 'No diag', 0, mx, nx,
594 $ 1, 1, descyr, ixseed, mem( ipmaty ),
595 $ descyr( lld_ ) )
596*
597 END IF
598*
599* Zero non referenced part of the matrices A
600*
601 IF( ( k.EQ.2 .OR. k.EQ.6 .OR. k.EQ.7 ).AND.
602 $ ( max( nrowa, ncola ).GT.1 ) ) THEN
603*
604* The distributed matrix A is symmetric
605*
606 IF( lsame( uplo, 'L' ) ) THEN
607*
608* Zeros the strict upper triangular part of A.
609*
610 CALL pdlaset( 'Upper', nrowa-1, ncola-1, rogue,
611 $ rogue, mem( ipa ), ia, ja+1, desca )
612 IF( k.NE.2 ) THEN
613 CALL pb_dlaset( 'Upper', nrowa-1, ncola-1, 0,
614 $ rogue, rogue,
615 $ mem( ipmata+ia-1+ja*lda ), lda )
616 END IF
617*
618 ELSE IF( lsame( uplo, 'U' ) ) THEN
619*
620* Zeros the strict lower triangular part of A.
621*
622 CALL pdlaset( 'Lower', nrowa-1, ncola-1, rogue,
623 $ rogue, mem( ipa ), ia+1, ja, desca )
624 IF( k.NE.2 ) THEN
625 CALL pb_dlaset( 'Lower', nrowa-1, ncola-1, 0,
626 $ rogue, rogue,
627 $ mem( ipmata+ia+(ja-1)*lda ),
628 $ lda )
629 END IF
630*
631 END IF
632*
633 ELSE IF( k.EQ.3 .OR. k.EQ.4 ) THEN
634*
635 IF( lsame( uplo, 'L' ) ) THEN
636*
637* The distributed matrix A is lower triangular
638*
639 IF( lsame( diag, 'N' ) ) THEN
640*
641 IF( max( nrowa, ncola ).GT.1 ) THEN
642 CALL pdlaset( 'Upper', nrowa-1, ncola-1,
643 $ rogue, rogue, mem( ipa ), ia,
644 $ ja+1, desca )
645 CALL pb_dlaset( 'Upper', nrowa-1, ncola-1, 0,
646 $ zero, zero,
647 $ mem( ipmata+ia-1+ja*lda ),
648 $ lda )
649 END IF
650*
651 ELSE
652*
653 CALL pdlaset( 'Upper', nrowa, ncola, rogue, one,
654 $ mem( ipa ), ia, ja, desca )
655 CALL pb_dlaset( 'Upper', nrowa, ncola, 0, zero,
656 $ one,
657 $ mem( ipmata+ia-1+(ja-1)*lda ),
658 $ lda )
659 IF( ( k.EQ.4 ).AND.
660 $ ( max( nrowa, ncola ).GT.1 ) ) THEN
661 scale = one / dble( max( nrowa, ncola ) )
662 CALL pdlascal( 'Lower', nrowa-1, ncola-1,
663 $ scale, mem( ipa ), ia+1, ja,
664 $ desca )
665 CALL pb_dlascal( 'Lower', nrowa-1, ncola-1,
666 $ 0, scale,
667 $ mem( ipmata+ia+(ja-1)*lda ),
668 $ lda )
669 END IF
670*
671 END IF
672*
673 ELSE IF( lsame( uplo, 'U' ) ) THEN
674*
675* The distributed matrix A is upper triangular
676*
677 IF( lsame( diag, 'N' ) ) THEN
678*
679 IF( max( nrowa, ncola ).GT.1 ) THEN
680 CALL pdlaset( 'Lower', nrowa-1, ncola-1,
681 $ rogue, rogue, mem( ipa ), ia+1,
682 $ ja, desca )
683 CALL pb_dlaset( 'Lower', nrowa-1, ncola-1, 0,
684 $ zero, zero,
685 $ mem( ipmata+ia+(ja-1)*lda ),
686 $ lda )
687 END IF
688*
689 ELSE
690*
691 CALL pdlaset( 'Lower', nrowa, ncola, rogue, one,
692 $ mem( ipa ), ia, ja, desca )
693 CALL pb_dlaset( 'Lower', nrowa, ncola, 0, zero,
694 $ one,
695 $ mem( ipmata+ia-1+(ja-1)*lda ),
696 $ lda )
697 IF( ( k.EQ.4 ).AND.
698 $ ( max( nrowa, ncola ).GT.1 ) ) THEN
699 scale = one / dble( max( nrowa, ncola ) )
700 CALL pdlascal( 'Upper', nrowa-1, ncola-1,
701 $ scale, mem( ipa ), ia, ja+1,
702 $ desca )
703 CALL pb_dlascal( 'Upper', nrowa-1, ncola-1,
704 $ 0, scale,
705 $ mem( ipmata+ia-1+ja*lda ), lda )
706 END IF
707*
708 END IF
709*
710 END IF
711*
712 END IF
713*
714* Pad the guard zones of A, X and Y
715*
716 CALL pb_dfillpad( ictxt, mpa, nqa, mem( ipa-iprea ),
717 $ desca( lld_ ), iprea, iposta, padval )
718*
719 CALL pb_dfillpad( ictxt, mpx, nqx, mem( ipx-iprex ),
720 $ descx( lld_ ), iprex, ipostx, padval )
721*
722 IF( ycheck( k ) ) THEN
723 CALL pb_dfillpad( ictxt, mpy, nqy, mem( ipy-iprey ),
724 $ descy( lld_ ), iprey, iposty,
725 $ padval )
726 END IF
727*
728* Initialize the check for INPUT-only arguments.
729*
730 info = 0
731 CALL pdchkarg2( ictxt, nout, snames( k ), uplo, trans,
732 $ diag, m, n, alpha, ia, ja, desca, ix,
733 $ jx, descx, incx, beta, iy, jy, descy,
734 $ incy, info )
735*
736* Print initial parallel data if IVERB >= 2.
737*
738 IF( iverb.EQ.2 ) THEN
739 CALL pb_pdlaprnt( nrowa, ncola, mem( ipa ), ia, ja,
740 $ desca, 0, 0, 'PARALLEL_INITIAL_A',
741 $ nout, mem( ipg ) )
742 ELSE IF( iverb.GE.3 ) THEN
743 CALL pb_pdlaprnt( ma, na, mem( ipa ), 1, 1, desca, 0,
744 $ 0, 'PARALLEL_INITIAL_A', nout,
745 $ mem( ipg ) )
746 END IF
747*
748 IF( iverb.EQ.2 ) THEN
749 IF( incx.EQ.descx( m_ ) ) THEN
750 CALL pb_pdlaprnt( 1, nlx, mem( ipx ), ix, jx,
751 $ descx, 0, 0,
752 $ 'PARALLEL_INITIAL_X', nout,
753 $ mem( ipg ) )
754 ELSE
755 CALL pb_pdlaprnt( nlx, 1, mem( ipx ), ix, jx,
756 $ descx, 0, 0,
757 $ 'PARALLEL_INITIAL_X', nout,
758 $ mem( ipg ) )
759 END IF
760 ELSE IF( iverb.GE.3 ) THEN
761 CALL pb_pdlaprnt( mx, nx, mem( ipx ), 1, 1, descx, 0,
762 $ 0, 'PARALLEL_INITIAL_X', nout,
763 $ mem( ipg ) )
764 END IF
765*
766 IF( ycheck( k ) ) THEN
767 IF( iverb.EQ.2 ) THEN
768 IF( incy.EQ.descy( m_ ) ) THEN
769 CALL pb_pdlaprnt( 1, nly, mem( ipy ), iy, jy,
770 $ descy, 0, 0,
771 $ 'PARALLEL_INITIAL_Y', nout,
772 $ mem( ipg ) )
773 ELSE
774 CALL pb_pdlaprnt( nly, 1, mem( ipy ), iy, jy,
775 $ descy, 0, 0,
776 $ 'PARALLEL_INITIAL_Y', nout,
777 $ mem( ipg ) )
778 END IF
779 ELSE IF( iverb.GE.3 ) THEN
780 CALL pb_pdlaprnt( my, ny, mem( ipy ), 1, 1, descy,
781 $ 0, 0, 'PARALLEL_INITIAL_Y', nout,
782 $ mem( ipg ) )
783 END IF
784 END IF
785*
786* Call the Level 2 PBLAS routine
787*
788 info = 0
789 IF( k.EQ.1 ) THEN
790*
791* Test PDGEMV
792*
793 CALL pdgemv( trans, m, n, alpha, mem( ipa ), ia, ja,
794 $ desca, mem( ipx ), ix, jx, descx, incx,
795 $ beta, mem( ipy ), iy, jy, descy, incy )
796*
797 ELSE IF( k.EQ.2 ) THEN
798*
799* Test PDSYMV
800*
801 CALL pdsymv( uplo, n, alpha, mem( ipa ), ia, ja,
802 $ desca, mem( ipx ), ix, jx, descx, incx,
803 $ beta, mem( ipy ), iy, jy, descy, incy )
804*
805 ELSE IF( k.EQ.3 ) THEN
806*
807* Test PDTRMV
808*
809 CALL pdtrmv( uplo, trans, diag, n, mem( ipa ), ia, ja,
810 $ desca, mem( ipx ), ix, jx, descx, incx )
811*
812 ELSE IF( k.EQ.4 ) THEN
813*
814* Test PDTRSV
815*
816 CALL pdtrsv( uplo, trans, diag, n, mem( ipa ), ia, ja,
817 $ desca, mem( ipx ), ix, jx, descx, incx )
818*
819 ELSE IF( k.EQ.5 ) THEN
820*
821* Test PDGER
822*
823 CALL pdger( m, n, alpha, mem( ipx ), ix, jx, descx,
824 $ incx, mem( ipy ), iy, jy, descy, incy,
825 $ mem( ipa ), ia, ja, desca )
826*
827 ELSE IF( k.EQ.6 ) THEN
828*
829* Test PDSYR
830*
831 CALL pdsyr( uplo, n, alpha, mem( ipx ), ix, jx, descx,
832 $ incx, mem( ipa ), ia, ja, desca )
833*
834 ELSE IF( k.EQ.7 ) THEN
835*
836* Test PDSYR2
837*
838 CALL pdsyr2( uplo, n, alpha, mem( ipx ), ix, jx,
839 $ descx, incx, mem( ipy ), iy, jy, descy,
840 $ incy, mem( ipa ), ia, ja, desca )
841*
842 END IF
843*
844* Check if the operation has been performed.
845*
846 IF( info.NE.0 ) THEN
847 kskip( k ) = kskip( k ) + 1
848 IF( iam.EQ.0 )
849 $ WRITE( nout, fmt = 9974 ) info
850 GO TO 30
851 END IF
852*
853* Check padding
854*
855 CALL pb_dchekpad( ictxt, snames( k ), mpa, nqa,
856 $ mem( ipa-iprea ), desca( lld_ ), iprea,
857 $ iposta, padval )
858*
859 CALL pb_dchekpad( ictxt, snames( k ), mpx, nqx,
860 $ mem( ipx-iprex ), descx( lld_ ), iprex,
861 $ ipostx, padval )
862*
863 IF( ycheck( k ) ) THEN
864 CALL pb_dchekpad( ictxt, snames( k ), mpy, nqy,
865 $ mem( ipy-iprey ), descy( lld_ ),
866 $ iprey, iposty, padval )
867 END IF
868*
869* Check the computations
870*
871 CALL pdblas2tstchk( ictxt, nout, k, uplo, trans, diag, m,
872 $ n, alpha, mem( ipmata ), mem( ipa ),
873 $ ia, ja, desca, mem( ipmatx ),
874 $ mem( ipx ), ix, jx, descx, incx,
875 $ beta, mem( ipmaty ), mem( ipy ), iy,
876 $ jy, descy, incy, thresh, rogue,
877 $ mem( ipg ), info )
878 IF( mod( info, 2 ).EQ.1 ) THEN
879 ierr( 1 ) = 1
880 ELSE IF( mod( info / 2, 2 ).EQ.1 ) THEN
881 ierr( 2 ) = 1
882 ELSE IF( mod( info / 4, 2 ).EQ.1 ) THEN
883 ierr( 3 ) = 1
884 ELSE IF( info.NE.0 ) THEN
885 ierr( 1 ) = 1
886 ierr( 2 ) = 1
887 ierr( 3 ) = 1
888 END IF
889*
890* Check input-only scalar arguments
891*
892 info = 1
893 CALL pdchkarg2( ictxt, nout, snames( k ), uplo, trans,
894 $ diag, m, n, alpha, ia, ja, desca, ix,
895 $ jx, descx, incx, beta, iy, jy, descy,
896 $ incy, info )
897*
898* Check input-only array arguments
899*
900 CALL pdchkmout( nrowa, ncola, mem( ipmata ), mem( ipa ),
901 $ ia, ja, desca, ierr( 4 ) )
902 CALL pdchkvout( nlx, mem( ipmatx ), mem( ipx ), ix, jx,
903 $ descx, incx, ierr( 5 ) )
904*
905 IF( ierr( 4 ).NE.0 ) THEN
906 IF( iam.EQ.0 )
907 $ WRITE( nout, fmt = 9982 ) 'PARALLEL_A',
908 $ snames( k )
909 END IF
910*
911 IF( ierr( 5 ).NE.0 ) THEN
912 IF( iam.EQ.0 )
913 $ WRITE( nout, fmt = 9982 ) 'PARALLEL_X',
914 $ snames( k )
915 END IF
916*
917 IF( ycheck( k ) ) THEN
918 CALL pdchkvout( nly, mem( ipmaty ), mem( ipy ), iy,
919 $ jy, descy, incy, ierr( 6 ) )
920 IF( ierr( 6 ).NE.0 ) THEN
921 IF( iam.EQ.0 )
922 $ WRITE( nout, fmt = 9982 ) 'PARALLEL_Y',
923 $ snames( k )
924 END IF
925 END IF
926*
927* Only node 0 prints computational test result
928*
929 IF( info.NE.0 .OR. ierr( 1 ).NE.0 .OR.
930 $ ierr( 2 ).NE.0 .OR. ierr( 3 ).NE.0 .OR.
931 $ ierr( 4 ).NE.0 .OR. ierr( 5 ).NE.0 .OR.
932 $ ierr( 6 ).NE.0 ) THEN
933 IF( iam.EQ.0 )
934 $ WRITE( nout, fmt = 9984 ) snames( k )
935 kfail( k ) = kfail( k ) + 1
936 errflg = .true.
937 ELSE
938 IF( iam.EQ.0 )
939 $ WRITE( nout, fmt = 9983 ) snames( k )
940 kpass( k ) = kpass( k ) + 1
941 END IF
942*
943* Dump matrix if IVERB >= 1 and error.
944*
945 IF( iverb.GE.1 .AND. errflg ) THEN
946 IF( ierr( 4 ).NE.0 .OR. iverb.GE.3 ) THEN
947 CALL pdmprnt( ictxt, nout, ma, na, mem( ipmata ),
948 $ lda, 0, 0, 'SERIAL_A' )
949 CALL pb_pdlaprnt( ma, na, mem( ipa ), 1, 1, desca,
950 $ 0, 0, 'PARALLEL_A', nout,
951 $ mem( ipmata ) )
952 ELSE IF( ierr( 1 ).NE.0 ) THEN
953 IF( ( nrowa.GT.0 ).AND.( ncola.GT.0 ) )
954 $ CALL pdmprnt( ictxt, nout, nrowa, ncola,
955 $ mem( ipmata+ia-1+(ja-1)*lda ),
956 $ lda, 0, 0, 'SERIAL_A' )
957 CALL pb_pdlaprnt( nrowa, ncola, mem( ipa ), ia, ja,
958 $ desca, 0, 0, 'PARALLEL_A',
959 $ nout, mem( ipmata ) )
960 END IF
961 IF( ierr( 5 ).NE.0 .OR. iverb.GE.3 ) THEN
962 CALL pdmprnt( ictxt, nout, mx, nx, mem( ipmatx ),
963 $ ldx, 0, 0, 'SERIAL_X' )
964 CALL pb_pdlaprnt( mx, nx, mem( ipx ), 1, 1, descx,
965 $ 0, 0, 'PARALLEL_X', nout,
966 $ mem( ipmatx ) )
967 ELSE IF( ierr( 2 ).NE.0 ) THEN
968 IF( nlx.GT.0 )
969 $ CALL pdvprnt( ictxt, nout, nlx,
970 $ mem( ipmatx+ix-1+(jx-1)*ldx ),
971 $ incx, 0, 0, 'SERIAL_X' )
972 IF( incx.EQ.descx( m_ ) ) THEN
973 CALL pb_pdlaprnt( 1, nlx, mem( ipx ), ix, jx,
974 $ descx, 0, 0, 'PARALLEL_X',
975 $ nout, mem( ipmatx ) )
976 ELSE
977 CALL pb_pdlaprnt( nlx, 1, mem( ipx ), ix, jx,
978 $ descx, 0, 0, 'PARALLEL_X',
979 $ nout, mem( ipmatx ) )
980 END IF
981 END IF
982 IF( ycheck( k ) ) THEN
983 IF( ierr( 6 ).NE.0 .OR. iverb.GE.3 ) THEN
984 CALL pdmprnt( ictxt, nout, my, ny,
985 $ mem( ipmaty ), ldy, 0, 0,
986 $ 'SERIAL_Y' )
987 CALL pb_pdlaprnt( my, ny, mem( ipy ), 1, 1,
988 $ descy, 0, 0, 'PARALLEL_Y',
989 $ nout, mem( ipmatx ) )
990 ELSE IF( ierr( 3 ).NE.0 ) THEN
991 IF( nly.GT.0 )
992 $ CALL pdvprnt( ictxt, nout, nly,
993 $ mem( ipmaty+iy-1+(jy-1)*ldy ),
994 $ incy, 0, 0, 'SERIAL_Y' )
995 IF( incy.EQ.descy( m_ ) ) THEN
996 CALL pb_pdlaprnt( 1, nly, mem( ipy ), iy, jy,
997 $ descy, 0, 0, 'PARALLEL_Y',
998 $ nout, mem( ipmatx ) )
999 ELSE
1000 CALL pb_pdlaprnt( nly, 1, mem( ipy ), iy, jy,
1001 $ descy, 0, 0, 'PARALLEL_Y',
1002 $ nout, mem( ipmatx ) )
1003 END IF
1004 END IF
1005 END IF
1006 END IF
1007*
1008* Leave if error and "Stop On Failure"
1009*
1010 IF( sof.AND.errflg )
1011 $ GO TO 70
1012*
1013 30 CONTINUE
1014*
1015 40 IF( iam.EQ.0 ) THEN
1016 WRITE( nout, fmt = * )
1017 WRITE( nout, fmt = 9981 ) j
1018 END IF
1019*
1020 50 CONTINUE
1021*
1022 CALL blacs_gridexit( ictxt )
1023*
1024 60 CONTINUE
1025*
1026* Come here, if error and "Stop On Failure"
1027*
1028 70 CONTINUE
1029*
1030* Before printing out final stats, add TSKIP to all skips
1031*
1032 DO 80 i = 1, nsubs
1033 IF( ltest( i ) ) THEN
1034 kskip( i ) = kskip( i ) + tskip
1035 ktests( i ) = kskip( i ) + kfail( i ) + kpass( i )
1036 END IF
1037 80 CONTINUE
1038*
1039* Print results
1040*
1041 IF( iam.EQ.0 ) THEN
1042 WRITE( nout, fmt = * )
1043 WRITE( nout, fmt = 9977 )
1044 WRITE( nout, fmt = * )
1045 WRITE( nout, fmt = 9979 )
1046 WRITE( nout, fmt = 9978 )
1047*
1048 DO 90 i = 1, nsubs
1049 WRITE( nout, fmt = 9980 ) '|', snames( i ), ktests( i ),
1050 $ kpass( i ), kfail( i ), kskip( i )
1051 90 CONTINUE
1052 WRITE( nout, fmt = * )
1053 WRITE( nout, fmt = 9976 )
1054 WRITE( nout, fmt = * )
1055*
1056 END IF
1057*
1058 CALL blacs_exit( 0 )
1059*
1060 9999 FORMAT( 'ILLEGAL ', a, ': ', a, ' = ', i10,
1061 $ ' should be at least 1' )
1062 9998 FORMAT( 'ILLEGAL GRID: NPROW*NPCOL = ', i4,
1063 $ '. It can be at most', i4 )
1064 9997 FORMAT( 'Bad ', a, ' parameters: going on to next test case.' )
1065 9996 FORMAT( 2x, 'Test number ', i4 , ' started on a ', i6, ' x ',
1066 $ i6, ' process grid.' )
1067 9995 FORMAT( 2x, ' ------------------------------------------------',
1068 $ '--------------------------' )
1069 9994 FORMAT( 2x, ' M N UPLO TRANS DIAG' )
1070 9993 FORMAT( 5x,i6,1x,i6,9x,a1,11x,a1,10x,a1 )
1071 9992 FORMAT( 2x, ' IA JA MA NA IMBA INBA',
1072 $ ' MBA NBA RSRCA CSRCA' )
1073 9991 FORMAT( 5x,i6,1x,i6,1x,i6,1x,i6,1x,i6,1x,i6,1x,i6,1x,i6,
1074 $ 1x,i5,1x,i5 )
1075 9990 FORMAT( 2x, ' IX JX MX NX IMBX INBX',
1076 $ ' MBX NBX RSRCX CSRCX INCX' )
1077 9989 FORMAT( 5x,i6,1x,i6,1x,i6,1x,i6,1x,i6,1x,i6,1x,i6,1x,i6,
1078 $ 1x,i5,1x,i5,1x,i6 )
1079 9988 FORMAT( 2x, ' IY JY MY NY IMBY INBY',
1080 $ ' MBY NBY RSRCY CSRCY INCY' )
1081 9987 FORMAT( 'Not enough memory for this test: going on to',
1082 $ ' next test case.' )
1083 9986 FORMAT( 'Not enough memory. Need: ', i12 )
1084 9985 FORMAT( 2x, ' Tested Subroutine: ', a )
1085 9984 FORMAT( 2x, ' ***** Computational check: ', a, ' ',
1086 $ ' FAILED ',' *****' )
1087 9983 FORMAT( 2x, ' ***** Computational check: ', a, ' ',
1088 $ ' PASSED ',' *****' )
1089 9982 FORMAT( 2x, ' ***** ERROR ***** Matrix operand ', a,
1090 $ ' modified by ', a, ' *****' )
1091 9981 FORMAT( 2x, 'Test number ', i4, ' completed.' )
1092 9980 FORMAT( 2x,a1,2x,a7,8x,i4,6x,i4,5x,i4,4x,i4 )
1093 9979 FORMAT( 2x, ' SUBROUTINE TOTAL TESTS PASSED FAILED ',
1094 $ 'SKIPPED' )
1095 9978 FORMAT( 2x, ' ---------- ----------- ------ ------ ',
1096 $ '-------' )
1097 9977 FORMAT( 2x, 'Testing Summary')
1098 9976 FORMAT( 2x, 'End of Tests.' )
1099 9975 FORMAT( 2x, 'Tests started.' )
1100 9974 FORMAT( 2x, ' ***** Operation not supported, error code: ',
1101 $ i5, ' *****' )
1102*
1103 stop
1104*
1105* End of PDBLA2TST
1106*
1107 END
1108 SUBROUTINE pdbla2tstinfo( SUMMRY, NOUT, NMAT, DIAGVAL, TRANVAL,
1109 $ UPLOVAL, MVAL, NVAL, MAVAL, NAVAL,
1110 $ IMBAVAL, MBAVAL, INBAVAL, NBAVAL,
1111 $ RSCAVAL, CSCAVAL, IAVAL, JAVAL,
1112 $ MXVAL, NXVAL, IMBXVAL, MBXVAL,
1113 $ INBXVAL, NBXVAL, RSCXVAL, CSCXVAL,
1114 $ IXVAL, JXVAL, INCXVAL, MYVAL, NYVAL,
1115 $ IMBYVAL, MBYVAL, INBYVAL, NBYVAL,
1116 $ RSCYVAL, CSCYVAL, IYVAL, JYVAL,
1117 $ INCYVAL, LDVAL, NGRIDS, PVAL, LDPVAL,
1118 $ QVAL, LDQVAL, NBLOG, LTEST, SOF, TEE,
1119 $ IAM, IGAP, IVERB, NPROCS, THRESH, ALPHA,
1120 $ BETA, WORK )
1121*
1122* -- PBLAS test routine (version 2.0) --
1123* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
1124* and University of California, Berkeley.
1125* April 1, 1998
1126*
1127* .. Scalar Arguments ..
1128 LOGICAL SOF, TEE
1129 INTEGER IAM, IGAP, IVERB, LDPVAL, LDQVAL, LDVAL, NBLOG,
1130 $ NGRIDS, NMAT, NOUT, NPROCS
1131 REAL THRESH
1132 DOUBLE PRECISION ALPHA, BETA
1133* ..
1134* .. Array Arguments ..
1135 CHARACTER*( * ) SUMMRY
1136 CHARACTER*1 DIAGVAL( LDVAL ), TRANVAL( LDVAL ),
1137 $ UPLOVAL( LDVAL )
1138 LOGICAL LTEST( * )
1139 INTEGER CSCAVAL( LDVAL ), CSCXVAL( LDVAL ),
1140 $ cscyval( ldval ), iaval( ldval ),
1141 $ imbaval( ldval ), imbxval( ldval ),
1142 $ imbyval( ldval ), inbaval( ldval ),
1143 $ inbxval( ldval ), inbyval( ldval ),
1144 $ incxval( ldval ), incyval( ldval ),
1145 $ ixval( ldval ), iyval( ldval ), javal( ldval ),
1146 $ jxval( ldval ), jyval( ldval ), maval( ldval ),
1147 $ mbaval( ldval ), mbxval( ldval ),
1148 $ mbyval( ldval ), mval( ldval ), mxval( ldval ),
1149 $ myval( ldval ), naval( ldval ),
1150 $ nbaval( ldval ), nbxval( ldval ),
1151 $ nbyval( ldval ), nval( ldval ), nxval( ldval ),
1152 $ nyval( ldval ), pval( ldpval ), qval( ldqval ),
1153 $ rscaval( ldval ), rscxval( ldval ),
1154 $ rscyval( ldval ), work( * )
1155* ..
1156*
1157* Purpose
1158* =======
1159*
1160* PDBLA2TSTINFO get the needed startup information for testing various
1161* Level 2 PBLAS routines, and transmits it to all processes.
1162*
1163* Notes
1164* =====
1165*
1166* For packing the information we assumed that the length in bytes of an
1167* integer is equal to the length in bytes of a real single precision.
1168*
1169* Arguments
1170* =========
1171*
1172* SUMMRY (global output) CHARACTER*(*)
1173* On exit, SUMMRY is the name of output (summary) file (if
1174* any). SUMMRY is only defined for process 0.
1175*
1176* NOUT (global output) INTEGER
1177* On exit, NOUT specifies the unit number for the output file.
1178* When NOUT is 6, output to screen, when NOUT is 0, output to
1179* stderr. NOUT is only defined for process 0.
1180*
1181* NMAT (global output) INTEGER
1182* On exit, NMAT specifies the number of different test cases.
1183*
1184* DIAGVAL (global output) CHARACTER array
1185* On entry, DIAGVAL is an array of dimension LDVAL. On exit,
1186* this array contains the values of DIAG to run the code with.
1187*
1188* TRANVAL (global output) CHARACTER array
1189* On entry, TRANVAL is an array of dimension LDVAL. On exit,
1190* this array contains the values of TRANS to run the code
1191* with.
1192*
1193* UPLOVAL (global output) CHARACTER array
1194* On entry, UPLOVAL is an array of dimension LDVAL. On exit,
1195* this array contains the values of UPLO to run the code with.
1196*
1197* MVAL (global output) INTEGER array
1198* On entry, MVAL is an array of dimension LDVAL. On exit, this
1199* array contains the values of M to run the code with.
1200*
1201* NVAL (global output) INTEGER array
1202* On entry, NVAL is an array of dimension LDVAL. On exit, this
1203* array contains the values of N to run the code with.
1204*
1205* MAVAL (global output) INTEGER array
1206* On entry, MAVAL is an array of dimension LDVAL. On exit, this
1207* array contains the values of DESCA( M_ ) to run the code
1208* with.
1209*
1210* NAVAL (global output) INTEGER array
1211* On entry, NAVAL is an array of dimension LDVAL. On exit, this
1212* array contains the values of DESCA( N_ ) to run the code
1213* with.
1214*
1215* IMBAVAL (global output) INTEGER array
1216* On entry, IMBAVAL is an array of dimension LDVAL. On exit,
1217* this array contains the values of DESCA( IMB_ ) to run the
1218* code with.
1219*
1220* MBAVAL (global output) INTEGER array
1221* On entry, MBAVAL is an array of dimension LDVAL. On exit,
1222* this array contains the values of DESCA( MB_ ) to run the
1223* code with.
1224*
1225* INBAVAL (global output) INTEGER array
1226* On entry, INBAVAL is an array of dimension LDVAL. On exit,
1227* this array contains the values of DESCA( INB_ ) to run the
1228* code with.
1229*
1230* NBAVAL (global output) INTEGER array
1231* On entry, NBAVAL is an array of dimension LDVAL. On exit,
1232* this array contains the values of DESCA( NB_ ) to run the
1233* code with.
1234*
1235* RSCAVAL (global output) INTEGER array
1236* On entry, RSCAVAL is an array of dimension LDVAL. On exit,
1237* this array contains the values of DESCA( RSRC_ ) to run the
1238* code with.
1239*
1240* CSCAVAL (global output) INTEGER array
1241* On entry, CSCAVAL is an array of dimension LDVAL. On exit,
1242* this array contains the values of DESCA( CSRC_ ) to run the
1243* code with.
1244*
1245* IAVAL (global output) INTEGER array
1246* On entry, IAVAL is an array of dimension LDVAL. On exit, this
1247* array contains the values of IA to run the code with.
1248*
1249* JAVAL (global output) INTEGER array
1250* On entry, JAVAL is an array of dimension LDVAL. On exit, this
1251* array contains the values of JA to run the code with.
1252*
1253* MXVAL (global output) INTEGER array
1254* On entry, MXVAL is an array of dimension LDVAL. On exit, this
1255* array contains the values of DESCX( M_ ) to run the code
1256* with.
1257*
1258* NXVAL (global output) INTEGER array
1259* On entry, NXVAL is an array of dimension LDVAL. On exit, this
1260* array contains the values of DESCX( N_ ) to run the code
1261* with.
1262*
1263* IMBXVAL (global output) INTEGER array
1264* On entry, IMBXVAL is an array of dimension LDVAL. On exit,
1265* this array contains the values of DESCX( IMB_ ) to run the
1266* code with.
1267*
1268* MBXVAL (global output) INTEGER array
1269* On entry, MBXVAL is an array of dimension LDVAL. On exit,
1270* this array contains the values of DESCX( MB_ ) to run the
1271* code with.
1272*
1273* INBXVAL (global output) INTEGER array
1274* On entry, INBXVAL is an array of dimension LDVAL. On exit,
1275* this array contains the values of DESCX( INB_ ) to run the
1276* code with.
1277*
1278* NBXVAL (global output) INTEGER array
1279* On entry, NBXVAL is an array of dimension LDVAL. On exit,
1280* this array contains the values of DESCX( NB_ ) to run the
1281* code with.
1282*
1283* RSCXVAL (global output) INTEGER array
1284* On entry, RSCXVAL is an array of dimension LDVAL. On exit,
1285* this array contains the values of DESCX( RSRC_ ) to run the
1286* code with.
1287*
1288* CSCXVAL (global output) INTEGER array
1289* On entry, CSCXVAL is an array of dimension LDVAL. On exit,
1290* this array contains the values of DESCX( CSRC_ ) to run the
1291* code with.
1292*
1293* IXVAL (global output) INTEGER array
1294* On entry, IXVAL is an array of dimension LDVAL. On exit, this
1295* array contains the values of IX to run the code with.
1296*
1297* JXVAL (global output) INTEGER array
1298* On entry, JXVAL is an array of dimension LDVAL. On exit, this
1299* array contains the values of JX to run the code with.
1300*
1301* INCXVAL (global output) INTEGER array
1302* On entry, INCXVAL is an array of dimension LDVAL. On exit,
1303* this array contains the values of INCX to run the code with.
1304*
1305* MYVAL (global output) INTEGER array
1306* On entry, MYVAL is an array of dimension LDVAL. On exit, this
1307* array contains the values of DESCY( M_ ) to run the code
1308* with.
1309*
1310* NYVAL (global output) INTEGER array
1311* On entry, NYVAL is an array of dimension LDVAL. On exit, this
1312* array contains the values of DESCY( N_ ) to run the code
1313* with.
1314*
1315* IMBYVAL (global output) INTEGER array
1316* On entry, IMBYVAL is an array of dimension LDVAL. On exit,
1317* this array contains the values of DESCY( IMB_ ) to run the
1318* code with.
1319*
1320* MBYVAL (global output) INTEGER array
1321* On entry, MBYVAL is an array of dimension LDVAL. On exit,
1322* this array contains the values of DESCY( MB_ ) to run the
1323* code with.
1324*
1325* INBYVAL (global output) INTEGER array
1326* On entry, INBYVAL is an array of dimension LDVAL. On exit,
1327* this array contains the values of DESCY( INB_ ) to run the
1328* code with.
1329*
1330* NBYVAL (global output) INTEGER array
1331* On entry, NBYVAL is an array of dimension LDVAL. On exit,
1332* this array contains the values of DESCY( NB_ ) to run the
1333* code with.
1334*
1335* RSCYVAL (global output) INTEGER array
1336* On entry, RSCYVAL is an array of dimension LDVAL. On exit,
1337* this array contains the values of DESCY( RSRC_ ) to run the
1338* code with.
1339*
1340* CSCYVAL (global output) INTEGER array
1341* On entry, CSCYVAL is an array of dimension LDVAL. On exit,
1342* this array contains the values of DESCY( CSRC_ ) to run the
1343* code with.
1344*
1345* IYVAL (global output) INTEGER array
1346* On entry, IYVAL is an array of dimension LDVAL. On exit, this
1347* array contains the values of IY to run the code with.
1348*
1349* JYVAL (global output) INTEGER array
1350* On entry, JYVAL is an array of dimension LDVAL. On exit, this
1351* array contains the values of JY to run the code with.
1352*
1353* INCYVAL (global output) INTEGER array
1354* On entry, INCYVAL is an array of dimension LDVAL. On exit,
1355* this array contains the values of INCY to run the code with.
1356*
1357* LDVAL (global input) INTEGER
1358* On entry, LDVAL specifies the maximum number of different va-
1359* lues that can be used for DIAG, TRANS, UPLO, M, N, DESCA(:),
1360* IA, JA, DESCX(:), IX, JX, INCX, DESCY(:), IY, JY and INCY.
1361* This is also the maximum number of test cases.
1362*
1363* NGRIDS (global output) INTEGER
1364* On exit, NGRIDS specifies the number of different values that
1365* can be used for P and Q.
1366*
1367* PVAL (global output) INTEGER array
1368* On entry, PVAL is an array of dimension LDPVAL. On exit, this
1369* array contains the values of P to run the code with.
1370*
1371* LDPVAL (global input) INTEGER
1372* On entry, LDPVAL specifies the maximum number of different
1373* values that can be used for P.
1374*
1375* QVAL (global output) INTEGER array
1376* On entry, QVAL is an array of dimension LDQVAL. On exit, this
1377* array contains the values of Q to run the code with.
1378*
1379* LDQVAL (global input) INTEGER
1380* On entry, LDQVAL specifies the maximum number of different
1381* values that can be used for Q.
1382*
1383* NBLOG (global output) INTEGER
1384* On exit, NBLOG specifies the logical computational block size
1385* to run the tests with. NBLOG must be at least one.
1386*
1387* LTEST (global output) LOGICAL array
1388* On entry, LTEST is an array of dimension at least seven. On
1389* exit, if LTEST( i ) is .TRUE., the i-th Level 2 PBLAS routine
1390* will be tested. See the input file for the ordering of the
1391* routines.
1392*
1393* SOF (global output) LOGICAL
1394* On exit, if SOF is .TRUE., the tester will stop on the first
1395* detected failure. Otherwise, it won't.
1396*
1397* TEE (global output) LOGICAL
1398* On exit, if TEE is .TRUE., the tester will perform the error
1399* exit tests. These tests won't be performed otherwise.
1400*
1401* IAM (local input) INTEGER
1402* On entry, IAM specifies the number of the process executing
1403* this routine.
1404*
1405* IGAP (global output) INTEGER
1406* On exit, IGAP specifies the user-specified gap used for pad-
1407* ding. IGAP must be at least zero.
1408*
1409* IVERB (global output) INTEGER
1410* On exit, IVERB specifies the output verbosity level: 0 for
1411* pass/fail, 1, 2 or 3 for matrix dump on errors.
1412*
1413* NPROCS (global input) INTEGER
1414* On entry, NPROCS specifies the total number of processes.
1415*
1416* THRESH (global output) REAL
1417* On exit, THRESH specifies the threshhold value for the test
1418* ratio.
1419*
1420* ALPHA (global output) DOUBLE PRECISION
1421* On exit, ALPHA specifies the value of alpha to be used in all
1422* the test cases.
1423*
1424* BETA (global output) DOUBLE PRECISION
1425* On exit, BETA specifies the value of beta to be used in all
1426* the test cases.
1427*
1428* WORK (local workspace) INTEGER array
1429* On entry, WORK is an array of dimension at least
1430* MAX( 3, 2*NGRIDS+37*NMAT+NSUBS+4 ) with NSUBS equal to 7.
1431* This array is used to pack all output arrays in order to send
1432* the information in one message.
1433*
1434* -- Written on April 1, 1998 by
1435* Antoine Petitet, University of Tennessee, Knoxville 37996, USA.
1436*
1437* =====================================================================
1438*
1439* .. Parameters ..
1440 INTEGER NIN, NSUBS
1441 PARAMETER ( NIN = 11, nsubs = 7 )
1442* ..
1443* .. Local Scalars ..
1444 LOGICAL LTESTT
1445 INTEGER I, ICTXT, J
1446 DOUBLE PRECISION EPS
1447* ..
1448* .. Local Arrays ..
1449 CHARACTER*7 SNAMET
1450 CHARACTER*79 USRINFO
1451* ..
1452* .. External Subroutines ..
1453 EXTERNAL blacs_abort, blacs_get, blacs_gridexit,
1454 $ blacs_gridinit, blacs_setup, dgebr2d, dgebs2d,
1455 $ icopy, igebr2d, igebs2d, sgebr2d, sgebs2d
1456*ype real dble cplx zplx
1457* ..
1458* .. External Functions ..
1459 DOUBLE PRECISION PDLAMCH
1460 EXTERNAL PDLAMCH
1461* ..
1462* .. Intrinsic Functions ..
1463 INTRINSIC char, ichar, max, min
1464* ..
1465* .. Common Blocks ..
1466 CHARACTER*7 SNAMES( NSUBS )
1467 COMMON /SNAMEC/SNAMES
1468* ..
1469* .. Executable Statements ..
1470*
1471* Process 0 reads the input data, broadcasts to other processes and
1472* writes needed information to NOUT
1473*
1474 IF( iam.EQ.0 ) THEN
1475*
1476* Open file and skip data file header
1477*
1478 OPEN( nin, file='PDBLAS2TST.dat', status='OLD' )
1479 READ( nin, fmt = * ) summry
1480 summry = ' '
1481*
1482* Read in user-supplied info about machine type, compiler, etc.
1483*
1484 READ( nin, fmt = 9999 ) usrinfo
1485*
1486* Read name and unit number for summary output file
1487*
1488 READ( nin, fmt = * ) summry
1489 READ( nin, fmt = * ) nout
1490 IF( nout.NE.0 .AND. nout.NE.6 )
1491 $ OPEN( nout, file = summry, status = 'UNKNOWN' )
1492*
1493* Read and check the parameter values for the tests.
1494*
1495* Read the flag that indicates if Stop on Failure
1496*
1497 READ( nin, fmt = * ) sof
1498*
1499* Read the flag that indicates if Test Error Exits
1500*
1501 READ( nin, fmt = * ) tee
1502*
1503* Read the verbosity level
1504*
1505 READ( nin, fmt = * ) iverb
1506 IF( iverb.LT.0 .OR. iverb.GT.3 )
1507 $ iverb = 0
1508*
1509* Read the leading dimension gap
1510*
1511 READ( nin, fmt = * ) igap
1512 IF( igap.LT.0 )
1513 $ igap = 0
1514*
1515* Read the threshold value for test ratio
1516*
1517 READ( nin, fmt = * ) thresh
1518 IF( thresh.LT.0.0 )
1519 $ thresh = 16.0
1520*
1521* Get logical computational block size
1522*
1523 READ( nin, fmt = * ) nblog
1524 IF( nblog.LT.1 )
1525 $ nblog = 32
1526*
1527* Get number of grids
1528*
1529 READ( nin, fmt = * ) ngrids
1530 IF( ngrids.LT.1 .OR. ngrids.GT.ldpval ) THEN
1531 WRITE( nout, fmt = 9998 ) 'Grids', ldpval
1532 GO TO 120
1533 ELSE IF( ngrids.GT.ldqval ) THEN
1534 WRITE( nout, fmt = 9998 ) 'Grids', ldqval
1535 GO TO 120
1536 END IF
1537*
1538* Get values of P and Q
1539*
1540 READ( nin, fmt = * ) ( pval( i ), i = 1, ngrids )
1541 READ( nin, fmt = * ) ( qval( i ), i = 1, ngrids )
1542*
1543* Read ALPHA, BETA
1544*
1545 READ( nin, fmt = * ) alpha
1546 READ( nin, fmt = * ) beta
1547*
1548* Read number of tests.
1549*
1550 READ( nin, fmt = * ) nmat
1551 IF( nmat.LT.1 .OR. nmat.GT.ldval ) THEN
1552 WRITE( nout, fmt = 9998 ) 'Tests', ldval
1553 GO TO 120
1554 ENDIF
1555*
1556* Read in input data into arrays.
1557*
1558 READ( nin, fmt = * ) ( uploval( i ), i = 1, nmat )
1559 READ( nin, fmt = * ) ( tranval( i ), i = 1, nmat )
1560 READ( nin, fmt = * ) ( diagval( i ), i = 1, nmat )
1561 READ( nin, fmt = * ) ( mval( i ), i = 1, nmat )
1562 READ( nin, fmt = * ) ( nval( i ), i = 1, nmat )
1563 READ( nin, fmt = * ) ( maval( i ), i = 1, nmat )
1564 READ( nin, fmt = * ) ( naval( i ), i = 1, nmat )
1565 READ( nin, fmt = * ) ( imbaval( i ), i = 1, nmat )
1566 READ( nin, fmt = * ) ( inbaval( i ), i = 1, nmat )
1567 READ( nin, fmt = * ) ( mbaval( i ), i = 1, nmat )
1568 READ( nin, fmt = * ) ( nbaval( i ), i = 1, nmat )
1569 READ( nin, fmt = * ) ( rscaval( i ), i = 1, nmat )
1570 READ( nin, fmt = * ) ( cscaval( i ), i = 1, nmat )
1571 READ( nin, fmt = * ) ( iaval( i ), i = 1, nmat )
1572 READ( nin, fmt = * ) ( javal( i ), i = 1, nmat )
1573 READ( nin, fmt = * ) ( mxval( i ), i = 1, nmat )
1574 READ( nin, fmt = * ) ( nxval( i ), i = 1, nmat )
1575 READ( nin, fmt = * ) ( imbxval( i ), i = 1, nmat )
1576 READ( nin, fmt = * ) ( inbxval( i ), i = 1, nmat )
1577 READ( nin, fmt = * ) ( mbxval( i ), i = 1, nmat )
1578 READ( nin, fmt = * ) ( nbxval( i ), i = 1, nmat )
1579 READ( nin, fmt = * ) ( rscxval( i ), i = 1, nmat )
1580 READ( nin, fmt = * ) ( cscxval( i ), i = 1, nmat )
1581 READ( nin, fmt = * ) ( ixval( i ), i = 1, nmat )
1582 READ( nin, fmt = * ) ( jxval( i ), i = 1, nmat )
1583 READ( nin, fmt = * ) ( incxval( i ), i = 1, nmat )
1584 READ( nin, fmt = * ) ( myval( i ), i = 1, nmat )
1585 READ( nin, fmt = * ) ( nyval( i ), i = 1, nmat )
1586 READ( nin, fmt = * ) ( imbyval( i ), i = 1, nmat )
1587 READ( nin, fmt = * ) ( inbyval( i ), i = 1, nmat )
1588 READ( nin, fmt = * ) ( mbyval( i ), i = 1, nmat )
1589 READ( nin, fmt = * ) ( nbyval( i ), i = 1, nmat )
1590 READ( nin, fmt = * ) ( rscyval( i ), i = 1, nmat )
1591 READ( nin, fmt = * ) ( cscyval( i ), i = 1, nmat )
1592 READ( nin, fmt = * ) ( iyval( i ), i = 1, nmat )
1593 READ( nin, fmt = * ) ( jyval( i ), i = 1, nmat )
1594 READ( nin, fmt = * ) ( incyval( i ), i = 1, nmat )
1595*
1596* Read names of subroutines and flags which indicate
1597* whether they are to be tested.
1598*
1599 DO 10 i = 1, nsubs
1600 ltest( i ) = .false.
1601 10 CONTINUE
1602 20 CONTINUE
1603 READ( nin, fmt = 9996, END = 50 ) SNAMET, ltestt
1604 DO 30 i = 1, nsubs
1605 IF( snamet.EQ.snames( i ) )
1606 $ GO TO 40
1607 30 CONTINUE
1608*
1609 WRITE( nout, fmt = 9995 )snamet
1610 GO TO 120
1611*
1612 40 CONTINUE
1613 ltest( i ) = ltestt
1614 GO TO 20
1615*
1616 50 CONTINUE
1617*
1618* Close input file
1619*
1620 CLOSE ( nin )
1621*
1622* For pvm only: if virtual machine not set up, allocate it and
1623* spawn the correct number of processes.
1624*
1625 IF( nprocs.LT.1 ) THEN
1626 nprocs = 0
1627 DO 60 i = 1, ngrids
1628 nprocs = max( nprocs, pval( i )*qval( i ) )
1629 60 CONTINUE
1630 CALL blacs_setup( iam, nprocs )
1631 END IF
1632*
1633* Temporarily define blacs grid to include all processes so
1634* information can be broadcast to all processes
1635*
1636 CALL blacs_get( -1, 0, ictxt )
1637 CALL blacs_gridinit( ictxt, 'Row-major', 1, nprocs )
1638*
1639* Compute machine epsilon
1640*
1641 eps = pdlamch( ictxt, 'eps' )
1642*
1643* Pack information arrays and broadcast
1644*
1645 CALL sgebs2d( ictxt, 'All', ' ', 1, 1, thresh, 1 )
1646 CALL dgebs2d( ictxt, 'All', ' ', 1, 1, alpha, 1 )
1647 CALL dgebs2d( ictxt, 'All', ' ', 1, 1, beta, 1 )
1648*
1649 work( 1 ) = ngrids
1650 work( 2 ) = nmat
1651 work( 3 ) = nblog
1652 CALL igebs2d( ictxt, 'All', ' ', 3, 1, work, 3 )
1653*
1654 i = 1
1655 IF( sof ) THEN
1656 work( i ) = 1
1657 ELSE
1658 work( i ) = 0
1659 END IF
1660 i = i + 1
1661 IF( tee ) THEN
1662 work( i ) = 1
1663 ELSE
1664 work( i ) = 0
1665 END IF
1666 i = i + 1
1667 work( i ) = iverb
1668 i = i + 1
1669 work( i ) = igap
1670 i = i + 1
1671 DO 70 j = 1, nmat
1672 work( i ) = ichar( diagval( j ) )
1673 work( i+1 ) = ichar( tranval( j ) )
1674 work( i+2 ) = ichar( uploval( j ) )
1675 i = i + 3
1676 70 CONTINUE
1677 CALL icopy( ngrids, pval, 1, work( i ), 1 )
1678 i = i + ngrids
1679 CALL icopy( ngrids, qval, 1, work( i ), 1 )
1680 i = i + ngrids
1681 CALL icopy( nmat, mval, 1, work( i ), 1 )
1682 i = i + nmat
1683 CALL icopy( nmat, nval, 1, work( i ), 1 )
1684 i = i + nmat
1685 CALL icopy( nmat, maval, 1, work( i ), 1 )
1686 i = i + nmat
1687 CALL icopy( nmat, naval, 1, work( i ), 1 )
1688 i = i + nmat
1689 CALL icopy( nmat, imbaval, 1, work( i ), 1 )
1690 i = i + nmat
1691 CALL icopy( nmat, inbaval, 1, work( i ), 1 )
1692 i = i + nmat
1693 CALL icopy( nmat, mbaval, 1, work( i ), 1 )
1694 i = i + nmat
1695 CALL icopy( nmat, nbaval, 1, work( i ), 1 )
1696 i = i + nmat
1697 CALL icopy( nmat, rscaval, 1, work( i ), 1 )
1698 i = i + nmat
1699 CALL icopy( nmat, cscaval, 1, work( i ), 1 )
1700 i = i + nmat
1701 CALL icopy( nmat, iaval, 1, work( i ), 1 )
1702 i = i + nmat
1703 CALL icopy( nmat, javal, 1, work( i ), 1 )
1704 i = i + nmat
1705 CALL icopy( nmat, mxval, 1, work( i ), 1 )
1706 i = i + nmat
1707 CALL icopy( nmat, nxval, 1, work( i ), 1 )
1708 i = i + nmat
1709 CALL icopy( nmat, imbxval, 1, work( i ), 1 )
1710 i = i + nmat
1711 CALL icopy( nmat, inbxval, 1, work( i ), 1 )
1712 i = i + nmat
1713 CALL icopy( nmat, mbxval, 1, work( i ), 1 )
1714 i = i + nmat
1715 CALL icopy( nmat, nbxval, 1, work( i ), 1 )
1716 i = i + nmat
1717 CALL icopy( nmat, rscxval, 1, work( i ), 1 )
1718 i = i + nmat
1719 CALL icopy( nmat, cscxval, 1, work( i ), 1 )
1720 i = i + nmat
1721 CALL icopy( nmat, ixval, 1, work( i ), 1 )
1722 i = i + nmat
1723 CALL icopy( nmat, jxval, 1, work( i ), 1 )
1724 i = i + nmat
1725 CALL icopy( nmat, incxval, 1, work( i ), 1 )
1726 i = i + nmat
1727 CALL icopy( nmat, myval, 1, work( i ), 1 )
1728 i = i + nmat
1729 CALL icopy( nmat, nyval, 1, work( i ), 1 )
1730 i = i + nmat
1731 CALL icopy( nmat, imbyval, 1, work( i ), 1 )
1732 i = i + nmat
1733 CALL icopy( nmat, inbyval, 1, work( i ), 1 )
1734 i = i + nmat
1735 CALL icopy( nmat, mbyval, 1, work( i ), 1 )
1736 i = i + nmat
1737 CALL icopy( nmat, nbyval, 1, work( i ), 1 )
1738 i = i + nmat
1739 CALL icopy( nmat, rscyval, 1, work( i ), 1 )
1740 i = i + nmat
1741 CALL icopy( nmat, cscyval, 1, work( i ), 1 )
1742 i = i + nmat
1743 CALL icopy( nmat, iyval, 1, work( i ), 1 )
1744 i = i + nmat
1745 CALL icopy( nmat, jyval, 1, work( i ), 1 )
1746 i = i + nmat
1747 CALL icopy( nmat, incyval, 1, work( i ), 1 )
1748 i = i + nmat
1749*
1750 DO 80 j = 1, nsubs
1751 IF( ltest( j ) ) THEN
1752 work( i ) = 1
1753 ELSE
1754 work( i ) = 0
1755 END IF
1756 i = i + 1
1757 80 CONTINUE
1758 i = i - 1
1759 CALL igebs2d( ictxt, 'All', ' ', i, 1, work, i )
1760*
1761* regurgitate input
1762*
1763 WRITE( nout, fmt = 9999 ) 'Level 2 PBLAS testing program.'
1764 WRITE( nout, fmt = 9999 ) usrinfo
1765 WRITE( nout, fmt = * )
1766 WRITE( nout, fmt = 9999 )
1767 $ 'Tests of the real double precision '//
1768 $ 'Level 2 PBLAS'
1769 WRITE( nout, fmt = * )
1770 WRITE( nout, fmt = 9993 ) nmat
1771 WRITE( nout, fmt = 9979 ) nblog
1772 WRITE( nout, fmt = 9992 ) ngrids
1773 WRITE( nout, fmt = 9990 )
1774 $ 'P', ( pval(i), i = 1, min(ngrids, 5) )
1775 IF( ngrids.GT.5 )
1776 $ WRITE( nout, fmt = 9991 ) ( pval(i), i = 6,
1777 $ min( 10, ngrids ) )
1778 IF( ngrids.GT.10 )
1779 $ WRITE( nout, fmt = 9991 ) ( pval(i), i = 11,
1780 $ min( 15, ngrids ) )
1781 IF( ngrids.GT.15 )
1782 $ WRITE( nout, fmt = 9991 ) ( pval(i), i = 16, ngrids )
1783 WRITE( nout, fmt = 9990 )
1784 $ 'Q', ( qval(i), i = 1, min(ngrids, 5) )
1785 IF( ngrids.GT.5 )
1786 $ WRITE( nout, fmt = 9991 ) ( qval(i), i = 6,
1787 $ min( 10, ngrids ) )
1788 IF( ngrids.GT.10 )
1789 $ WRITE( nout, fmt = 9991 ) ( qval(i), i = 11,
1790 $ min( 15, ngrids ) )
1791 IF( ngrids.GT.15 )
1792 $ WRITE( nout, fmt = 9991 ) ( qval(i), i = 16, ngrids )
1793 WRITE( nout, fmt = 9988 ) sof
1794 WRITE( nout, fmt = 9987 ) tee
1795 WRITE( nout, fmt = 9983 ) igap
1796 WRITE( nout, fmt = 9986 ) iverb
1797 WRITE( nout, fmt = 9980 ) thresh
1798 WRITE( nout, fmt = 9982 ) alpha
1799 WRITE( nout, fmt = 9981 ) beta
1800 IF( ltest( 1 ) ) THEN
1801 WRITE( nout, fmt = 9985 ) snames( 1 ), ' ... Yes'
1802 ELSE
1803 WRITE( nout, fmt = 9985 ) snames( 1 ), ' ... No '
1804 END IF
1805 DO 90 i = 2, nsubs
1806 IF( ltest( i ) ) THEN
1807 WRITE( nout, fmt = 9984 ) snames( i ), ' ... Yes'
1808 ELSE
1809 WRITE( nout, fmt = 9984 ) snames( i ), ' ... No '
1810 END IF
1811 90 CONTINUE
1812 WRITE( nout, fmt = 9994 ) eps
1813 WRITE( nout, fmt = * )
1814*
1815 ELSE
1816*
1817* If in pvm, must participate setting up virtual machine
1818*
1819 IF( nprocs.LT.1 )
1820 $ CALL blacs_setup( iam, nprocs )
1821*
1822* Temporarily define blacs grid to include all processes so
1823* information can be broadcast to all processes
1824*
1825 CALL blacs_get( -1, 0, ictxt )
1826 CALL blacs_gridinit( ictxt, 'Row-major', 1, nprocs )
1827*
1828* Compute machine epsilon
1829*
1830 eps = pdlamch( ictxt, 'eps' )
1831*
1832 CALL sgebr2d( ictxt, 'All', ' ', 1, 1, thresh, 1, 0, 0 )
1833 CALL dgebr2d( ictxt, 'All', ' ', 1, 1, alpha, 1, 0, 0 )
1834 CALL dgebr2d( ictxt, 'All', ' ', 1, 1, beta, 1, 0, 0 )
1835*
1836 CALL igebr2d( ictxt, 'All', ' ', 3, 1, work, 3, 0, 0 )
1837 ngrids = work( 1 )
1838 nmat = work( 2 )
1839 nblog = work( 3 )
1840*
1841 i = 2*ngrids + 37*nmat + nsubs + 4
1842 CALL igebr2d( ictxt, 'All', ' ', i, 1, work, i, 0, 0 )
1843*
1844 i = 1
1845 IF( work( i ).EQ.1 ) THEN
1846 sof = .true.
1847 ELSE
1848 sof = .false.
1849 END IF
1850 i = i + 1
1851 IF( work( i ).EQ.1 ) THEN
1852 tee = .true.
1853 ELSE
1854 tee = .false.
1855 END IF
1856 i = i + 1
1857 iverb = work( i )
1858 i = i + 1
1859 igap = work( i )
1860 i = i + 1
1861 DO 100 j = 1, nmat
1862 diagval( j ) = char( work( i ) )
1863 tranval( j ) = char( work( i+1 ) )
1864 uploval( j ) = char( work( i+2 ) )
1865 i = i + 3
1866 100 CONTINUE
1867 CALL icopy( ngrids, work( i ), 1, pval, 1 )
1868 i = i + ngrids
1869 CALL icopy( ngrids, work( i ), 1, qval, 1 )
1870 i = i + ngrids
1871 CALL icopy( nmat, work( i ), 1, mval, 1 )
1872 i = i + nmat
1873 CALL icopy( nmat, work( i ), 1, nval, 1 )
1874 i = i + nmat
1875 CALL icopy( nmat, work( i ), 1, maval, 1 )
1876 i = i + nmat
1877 CALL icopy( nmat, work( i ), 1, naval, 1 )
1878 i = i + nmat
1879 CALL icopy( nmat, work( i ), 1, imbaval, 1 )
1880 i = i + nmat
1881 CALL icopy( nmat, work( i ), 1, inbaval, 1 )
1882 i = i + nmat
1883 CALL icopy( nmat, work( i ), 1, mbaval, 1 )
1884 i = i + nmat
1885 CALL icopy( nmat, work( i ), 1, nbaval, 1 )
1886 i = i + nmat
1887 CALL icopy( nmat, work( i ), 1, rscaval, 1 )
1888 i = i + nmat
1889 CALL icopy( nmat, work( i ), 1, cscaval, 1 )
1890 i = i + nmat
1891 CALL icopy( nmat, work( i ), 1, iaval, 1 )
1892 i = i + nmat
1893 CALL icopy( nmat, work( i ), 1, javal, 1 )
1894 i = i + nmat
1895 CALL icopy( nmat, work( i ), 1, mxval, 1 )
1896 i = i + nmat
1897 CALL icopy( nmat, work( i ), 1, nxval, 1 )
1898 i = i + nmat
1899 CALL icopy( nmat, work( i ), 1, imbxval, 1 )
1900 i = i + nmat
1901 CALL icopy( nmat, work( i ), 1, inbxval, 1 )
1902 i = i + nmat
1903 CALL icopy( nmat, work( i ), 1, mbxval, 1 )
1904 i = i + nmat
1905 CALL icopy( nmat, work( i ), 1, nbxval, 1 )
1906 i = i + nmat
1907 CALL icopy( nmat, work( i ), 1, rscxval, 1 )
1908 i = i + nmat
1909 CALL icopy( nmat, work( i ), 1, cscxval, 1 )
1910 i = i + nmat
1911 CALL icopy( nmat, work( i ), 1, ixval, 1 )
1912 i = i + nmat
1913 CALL icopy( nmat, work( i ), 1, jxval, 1 )
1914 i = i + nmat
1915 CALL icopy( nmat, work( i ), 1, incxval, 1 )
1916 i = i + nmat
1917 CALL icopy( nmat, work( i ), 1, myval, 1 )
1918 i = i + nmat
1919 CALL icopy( nmat, work( i ), 1, nyval, 1 )
1920 i = i + nmat
1921 CALL icopy( nmat, work( i ), 1, imbyval, 1 )
1922 i = i + nmat
1923 CALL icopy( nmat, work( i ), 1, inbyval, 1 )
1924 i = i + nmat
1925 CALL icopy( nmat, work( i ), 1, mbyval, 1 )
1926 i = i + nmat
1927 CALL icopy( nmat, work( i ), 1, nbyval, 1 )
1928 i = i + nmat
1929 CALL icopy( nmat, work( i ), 1, rscyval, 1 )
1930 i = i + nmat
1931 CALL icopy( nmat, work( i ), 1, cscyval, 1 )
1932 i = i + nmat
1933 CALL icopy( nmat, work( i ), 1, iyval, 1 )
1934 i = i + nmat
1935 CALL icopy( nmat, work( i ), 1, jyval, 1 )
1936 i = i + nmat
1937 CALL icopy( nmat, work( i ), 1, incyval, 1 )
1938 i = i + nmat
1939*
1940 DO 110 j = 1, nsubs
1941 IF( work( i ).EQ.1 ) THEN
1942 ltest( j ) = .true.
1943 ELSE
1944 ltest( j ) = .false.
1945 END IF
1946 i = i + 1
1947 110 CONTINUE
1948*
1949 END IF
1950*
1951 CALL blacs_gridexit( ictxt )
1952*
1953 RETURN
1954*
1955 120 WRITE( nout, fmt = 9997 )
1956 CLOSE( nin )
1957 IF( nout.NE.6 .AND. nout.NE.0 )
1958 $ CLOSE( nout )
1959 CALL blacs_abort( ictxt, 1 )
1960*
1961 stop
1962*
1963 9999 FORMAT( a )
1964 9998 FORMAT( ' Number of values of ',5a, ' is less than 1 or greater ',
1965 $ 'than ', i2 )
1966 9997 FORMAT( ' Illegal input in file ',40a,'. Aborting run.' )
1967 9996 FORMAT( a7, l2 )
1968 9995 FORMAT( ' Subprogram name ', a7, ' not recognized',
1969 $ /' ******* TESTS ABANDONED *******' )
1970 9994 FORMAT( 2x, 'Relative machine precision (eps) is taken to be ',
1971 $ e18.6 )
1972 9993 FORMAT( 2x, 'Number of Tests : ', i6 )
1973 9992 FORMAT( 2x, 'Number of process grids : ', i6 )
1974 9991 FORMAT( 2x, ' : ', 5i6 )
1975 9990 FORMAT( 2x, a1, ' : ', 5i6 )
1976 9988 FORMAT( 2x, 'Stop on failure flag : ', l6 )
1977 9987 FORMAT( 2x, 'Test for error exits flag : ', l6 )
1978 9986 FORMAT( 2x, 'Verbosity level : ', i6 )
1979 9985 FORMAT( 2x, 'Routines to be tested : ', a, a8 )
1980 9984 FORMAT( 2x, ' ', a, a8 )
1981 9983 FORMAT( 2x, 'Leading dimension gap : ', i6 )
1982 9982 FORMAT( 2x, 'Alpha : ', g16.6 )
1983 9981 FORMAT( 2x, 'Beta : ', g16.6 )
1984 9980 FORMAT( 2x, 'Threshold value : ', g16.6 )
1985 9979 FORMAT( 2x, 'Logical block size : ', i6 )
1986*
1987* End of PDBLA2TSTINFO
1988*
1989 END
1990 SUBROUTINE pdblas2tstchke( LTEST, INOUT, NPROCS )
1991*
1992* -- PBLAS test routine (version 2.0) --
1993* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
1994* and University of California, Berkeley.
1995* April 1, 1998
1996*
1997* .. Scalar Arguments ..
1998 INTEGER INOUT, NPROCS
1999* ..
2000* .. Array Arguments ..
2001 LOGICAL LTEST( * )
2002* ..
2003*
2004* Purpose
2005* =======
2006*
2007* PDBLAS2TSTCHKE tests the error exits of the Level 2 PBLAS.
2008*
2009* Arguments
2010* =========
2011*
2012* LTEST (global input) LOGICAL array
2013* On entry, LTEST is an array of dimension at least 7 (NSUBS).
2014* If LTEST( 1 ) is .TRUE., PDGEMV will be tested;
2015* If LTEST( 2 ) is .TRUE., PDSYMV will be tested;
2016* If LTEST( 3 ) is .TRUE., PDTRMV will be tested;
2017* If LTEST( 4 ) is .TRUE., PDTRSV will be tested;
2018* If LTEST( 5 ) is .TRUE., PDGER will be tested;
2019* If LTEST( 6 ) is .TRUE., PDSYR will be tested;
2020* If LTEST( 7 ) is .TRUE., PDSYR2 will be tested;
2021*
2022* INOUT (global input) INTEGER
2023* On entry, INOUT specifies the unit number for output file.
2024* When INOUT is 6, output to screen, when INOUT = 0, output to
2025* stderr. INOUT is only defined in process 0.
2026*
2027* NPROCS (global input) INTEGER
2028* On entry, NPROCS specifies the total number of processes cal-
2029* ling this routine.
2030*
2031* Calling sequence encodings
2032* ==========================
2033*
2034* code Formal argument list Examples
2035*
2036* 11 (n, v1,v2) _SWAP, _COPY
2037* 12 (n,s1, v1 ) _SCAL, _SCAL
2038* 13 (n,s1, v1,v2) _AXPY, _DOT_
2039* 14 (n,s1,i1,v1 ) _AMAX
2040* 15 (n,u1, v1 ) _ASUM, _NRM2
2041*
2042* 21 ( trans, m,n,s1,m1,v1,s2,v2) _GEMV
2043* 22 (uplo, n,s1,m1,v1,s2,v2) _SYMV, _HEMV
2044* 23 (uplo,trans,diag, n, m1,v1 ) _TRMV, _TRSV
2045* 24 ( m,n,s1,v1,v2,m1) _GER_
2046* 25 (uplo, n,s1,v1, m1) _SYR
2047* 26 (uplo, n,u1,v1, m1) _HER
2048* 27 (uplo, n,s1,v1,v2,m1) _SYR2, _HER2
2049*
2050* 31 ( transa,transb, m,n,k,s1,m1,m2,s2,m3) _GEMM
2051* 32 (side,uplo, m,n, s1,m1,m2,s2,m3) _SYMM, _HEMM
2052* 33 ( uplo,trans, n,k,s1,m1, s2,m3) _SYRK
2053* 34 ( uplo,trans, n,k,u1,m1, u2,m3) _HERK
2054* 35 ( uplo,trans, n,k,s1,m1,m2,s2,m3) _SYR2K
2055* 36 ( uplo,trans, n,k,s1,m1,m2,u2,m3) _HER2K
2056* 37 ( m,n, s1,m1, s2,m3) _TRAN_
2057* 38 (side,uplo,transa, diag,m,n, s1,m1,m2 ) _TRMM, _TRSM
2058* 39 ( trans, m,n, s1,m1, s2,m3) _GEADD
2059* 40 ( uplo,trans, m,n, s1,m1, s2,m3) _TRADD
2060*
2061* -- Written on April 1, 1998 by
2062* Antoine Petitet, University of Tennessee, Knoxville 37996, USA.
2063*
2064* =====================================================================
2065*
2066* .. Parameters ..
2067 INTEGER NSUBS
2068 PARAMETER ( NSUBS = 7 )
2069* ..
2070* .. Local Scalars ..
2071 logical abrtsav
2072 INTEGER I, ICTXT, MYCOL, MYROW, NPCOL, NPROW
2073* ..
2074* .. Local Arrays ..
2075 INTEGER SCODE( NSUBS )
2076* ..
2077* .. External Subroutines ..
2078 EXTERNAL BLACS_GET, BLACS_GRIDEXIT, BLACS_GRIDINFO,
2079 $ blacs_gridinit, pddimee, pdgemv, pdger,
2080 $ pdmatee, pdoptee, pdsymv, pdsyr, pdsyr2,
2081 $ pdtrmv, pdtrsv, pdvecee
2082* ..
2083* .. Common Blocks ..
2084 LOGICAL ABRTFLG
2085 INTEGER NOUT
2086 CHARACTER*7 SNAMES( NSUBS )
2087 COMMON /snamec/snames
2088 COMMON /pberrorc/nout, abrtflg
2089* ..
2090* .. Data Statements ..
2091 DATA scode/21, 22, 23, 23, 24, 25, 27/
2092* ..
2093* .. Executable Statements ..
2094*
2095* Temporarily define blacs grid to include all processes so
2096* information can be broadcast to all processes.
2097*
2098 CALL blacs_get( -1, 0, ictxt )
2099 CALL blacs_gridinit( ictxt, 'Row-major', 1, nprocs )
2100 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
2101*
2102* Set ABRTFLG to FALSE so that the PBLAS error handler won't abort
2103* on errors during these tests and set the output device unit for
2104* it.
2105*
2106 abrtsav = abrtflg
2107 abrtflg = .false.
2108 nout = inout
2109*
2110* Test PDGEMV
2111*
2112 i = 1
2113 IF( ltest( i ) ) THEN
2114 CALL pdoptee( ictxt, nout, pdgemv, scode( i ), snames( i ) )
2115 CALL pddimee( ictxt, nout, pdgemv, scode( i ), snames( i ) )
2116 CALL pdmatee( ictxt, nout, pdgemv, scode( i ), snames( i ) )
2117 CALL pdvecee( ictxt, nout, pdgemv, scode( i ), snames( i ) )
2118 END IF
2119*
2120* Test PDSYMV
2121*
2122 i = i + 1
2123 IF( ltest( i ) ) THEN
2124 CALL pdoptee( ictxt, nout, pdsymv, scode( i ), snames( i ) )
2125 CALL pddimee( ictxt, nout, pdsymv, scode( i ), snames( i ) )
2126 CALL pdmatee( ictxt, nout, pdsymv, scode( i ), snames( i ) )
2127 CALL pdvecee( ictxt, nout, pdsymv, scode( i ), snames( i ) )
2128 END IF
2129*
2130* Test PDTRMV
2131*
2132 i = i + 1
2133 IF( ltest( i ) ) THEN
2134 CALL pdoptee( ictxt, nout, pdtrmv, scode( i ), snames( i ) )
2135 CALL pddimee( ictxt, nout, pdtrmv, scode( i ), snames( i ) )
2136 CALL pdmatee( ictxt, nout, pdtrmv, scode( i ), snames( i ) )
2137 CALL pdvecee( ictxt, nout, pdtrmv, scode( i ), snames( i ) )
2138 END IF
2139*
2140* Test PDTRSV
2141*
2142 i = i + 1
2143 IF( ltest( i ) ) THEN
2144 CALL pdoptee( ictxt, nout, pdtrsv, scode( i ), snames( i ) )
2145 CALL pddimee( ictxt, nout, pdtrsv, scode( i ), snames( i ) )
2146 CALL pdmatee( ictxt, nout, pdtrsv, scode( i ), snames( i ) )
2147 CALL pdvecee( ictxt, nout, pdtrsv, scode( i ), snames( i ) )
2148 END IF
2149*
2150* Test PDGER
2151*
2152 i = i + 1
2153 IF( ltest( i ) ) THEN
2154 CALL pddimee( ictxt, nout, pdger, scode( i ), snames( i ) )
2155 CALL pdvecee( ictxt, nout, pdger, scode( i ), snames( i ) )
2156 CALL pdmatee( ictxt, nout, pdger, scode( i ), snames( i ) )
2157 END IF
2158*
2159* Test PDSYR
2160*
2161 i = i + 1
2162 IF( ltest( i ) ) THEN
2163 CALL pdoptee( ictxt, nout, pdsyr, scode( i ), snames( i ) )
2164 CALL pddimee( ictxt, nout, pdsyr, scode( i ), snames( i ) )
2165 CALL pdvecee( ictxt, nout, pdsyr, scode( i ), snames( i ) )
2166 CALL pdmatee( ictxt, nout, pdsyr, scode( i ), snames( i ) )
2167 END IF
2168*
2169* Test PDSYR2
2170*
2171 i = i + 1
2172 IF( ltest( i ) ) THEN
2173 CALL pdoptee( ictxt, nout, pdsyr2, scode( i ), snames( i ) )
2174 CALL pddimee( ictxt, nout, pdsyr2, scode( i ), snames( i ) )
2175 CALL pdvecee( ictxt, nout, pdsyr2, scode( i ), snames( i ) )
2176 CALL pdmatee( ictxt, nout, pdsyr2, scode( i ), snames( i ) )
2177 END IF
2178*
2179 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2180 $ WRITE( nout, fmt = 9999 )
2181*
2182 CALL blacs_gridexit( ictxt )
2183*
2184* Reset ABRTFLG to the value it had before calling this routine
2185*
2186 abrtflg = abrtsav
2187*
2188 9999 FORMAT( 2x, 'Error-exit tests completed.' )
2189*
2190 RETURN
2191*
2192* End of PDBLAS2TSTCHKE
2193*
2194 END
2195 SUBROUTINE pdchkarg2( ICTXT, NOUT, SNAME, UPLO, TRANS, DIAG, M,
2196 $ N, ALPHA, IA, JA, DESCA, IX, JX, DESCX,
2197 $ INCX, BETA, IY, JY, DESCY, INCY, INFO )
2198*
2199* -- PBLAS test routine (version 2.0) --
2200* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
2201* and University of California, Berkeley.
2202* April 1, 1998
2203*
2204* .. Scalar Arguments ..
2205 CHARACTER*1 DIAG, TRANS, UPLO
2206 INTEGER IA, ICTXT, INCX, INCY, INFO, IX, IY, JA, JX,
2207 $ JY, M, N, NOUT
2208 DOUBLE PRECISION ALPHA, BETA
2209* ..
2210* .. Array Arguments ..
2211 CHARACTER*(*) SNAME
2212 INTEGER DESCA( * ), DESCX( * ), DESCY( * )
2213* ..
2214*
2215* Purpose
2216* =======
2217*
2218* PDCHKARG2 checks the input-only arguments of the Level 2 PBLAS. When
2219* INFO = 0, this routine makes a copy of its arguments (which are INPUT
2220* only arguments to PBLAS routines). Otherwise, it verifies the values
2221* of these arguments against the saved copies.
2222*
2223* Arguments
2224* =========
2225*
2226* ICTXT (local input) INTEGER
2227* On entry, ICTXT specifies the BLACS context handle, indica-
2228* ting the global context of the operation. The context itself
2229* is global, but the value of ICTXT is local.
2230*
2231* NOUT (global input) INTEGER
2232* On entry, NOUT specifies the unit number for the output file.
2233* When NOUT is 6, output to screen, when NOUT is 0, output to
2234* stderr. NOUT is only defined for process 0.
2235*
2236* SNAME (global input) CHARACTER*(*)
2237* On entry, SNAME specifies the subroutine name calling this
2238* subprogram.
2239*
2240* UPLO (global input) CHARACTER*1
2241* On entry, UPLO specifies the UPLO option in the Level 2 PBLAS
2242* operation.
2243*
2244* TRANS (global input) CHARACTER*1
2245* On entry, TRANS specifies the TRANS option in the Level 2
2246* PBLAS operation.
2247*
2248* DIAG (global input) CHARACTER*1
2249* On entry, DIAG specifies the DIAG option in the Level 2 PBLAS
2250* operation.
2251*
2252* M (global input) INTEGER
2253* On entry, M specifies the dimension of the submatrix ope-
2254* rands.
2255*
2256* N (global input) INTEGER
2257* On entry, N specifies the dimension of the submatrix ope-
2258* rands.
2259*
2260* ALPHA (global input) DOUBLE PRECISION
2261* On entry, ALPHA specifies the scalar alpha.
2262*
2263* IA (global input) INTEGER
2264* On entry, IA specifies A's global row index, which points to
2265* the beginning of the submatrix sub( A ).
2266*
2267* JA (global input) INTEGER
2268* On entry, JA specifies A's global column index, which points
2269* to the beginning of the submatrix sub( A ).
2270*
2271* DESCA (global and local input) INTEGER array
2272* On entry, DESCA is an integer array of dimension DLEN_. This
2273* is the array descriptor for the matrix A.
2274*
2275* IX (global input) INTEGER
2276* On entry, IX specifies X's global row index, which points to
2277* the beginning of the submatrix sub( X ).
2278*
2279* JX (global input) INTEGER
2280* On entry, JX specifies X's global column index, which points
2281* to the beginning of the submatrix sub( X ).
2282*
2283* DESCX (global and local input) INTEGER array
2284* On entry, DESCX is an integer array of dimension DLEN_. This
2285* is the array descriptor for the matrix X.
2286*
2287* INCX (global input) INTEGER
2288* On entry, INCX specifies the global increment for the
2289* elements of X. Only two values of INCX are supported in
2290* this version, namely 1 and M_X. INCX must not be zero.
2291*
2292* BETA (global input) DOUBLE PRECISION
2293* On entry, BETA specifies the scalar beta.
2294*
2295* IY (global input) INTEGER
2296* On entry, IY specifies Y's global row index, which points to
2297* the beginning of the submatrix sub( Y ).
2298*
2299* JY (global input) INTEGER
2300* On entry, JY specifies Y's global column index, which points
2301* to the beginning of the submatrix sub( Y ).
2302*
2303* DESCY (global and local input) INTEGER array
2304* On entry, DESCY is an integer array of dimension DLEN_. This
2305* is the array descriptor for the matrix Y.
2306*
2307* INCY (global input) INTEGER
2308* On entry, INCY specifies the global increment for the
2309* elements of Y. Only two values of INCY are supported in
2310* this version, namely 1 and M_Y. INCY must not be zero.
2311*
2312* INFO (global input/global output) INTEGER
2313* When INFO = 0 on entry, the values of the arguments which are
2314* INPUT only arguments to a PBLAS routine are copied into sta-
2315* tic variables and INFO is unchanged on exit. Otherwise, the
2316* values of the arguments are compared against the saved co-
2317* pies. In case no error has been found INFO is zero on return,
2318* otherwise it is non zero.
2319*
2320* -- Written on April 1, 1998 by
2321* Antoine Petitet, University of Tennessee, Knoxville 37996, USA.
2322*
2323* =====================================================================
2324*
2325* .. Parameters ..
2326 INTEGER BLOCK_CYCLIC_2D_INB, CSRC_, CTXT_, DLEN_,
2327 $ DTYPE_, IMB_, INB_, LLD_, MB_, M_, NB_, N_,
2328 $ RSRC_
2329 PARAMETER ( BLOCK_CYCLIC_2D_INB = 2, dlen_ = 11,
2330 $ dtype_ = 1, ctxt_ = 2, m_ = 3, n_ = 4,
2331 $ imb_ = 5, inb_ = 6, mb_ = 7, nb_ = 8,
2332 $ rsrc_ = 9, csrc_ = 10, lld_ = 11 )
2333* ..
2334* .. Local Scalars ..
2335 CHARACTER*1 DIAGREF, TRANSREF, UPLOREF
2336 INTEGER I, IAREF, INCXREF, INCYREF, IXREF, IYREF,
2337 $ JAREF, JXREF, JYREF, MREF, MYCOL, MYROW, NPCOL,
2338 $ NPROW, NREF
2339 DOUBLE PRECISION ALPHAREF, BETAREF
2340* ..
2341* .. Local Arrays ..
2342 CHARACTER*15 ARGNAME
2343 INTEGER DESCAREF( DLEN_ ), DESCXREF( DLEN_ ),
2344 $ DESCYREF( DLEN_ )
2345* ..
2346* .. External Subroutines ..
2347 EXTERNAL BLACS_GRIDINFO, IGSUM2D
2348* ..
2349* .. External Functions ..
2350 LOGICAL LSAME
2351 EXTERNAL LSAME
2352* ..
2353* .. Save Statements ..
2354 SAVE
2355* ..
2356* .. Executable Statements ..
2357*
2358* Get grid parameters
2359*
2360 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
2361*
2362* Check if first call. If yes, then save.
2363*
2364 IF( info.EQ.0 ) THEN
2365*
2366 diagref = diag
2367 transref = trans
2368 uploref = uplo
2369 mref = m
2370 nref = n
2371 alpharef = alpha
2372 iaref = ia
2373 jaref = ja
2374 DO 10 i = 1, dlen_
2375 descaref( i ) = desca( i )
2376 10 CONTINUE
2377 ixref = ix
2378 jxref = jx
2379 DO 20 i = 1, dlen_
2380 descxref( i ) = descx( i )
2381 20 CONTINUE
2382 incxref = incx
2383 betaref = beta
2384 iyref = iy
2385 jyref = jy
2386 DO 30 i = 1, dlen_
2387 descyref( i ) = descy( i )
2388 30 CONTINUE
2389 incyref = incy
2390*
2391 ELSE
2392*
2393* Test saved args. Return with first mismatch.
2394*
2395 argname = ' '
2396 IF( .NOT. lsame( diag, diagref ) ) THEN
2397 WRITE( argname, fmt = '(A)' ) 'DIAG'
2398 ELSE IF( .NOT. lsame( trans, transref ) ) THEN
2399 WRITE( argname, fmt = '(A)' ) 'TRANS'
2400 ELSE IF( .NOT. lsame( uplo, uploref ) ) THEN
2401 WRITE( argname, fmt = '(A)' ) 'UPLO'
2402 ELSE IF( m.NE.mref ) THEN
2403 WRITE( argname, fmt = '(A)' ) 'M'
2404 ELSE IF( n.NE.nref ) THEN
2405 WRITE( argname, fmt = '(A)' ) 'N'
2406 ELSE IF( alpha.NE.alpharef ) THEN
2407 WRITE( argname, fmt = '(A)' ) 'ALPHA'
2408 ELSE IF( ia.NE.iaref ) THEN
2409 WRITE( argname, fmt = '(A)' ) 'IA'
2410 ELSE IF( ja.NE.jaref ) THEN
2411 WRITE( argname, fmt = '(A)' ) 'JA'
2412 ELSE IF( desca( dtype_ ).NE.descaref( dtype_ ) ) THEN
2413 WRITE( argname, fmt = '(A)' ) 'DESCA( DTYPE_ )'
2414 ELSE IF( desca( m_ ).NE.descaref( m_ ) ) THEN
2415 WRITE( argname, fmt = '(A)' ) 'DESCA( M_ )'
2416 ELSE IF( desca( n_ ).NE.descaref( n_ ) ) THEN
2417 WRITE( argname, fmt = '(A)' ) 'DESCA( N_ )'
2418 ELSE IF( desca( imb_ ).NE.descaref( imb_ ) ) THEN
2419 WRITE( argname, fmt = '(A)' ) 'DESCA( IMB_ )'
2420 ELSE IF( desca( inb_ ).NE.descaref( inb_ ) ) THEN
2421 WRITE( argname, fmt = '(A)' ) 'DESCA( INB_ )'
2422 ELSE IF( desca( mb_ ).NE.descaref( mb_ ) ) THEN
2423 WRITE( argname, fmt = '(A)' ) 'DESCA( MB_ )'
2424 ELSE IF( desca( nb_ ).NE.descaref( nb_ ) ) THEN
2425 WRITE( argname, fmt = '(A)' ) 'DESCA( NB_ )'
2426 ELSE IF( desca( rsrc_ ).NE.descaref( rsrc_ ) ) THEN
2427 WRITE( argname, fmt = '(A)' ) 'DESCA( RSRC_ )'
2428 ELSE IF( desca( csrc_ ).NE.descaref( csrc_ ) ) THEN
2429 WRITE( argname, fmt = '(A)' ) 'DESCA( CSRC_ )'
2430 ELSE IF( desca( ctxt_ ).NE.descaref( ctxt_ ) ) THEN
2431 WRITE( argname, fmt = '(A)' ) 'DESCA( CTXT_ )'
2432 ELSE IF( desca( lld_ ).NE.descaref( lld_ ) ) THEN
2433 WRITE( argname, fmt = '(A)' ) 'DESCA( LLD_ )'
2434 ELSE IF( ix.NE.ixref ) THEN
2435 WRITE( argname, fmt = '(A)' ) 'IX'
2436 ELSE IF( jx.NE.jxref ) THEN
2437 WRITE( argname, fmt = '(A)' ) 'JX'
2438 ELSE IF( descx( dtype_ ).NE.descxref( dtype_ ) ) THEN
2439 WRITE( argname, fmt = '(A)' ) 'DESCX( DTYPE_ )'
2440 ELSE IF( descx( m_ ).NE.descxref( m_ ) ) THEN
2441 WRITE( argname, fmt = '(A)' ) 'DESCX( M_ )'
2442 ELSE IF( descx( n_ ).NE.descxref( n_ ) ) THEN
2443 WRITE( argname, fmt = '(A)' ) 'DESCX( N_ )'
2444 ELSE IF( descx( imb_ ).NE.descxref( imb_ ) ) THEN
2445 WRITE( argname, fmt = '(A)' ) 'DESCX( IMB_ )'
2446 ELSE IF( descx( inb_ ).NE.descxref( inb_ ) ) THEN
2447 WRITE( argname, fmt = '(A)' ) 'DESCX( INB_ )'
2448 ELSE IF( descx( mb_ ).NE.descxref( mb_ ) ) THEN
2449 WRITE( argname, fmt = '(A)' ) 'DESCX( MB_ )'
2450 ELSE IF( descx( nb_ ).NE.descxref( nb_ ) ) THEN
2451 WRITE( argname, fmt = '(A)' ) 'DESCX( NB_ )'
2452 ELSE IF( descx( rsrc_ ).NE.descxref( rsrc_ ) ) THEN
2453 WRITE( argname, fmt = '(A)' ) 'DESCX( RSRC_ )'
2454 ELSE IF( descx( csrc_ ).NE.descxref( csrc_ ) ) THEN
2455 WRITE( argname, fmt = '(A)' ) 'DESCX( CSRC_ )'
2456 ELSE IF( descx( ctxt_ ).NE.descxref( ctxt_ ) ) THEN
2457 WRITE( argname, fmt = '(A)' ) 'DESCX( CTXT_ )'
2458 ELSE IF( descx( lld_ ).NE.descxref( lld_ ) ) THEN
2459 WRITE( argname, fmt = '(A)' ) 'DESCX( LLD_ )'
2460 ELSE IF( incx.NE.incxref ) THEN
2461 WRITE( argname, fmt = '(A)' ) 'INCX'
2462 ELSE IF( beta.NE.betaref ) THEN
2463 WRITE( argname, fmt = '(A)' ) 'BETA'
2464 ELSE IF( iy.NE.iyref ) THEN
2465 WRITE( argname, fmt = '(A)' ) 'IY'
2466 ELSE IF( jy.NE.jyref ) THEN
2467 WRITE( argname, fmt = '(A)' ) 'JY'
2468 ELSE IF( descy( dtype_ ).NE.descyref( dtype_ ) ) THEN
2469 WRITE( argname, fmt = '(A)' ) 'DESCY( DTYPE_ )'
2470 ELSE IF( descy( m_ ).NE.descyref( m_ ) ) THEN
2471 WRITE( argname, fmt = '(A)' ) 'DESCY( M_ )'
2472 ELSE IF( descy( n_ ).NE.descyref( n_ ) ) THEN
2473 WRITE( argname, fmt = '(A)' ) 'DESCY( N_ )'
2474 ELSE IF( descy( imb_ ).NE.descyref( imb_ ) ) THEN
2475 WRITE( argname, fmt = '(A)' ) 'DESCY( IMB_ )'
2476 ELSE IF( descy( inb_ ).NE.descyref( inb_ ) ) THEN
2477 WRITE( argname, fmt = '(A)' ) 'DESCY( INB_ )'
2478 ELSE IF( descy( mb_ ).NE.descyref( mb_ ) ) THEN
2479 WRITE( argname, fmt = '(A)' ) 'DESCY( MB_ )'
2480 ELSE IF( descy( nb_ ).NE.descyref( nb_ ) ) THEN
2481 WRITE( argname, fmt = '(A)' ) 'DESCY( NB_ )'
2482 ELSE IF( descy( rsrc_ ).NE.descyref( rsrc_ ) ) THEN
2483 WRITE( argname, fmt = '(A)' ) 'DESCY( RSRC_ )'
2484 ELSE IF( descy( csrc_ ).NE.descyref( csrc_ ) ) THEN
2485 WRITE( argname, fmt = '(A)' ) 'DESCY( CSRC_ )'
2486 ELSE IF( descy( ctxt_ ).NE.descyref( ctxt_ ) ) THEN
2487 WRITE( argname, fmt = '(A)' ) 'DESCY( CTXT_ )'
2488 ELSE IF( descy( lld_ ).NE.descyref( lld_ ) ) THEN
2489 WRITE( argname, fmt = '(A)' ) 'DESCY( LLD_ )'
2490 ELSE IF( incy.NE.incyref ) THEN
2491 WRITE( argname, fmt = '(A)' ) 'INCY'
2492 ELSE
2493 info = 0
2494 END IF
2495*
2496 CALL igsum2d( ictxt, 'All', ' ', 1, 1, info, 1, -1, 0 )
2497*
2498 IF( myrow.EQ.0 .AND. mycol.EQ.0 ) THEN
2499*
2500 IF( info.NE.0 ) THEN
2501 WRITE( nout, fmt = 9999 ) argname, sname
2502 ELSE
2503 WRITE( nout, fmt = 9998 ) sname
2504 END IF
2505*
2506 END IF
2507*
2508 END IF
2509*
2510 9999 FORMAT( 2x, ' ***** Input-only parameter check: ', a,
2511 $ ' FAILED changed ', a, ' *****' )
2512 9998 FORMAT( 2x, ' ***** Input-only parameter check: ', a,
2513 $ ' PASSED *****' )
2514*
2515 RETURN
2516*
2517* End of PDCHKARG2
2518*
2519 END
2520 SUBROUTINE pdblas2tstchk( ICTXT, NOUT, NROUT, UPLO, TRANS, DIAG,
2521 $ M, N, ALPHA, A, PA, IA, JA, DESCA, X,
2522 $ PX, IX, JX, DESCX, INCX, BETA, Y, PY,
2523 $ IY, JY, DESCY, INCY, THRESH, ROGUE,
2524 $ WORK, INFO )
2525*
2526* -- PBLAS test routine (version 2.0) --
2527* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
2528* and University of California, Berkeley.
2529* April 1, 1998
2530*
2531* .. Scalar Arguments ..
2532 CHARACTER*1 DIAG, TRANS, UPLO
2533 INTEGER IA, ICTXT, INCX, INCY, INFO, IX, IY, JA, JX,
2534 $ JY, M, N, NOUT, NROUT
2535 REAL THRESH
2536 DOUBLE PRECISION ALPHA, BETA, ROGUE
2537* ..
2538* .. Array Arguments ..
2539 INTEGER DESCA( * ), DESCX( * ), DESCY( * )
2540 DOUBLE PRECISION A( * ), PA( * ), PX( * ), PY( * ), WORK( * ),
2541 $ X( * ), Y( * )
2542* ..
2543*
2544* Purpose
2545* =======
2546*
2547* PDBLAS2TSTCHK performs the computational tests of the Level 2 PBLAS.
2548*
2549* Notes
2550* =====
2551*
2552* A description vector is associated with each 2D block-cyclicly dis-
2553* tributed matrix. This vector stores the information required to
2554* establish the mapping between a matrix entry and its corresponding
2555* process and memory location.
2556*
2557* In the following comments, the character _ should be read as
2558* "of the distributed matrix". Let A be a generic term for any 2D
2559* block cyclicly distributed matrix. Its description vector is DESCA:
2560*
2561* NOTATION STORED IN EXPLANATION
2562* ---------------- --------------- ------------------------------------
2563* DTYPE_A (global) DESCA( DTYPE_ ) The descriptor type.
2564* CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
2565* the NPROW x NPCOL BLACS process grid
2566* A is distributed over. The context
2567* itself is global, but the handle
2568* (the integer value) may vary.
2569* M_A (global) DESCA( M_ ) The number of rows in the distribu-
2570* ted matrix A, M_A >= 0.
2571* N_A (global) DESCA( N_ ) The number of columns in the distri-
2572* buted matrix A, N_A >= 0.
2573* IMB_A (global) DESCA( IMB_ ) The number of rows of the upper left
2574* block of the matrix A, IMB_A > 0.
2575* INB_A (global) DESCA( INB_ ) The number of columns of the upper
2576* left block of the matrix A,
2577* INB_A > 0.
2578* MB_A (global) DESCA( MB_ ) The blocking factor used to distri-
2579* bute the last M_A-IMB_A rows of A,
2580* MB_A > 0.
2581* NB_A (global) DESCA( NB_ ) The blocking factor used to distri-
2582* bute the last N_A-INB_A columns of
2583* A, NB_A > 0.
2584* RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
2585* row of the matrix A is distributed,
2586* NPROW > RSRC_A >= 0.
2587* CSRC_A (global) DESCA( CSRC_ ) The process column over which the
2588* first column of A is distributed.
2589* NPCOL > CSRC_A >= 0.
2590* LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
2591* array storing the local blocks of
2592* the distributed matrix A,
2593* IF( Lc( 1, N_A ) > 0 )
2594* LLD_A >= MAX( 1, Lr( 1, M_A ) )
2595* ELSE
2596* LLD_A >= 1.
2597*
2598* Let K be the number of rows of a matrix A starting at the global in-
2599* dex IA,i.e, A( IA:IA+K-1, : ). Lr( IA, K ) denotes the number of rows
2600* that the process of row coordinate MYROW ( 0 <= MYROW < NPROW ) would
2601* receive if these K rows were distributed over NPROW processes. If K
2602* is the number of columns of a matrix A starting at the global index
2603* JA, i.e, A( :, JA:JA+K-1, : ), Lc( JA, K ) denotes the number of co-
2604* lumns that the process MYCOL ( 0 <= MYCOL < NPCOL ) would receive if
2605* these K columns were distributed over NPCOL processes.
2606*
2607* The values of Lr() and Lc() may be determined via a call to the func-
2608* tion PB_NUMROC:
2609* Lr( IA, K ) = PB_NUMROC( K, IA, IMB_A, MB_A, MYROW, RSRC_A, NPROW )
2610* Lc( JA, K ) = PB_NUMROC( K, JA, INB_A, NB_A, MYCOL, CSRC_A, NPCOL )
2611*
2612* Arguments
2613* =========
2614*
2615* ICTXT (local input) INTEGER
2616* On entry, ICTXT specifies the BLACS context handle, indica-
2617* ting the global context of the operation. The context itself
2618* is global, but the value of ICTXT is local.
2619*
2620* NOUT (global input) INTEGER
2621* On entry, NOUT specifies the unit number for the output file.
2622* When NOUT is 6, output to screen, when NOUT is 0, output to
2623* stderr. NOUT is only defined for process 0.
2624*
2625* NROUT (global input) INTEGER
2626* On entry, NROUT specifies which routine will be tested as
2627* follows:
2628* If NROUT = 1, PDGEMV will be tested;
2629* else if NROUT = 2, PDSYMV will be tested;
2630* else if NROUT = 3, PDTRMV will be tested;
2631* else if NROUT = 4, PDTRSV will be tested;
2632* else if NROUT = 5, PDGER will be tested;
2633* else if NROUT = 6, PDSYR will be tested;
2634* else if NROUT = 7, PDSYR2 will be tested;
2635*
2636* UPLO (global input) CHARACTER*1
2637* On entry, UPLO specifies if the upper or lower part of the
2638* matrix operand is to be referenced.
2639*
2640* TRANS (global input) CHARACTER*1
2641* On entry, TRANS specifies if the matrix operand A is to be
2642* transposed.
2643*
2644* DIAG (global input) CHARACTER*1
2645* On entry, DIAG specifies if the triangular matrix operand is
2646* unit or non-unit.
2647*
2648* M (global input) INTEGER
2649* On entry, M specifies the number of rows of A.
2650*
2651* N (global input) INTEGER
2652* On entry, N specifies the number of columns of A.
2653*
2654* ALPHA (global input) DOUBLE PRECISION
2655* On entry, ALPHA specifies the scalar alpha.
2656*
2657* A (local input/local output) DOUBLE PRECISION array
2658* On entry, A is an array of dimension (DESCA( M_ ),*). This
2659* array contains a local copy of the initial entire matrix PA.
2660*
2661* PA (local input) DOUBLE PRECISION array
2662* On entry, PA is an array of dimension (DESCA( LLD_ ),*). This
2663* array contains the local entries of the matrix PA.
2664*
2665* IA (global input) INTEGER
2666* On entry, IA specifies A's global row index, which points to
2667* the beginning of the submatrix sub( A ).
2668*
2669* JA (global input) INTEGER
2670* On entry, JA specifies A's global column index, which points
2671* to the beginning of the submatrix sub( A ).
2672*
2673* DESCA (global and local input) INTEGER array
2674* On entry, DESCA is an integer array of dimension DLEN_. This
2675* is the array descriptor for the matrix A.
2676*
2677* X (local input/local output) DOUBLE PRECISION array
2678* On entry, X is an array of dimension (DESCX( M_ ),*). This
2679* array contains a local copy of the initial entire matrix PX.
2680*
2681* PX (local input) DOUBLE PRECISION array
2682* On entry, PX is an array of dimension (DESCX( LLD_ ),*). This
2683* array contains the local entries of the matrix PX.
2684*
2685* IX (global input) INTEGER
2686* On entry, IX specifies X's global row index, which points to
2687* the beginning of the submatrix sub( X ).
2688*
2689* JX (global input) INTEGER
2690* On entry, JX specifies X's global column index, which points
2691* to the beginning of the submatrix sub( X ).
2692*
2693* DESCX (global and local input) INTEGER array
2694* On entry, DESCX is an integer array of dimension DLEN_. This
2695* is the array descriptor for the matrix X.
2696*
2697* INCX (global input) INTEGER
2698* On entry, INCX specifies the global increment for the
2699* elements of X. Only two values of INCX are supported in
2700* this version, namely 1 and M_X. INCX must not be zero.
2701*
2702* BETA (global input) DOUBLE PRECISION
2703* On entry, BETA specifies the scalar beta.
2704*
2705* Y (local input/local output) DOUBLE PRECISION array
2706* On entry, Y is an array of dimension (DESCY( M_ ),*). This
2707* array contains a local copy of the initial entire matrix PY.
2708*
2709* PY (local input) DOUBLE PRECISION array
2710* On entry, PY is an array of dimension (DESCY( LLD_ ),*). This
2711* array contains the local entries of the matrix PY.
2712*
2713* IY (global input) INTEGER
2714* On entry, IY specifies Y's global row index, which points to
2715* the beginning of the submatrix sub( Y ).
2716*
2717* JY (global input) INTEGER
2718* On entry, JY specifies Y's global column index, which points
2719* to the beginning of the submatrix sub( Y ).
2720*
2721* DESCY (global and local input) INTEGER array
2722* On entry, DESCY is an integer array of dimension DLEN_. This
2723* is the array descriptor for the matrix Y.
2724*
2725* INCY (global input) INTEGER
2726* On entry, INCY specifies the global increment for the
2727* elements of Y. Only two values of INCY are supported in
2728* this version, namely 1 and M_Y. INCY must not be zero.
2729*
2730* THRESH (global input) REAL
2731* On entry, THRESH is the threshold value for the test ratio.
2732*
2733* ROGUE (global input) DOUBLE PRECISION
2734* On entry, ROGUE specifies the constant used to pad the
2735* non-referenced part of triangular or symmetric matrices.
2736*
2737* WORK (workspace) DOUBLE PRECISION array
2738* On entry, WORK is an array of dimension LWORK where LWORK is
2739* at least MAX( M, N ). This array is used to store the compu-
2740* ted gauges (see PDMVCH).
2741*
2742* INFO (global output) INTEGER
2743* On exit, if INFO = 0, no error has been found, otherwise
2744* if( MOD( INFO, 2 ) = 1 ) then an error on A has been found,
2745* if( MOD( INFO/2, 2 ) = 1 ) then an error on X has been found,
2746* if( MOD( INFO/4, 2 ) = 1 ) then an error on Y has been found.
2747*
2748* -- Written on April 1, 1998 by
2749* Antoine Petitet, University of Tennessee, Knoxville 37996, USA.
2750*
2751* =====================================================================
2752*
2753* .. Parameters ..
2754 DOUBLE PRECISION ONE, ZERO
2755 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
2756 INTEGER BLOCK_CYCLIC_2D_INB, CSRC_, CTXT_, DLEN_,
2757 $ DTYPE_, IMB_, INB_, LLD_, MB_, M_, NB_, N_,
2758 $ RSRC_
2759 PARAMETER ( BLOCK_CYCLIC_2D_INB = 2, dlen_ = 11,
2760 $ dtype_ = 1, ctxt_ = 2, m_ = 3, n_ = 4,
2761 $ imb_ = 5, inb_ = 6, mb_ = 7, nb_ = 8,
2762 $ rsrc_ = 9, csrc_ = 10, lld_ = 11 )
2763* ..
2764* .. Local Scalars ..
2765 INTEGER I, MYCOL, MYROW, NPCOL, NPROW
2766 DOUBLE PRECISION ERR
2767* ..
2768* .. Local Arrays ..
2769 INTEGER IERR( 3 )
2770* ..
2771* .. External Subroutines ..
2772 EXTERNAL blacs_gridinfo, dtrsv, pb_dlaset, pdchkmin,
2773 $ pdchkvin, pdmvch, pdtrmv, pdvmch, pdvmch2
2774* ..
2775* .. External Functions ..
2776 LOGICAL LSAME
2777 EXTERNAL lsame
2778* ..
2779* .. Intrinsic Functions ..
2780 INTRINSIC dble
2781* ..
2782* .. Executable Statements ..
2783*
2784 info = 0
2785*
2786* Quick return if possible
2787*
2788 IF( ( m.LE.0 ).OR.( n.LE.0 ) )
2789 $ RETURN
2790*
2791* Start the operations
2792*
2793 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
2794*
2795 DO 10 i = 1, 3
2796 ierr( i ) = 0
2797 10 CONTINUE
2798*
2799 IF( nrout.EQ.1 ) THEN
2800*
2801* Test PDGEMV
2802*
2803* Check the resulting vector Y
2804*
2805 CALL pdmvch( ictxt, trans, m, n, alpha, a, ia, ja, desca, x,
2806 $ ix, jx, descx, incx, beta, y, py, iy, jy, descy,
2807 $ incy, work, err, ierr( 3 ) )
2808*
2809 IF( ierr( 3 ).NE.0 ) THEN
2810 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2811 $ WRITE( nout, fmt = 9997 )
2812 ELSE IF( err.GT.dble( thresh ) ) THEN
2813 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2814 $ WRITE( nout, fmt = 9996 ) err
2815 END IF
2816*
2817* Check the input-only arguments
2818*
2819 CALL pdchkmin( err, m, n, a, pa, ia, ja, desca, ierr( 1 ) )
2820 IF( lsame( trans, 'N' ) ) THEN
2821 CALL pdchkvin( err, n, x, px, ix, jx, descx, incx,
2822 $ ierr( 2 ) )
2823 ELSE
2824 CALL pdchkvin( err, m, x, px, ix, jx, descx, incx,
2825 $ ierr( 2 ) )
2826 END IF
2827*
2828 ELSE IF( nrout.EQ.2 ) THEN
2829*
2830* Test PDSYMV
2831*
2832* Check the resulting vector Y
2833*
2834 CALL pdmvch( ictxt, 'No transpose', n, n, alpha, a, ia, ja,
2835 $ desca, x, ix, jx, descx, incx, beta, y, py, iy,
2836 $ jy, descy, incy, work, err, ierr( 3 ) )
2837*
2838 IF( ierr( 3 ).NE.0 ) THEN
2839 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2840 $ WRITE( nout, fmt = 9997 )
2841 ELSE IF( err.GT.dble( thresh ) ) THEN
2842 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2843 $ WRITE( nout, fmt = 9996 ) err
2844 END IF
2845*
2846* Check the input-only arguments
2847*
2848 IF( lsame( uplo, 'L' ) ) THEN
2849 CALL pb_dlaset( 'Upper', n-1, n-1, 0, rogue, rogue,
2850 $ a( ia+ja*desca( m_ ) ), desca( m_ ) )
2851 ELSE
2852 CALL pb_dlaset( 'Lower', n-1, n-1, 0, rogue, rogue,
2853 $ a( ia+1+(ja-1)*desca( m_ ) ), desca( m_ ) )
2854 END IF
2855 CALL pdchkmin( err, n, n, a, pa, ia, ja, desca, ierr( 1 ) )
2856 CALL pdchkvin( err, n, x, px, ix, jx, descx, incx, ierr( 2 ) )
2857*
2858 ELSE IF( nrout.EQ.3 ) THEN
2859*
2860* Test PDTRMV
2861*
2862* Check the resulting vector X
2863*
2864 CALL pdmvch( ictxt, trans, n, n, one, a, ia, ja, desca, y, ix,
2865 $ jx, descx, incx, zero, x, px, ix, jx, descx, incx,
2866 $ work, err, ierr( 2 ) )
2867*
2868 IF( ierr( 2 ).NE.0 ) THEN
2869 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2870 $ WRITE( nout, fmt = 9997 )
2871 ELSE IF( err.GT.dble( thresh ) ) THEN
2872 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2873 $ WRITE( nout, fmt = 9996 ) err
2874 END IF
2875*
2876* Check the input-only arguments
2877*
2878 IF( lsame( uplo, 'L' ) ) THEN
2879 IF( lsame( diag, 'N' ) ) THEN
2880 CALL pb_dlaset( 'Upper', n-1, n-1, 0, rogue, rogue,
2881 $ a( ia+ja*desca( m_ ) ), desca( m_ ) )
2882 ELSE
2883 CALL pb_dlaset( 'Upper', n, n, 0, rogue, one,
2884 $ a( ia+(ja-1)*desca( m_ ) ), desca( m_ ) )
2885 END IF
2886 ELSE
2887 IF( lsame( diag, 'N' ) ) THEN
2888 CALL pb_dlaset( 'Lower', n-1, n-1, 0, rogue, rogue,
2889 $ a( ia+1+(ja-1)*desca( m_ ) ),
2890 $ desca( m_ ) )
2891 ELSE
2892 CALL pb_dlaset( 'Lower', n, n, 0, rogue, one,
2893 $ a( ia+(ja-1)*desca( m_ ) ), desca( m_ ) )
2894 END IF
2895 END IF
2896 CALL pdchkmin( err, n, n, a, pa, ia, ja, desca, ierr( 1 ) )
2897*
2898 ELSE IF( nrout.EQ.4 ) THEN
2899*
2900* Test PDTRSV
2901*
2902* Check the resulting vector X
2903*
2904 CALL dtrsv( uplo, trans, diag, n, a( ia+(ja-1)*desca( m_ ) ),
2905 $ desca( m_ ), x( ix+(jx-1)*descx( m_ ) ), incx )
2906 CALL pdtrmv( uplo, trans, diag, n, pa, ia, ja, desca, px, ix,
2907 $ jx, descx, incx )
2908 CALL pdmvch( ictxt, trans, n, n, one, a, ia, ja, desca, x, ix,
2909 $ jx, descx, incx, zero, y, px, ix, jx, descx, incx,
2910 $ work, err, ierr( 2 ) )
2911*
2912 IF( ierr( 2 ).NE.0 ) THEN
2913 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2914 $ WRITE( nout, fmt = 9997 )
2915 ELSE IF( err.GT.dble( thresh ) ) THEN
2916 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2917 $ WRITE( nout, fmt = 9996 ) err
2918 END IF
2919*
2920* Check the input-only arguments
2921*
2922 IF( lsame( uplo, 'L' ) ) THEN
2923 IF( lsame( diag, 'N' ) ) THEN
2924 CALL pb_dlaset( 'Upper', n-1, n-1, 0, rogue, rogue,
2925 $ a( ia+ja*desca( m_ ) ), desca( m_ ) )
2926 ELSE
2927 CALL pb_dlaset( 'Upper', n, n, 0, rogue, one,
2928 $ a( ia+(ja-1)*desca( m_ ) ), desca( m_ ) )
2929 END IF
2930 ELSE
2931 IF( lsame( diag, 'N' ) ) THEN
2932 CALL pb_dlaset( 'Lower', n-1, n-1, 0, rogue, rogue,
2933 $ a( ia+1+(ja-1)*desca( m_ ) ),
2934 $ desca( m_ ) )
2935 ELSE
2936 CALL pb_dlaset( 'Lower', n, n, 0, rogue, one,
2937 $ a( ia+(ja-1)*desca( m_ ) ), desca( m_ ) )
2938 END IF
2939 END IF
2940 CALL pdchkmin( err, n, n, a, pa, ia, ja, desca, ierr( 1 ) )
2941*
2942 ELSE IF( nrout.EQ.5 ) THEN
2943*
2944* Test PDGER
2945*
2946* Check the resulting matrix A
2947*
2948 CALL pdvmch( ictxt, 'Ge', m, n, alpha, x, ix, jx, descx,
2949 $ incx, y, iy, jy, descy, incy, a, pa, ia, ja,
2950 $ desca, work, err, ierr( 1 ) )
2951 IF( ierr( 1 ).NE.0 ) THEN
2952 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2953 $ WRITE( nout, fmt = 9997 )
2954 ELSE IF( err.GT.dble( thresh ) ) THEN
2955 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2956 $ WRITE( nout, fmt = 9996 ) err
2957 END IF
2958*
2959* Check the input-only arguments
2960*
2961 CALL pdchkvin( err, m, x, px, ix, jx, descx, incx, ierr( 2 ) )
2962 CALL pdchkvin( err, n, y, py, iy, jy, descy, incy, ierr( 3 ) )
2963*
2964 ELSE IF( nrout.EQ.6 ) THEN
2965*
2966* Test PDSYR
2967*
2968* Check the resulting matrix A
2969*
2970 CALL pdvmch( ictxt, uplo, n, n, alpha, x, ix, jx, descx,
2971 $ incx, x, ix, jx, descx, incx, a, pa, ia, ja,
2972 $ desca, work, err, ierr( 1 ) )
2973 IF( ierr( 1 ).NE.0 ) THEN
2974 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2975 $ WRITE( nout, fmt = 9997 )
2976 ELSE IF( err.GT.dble( thresh ) ) THEN
2977 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2978 $ WRITE( nout, fmt = 9996 ) err
2979 END IF
2980*
2981* Check the input-only arguments
2982*
2983 CALL pdchkvin( err, n, x, px, ix, jx, descx, incx, ierr( 2 ) )
2984*
2985 ELSE IF( nrout.EQ.7 ) THEN
2986*
2987* Test PDSYR2
2988*
2989* Check the resulting matrix A
2990*
2991 CALL pdvmch2( ictxt, uplo, n, n, alpha, x, ix, jx, descx, incx,
2992 $ y, iy, jy, descy, incy, a, pa, ia, ja, desca,
2993 $ work, err, ierr( 1 ) )
2994 IF( ierr( 1 ).NE.0 ) THEN
2995 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2996 $ WRITE( nout, fmt = 9997 )
2997 ELSE IF( err.GT.dble( thresh ) ) THEN
2998 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
2999 $ WRITE( nout, fmt = 9996 ) err
3000 END IF
3001*
3002* Check the input-only arguments
3003*
3004 CALL pdchkvin( err, n, x, px, ix, jx, descx, incx, ierr( 2 ) )
3005 CALL pdchkvin( err, n, y, py, iy, jy, descy, incy, ierr( 3 ) )
3006*
3007 END IF
3008*
3009 IF( ierr( 1 ).NE.0 ) THEN
3010 info = info + 1
3011 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
3012 $ WRITE( nout, fmt = 9999 ) 'A'
3013 END IF
3014*
3015 IF( ierr( 2 ).NE.0 ) THEN
3016 info = info + 2
3017 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
3018 $ WRITE( nout, fmt = 9998 ) 'X'
3019 END IF
3020*
3021 IF( ierr( 3 ).NE.0 ) THEN
3022 info = info + 4
3023 IF( myrow.EQ.0 .AND. mycol.EQ.0 )
3024 $ WRITE( nout, fmt = 9998 ) 'Y'
3025 END IF
3026*
3027 9999 FORMAT( 2x, ' ***** ERROR: Matrix operand ', a,
3028 $ ' is incorrect.' )
3029 9998 FORMAT( 2x, ' ***** ERROR: Vector operand ', a,
3030 $ ' is incorrect.' )
3031 9997 FORMAT( 2x, ' ***** FATAL ERROR - Computed result is less ',
3032 $ 'than half accurate *****' )
3033 9996 FORMAT( 2x, ' ***** Test completed with maximum test ratio: ',
3034 $ f11.5, ' SUSPECT *****' )
3035*
3036 RETURN
3037*
3038* End of PDBLAS2TSTCHK
3039*
3040 END
subroutine pmdimchk(ictxt, nout, m, n, matrix, ia, ja, desca, info)
Definition pblastst.f:202
subroutine pvdimchk(ictxt, nout, n, matrix, ix, jx, descx, incx, info)
Definition pblastst.f:3
subroutine icopy(n, sx, incx, sy, incy)
Definition pblastst.f:1525
subroutine pb_descset2(desc, m, n, imb, inb, mb, nb, rsrc, csrc, ctxt, lld)
Definition pblastst.f:3172
subroutine pmdescchk(ictxt, nout, matrix, desca, dta, ma, na, imba, inba, mba, nba, rsrca, csrca, mpa, nqa, iprea, imida, iposta, igap, gapmul, info)
Definition pblastst.f:746
subroutine pvdescchk(ictxt, nout, matrix, descx, dtx, mx, nx, imbx, inbx, mbx, nbx, rsrcx, csrcx, incx, mpx, nqx, iprex, imidx, ipostx, igap, gapmul, info)
Definition pblastst.f:388
#define max(A, B)
Definition pcgemr.c:180
#define min(A, B)
Definition pcgemr.c:181
subroutine pdchkarg2(ictxt, nout, sname, uplo, trans, diag, m, n, alpha, ia, ja, desca, ix, jx, descx, incx, beta, iy, jy, descy, incy, info)
subroutine pdblas2tstchke(ltest, inout, nprocs)
program pdbla2tst
Definition pdblas2tst.f:11
subroutine pdbla2tstinfo(summry, nout, nmat, diagval, tranval, uploval, mval, nval, maval, naval, imbaval, mbaval, inbaval, nbaval, rscaval, cscaval, iaval, javal, mxval, nxval, imbxval, mbxval, inbxval, nbxval, rscxval, cscxval, ixval, jxval, incxval, myval, nyval, imbyval, mbyval, inbyval, nbyval, rscyval, cscyval, iyval, jyval, incyval, ldval, ngrids, pval, ldpval, qval, ldqval, nblog, ltest, sof, tee, iam, igap, iverb, nprocs, thresh, alpha, beta, work)
subroutine pdblas2tstchk(ictxt, nout, nrout, uplo, trans, diag, m, n, alpha, a, pa, ia, ja, desca, x, px, ix, jx, descx, incx, beta, y, py, iy, jy, descy, incy, thresh, rogue, work, info)
subroutine pdchkmout(m, n, a, pa, ia, ja, desca, info)
Definition pdblastst.f:3627
subroutine pb_pdlaprnt(m, n, a, ia, ja, desca, irprnt, icprnt, cmatnm, nout, work)
Definition pdblastst.f:8634
subroutine pb_dlaset(uplo, m, n, ioffd, alpha, beta, a, lda)
Definition pdblastst.f:9359
subroutine pdchkvout(n, x, px, ix, jx, descx, incx, info)
Definition pdblastst.f:2870
subroutine pdvmch2(ictxt, uplo, m, n, alpha, x, ix, jx, descx, incx, y, iy, jy, descy, incy, a, pa, ia, ja, desca, g, err, info)
Definition pdblastst.f:4919
subroutine pdlagen(inplace, aform, diag, offa, m, n, ia, ja, desca, iaseed, a, lda)
Definition pdblastst.f:7845
subroutine pdoptee(ictxt, nout, subptr, scode, sname)
Definition pdblastst.f:2
subroutine pdmprnt(ictxt, nout, m, n, a, lda, irprnt, icprnt, cmatnm)
Definition pdblastst.f:3949
subroutine pdlascal(type, m, n, alpha, a, ia, ja, desca)
Definition pdblastst.f:7337
subroutine pb_dfillpad(ictxt, m, n, a, lda, ipre, ipost, chkval)
Definition pdblastst.f:9079
subroutine pdlaset(uplo, m, n, alpha, beta, a, ia, ja, desca)
Definition pdblastst.f:6862
subroutine pdvmch(ictxt, uplo, m, n, alpha, x, ix, jx, descx, incx, y, iy, jy, descy, incy, a, pa, ia, ja, desca, g, err, info)
Definition pdblastst.f:4570
subroutine pb_dchekpad(ictxt, mess, m, n, a, lda, ipre, ipost, chkval)
Definition pdblastst.f:9192
subroutine pdvprnt(ictxt, nout, n, x, incx, irprnt, icprnt, cvecnm)
Definition pdblastst.f:4056
subroutine pdmvch(ictxt, trans, m, n, alpha, a, ia, ja, desca, x, ix, jx, descx, incx, beta, y, py, iy, jy, descy, incy, g, err, info)
Definition pdblastst.f:4157
subroutine pdchkvin(errmax, n, x, px, ix, jx, descx, incx, info)
Definition pdblastst.f:2576
subroutine pb_dlascal(uplo, m, n, ioffd, alpha, a, lda)
Definition pdblastst.f:9556
subroutine pdmatee(ictxt, nout, subptr, scode, sname)
Definition pdblastst.f:1190
subroutine pddimee(ictxt, nout, subptr, scode, sname)
Definition pdblastst.f:455
subroutine pdchkmin(errmax, m, n, a, pa, ia, ja, desca, info)
Definition pdblastst.f:3326
subroutine pdvecee(ictxt, nout, subptr, scode, sname)
Definition pdblastst.f:936
double precision function pdlamch(ictxt, cmach)
Definition pdblastst.f:6769
logical function lsame(ca, cb)
Definition tools.f:1724