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