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