ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
pdnepdriver.f
Go to the documentation of this file.
1  PROGRAM pdnepdriver
2 *
3 * -- ScaLAPACK testing driver (version 1.7) --
4 * University of Tennessee, Knoxville, Oak Ridge National Laboratory,
5 * and University of California, Berkeley.
6 * May 1, 1997
7 *
8 * Purpose
9 * =======
10 *
11 * PDNEPDRIVER is the main test program for the DOUBLE PRECISION
12 * SCALAPACK NEP routines. This test driver performs a Schur
13 * decomposition followed by residual check of a Hessenberg matrix.
14 *
15 * The program must be driven by a short data file. An annotated
16 * example of a data file can be obtained by deleting the first 3
17 * characters from the following 18 lines:
18 * 'SCALAPACK, Version 1.4, NEP (Nonsymmetric EigenProblem) input file'
19 * 'Intel iPSC/860 hypercube, gamma model.'
20 * 'NEP.out' output file name (if any)
21 * 6 device out
22 * 8 number of problems sizes
23 * 1 2 3 4 6 10 100 200 vales of N
24 * 3 number of NB's
25 * 6 20 40 values of NB
26 * 4 number of process grids (ordered pairs of P & Q)
27 * 1 2 1 4 values of P
28 * 1 2 4 1 values of Q
29 * 20.0 threshold
30 *
31 * Internal Parameters
32 * ===================
33 *
34 * TOTMEM INTEGER, default = 2000000
35 * TOTMEM is a machine-specific parameter indicating the
36 * maximum amount of available memory in bytes.
37 * The user should customize TOTMEM to his platform. Remember
38 * to leave room in memory for the operating system, the BLACS
39 * buffer, etc. For example, on a system with 8 MB of memory
40 * per process (e.g., one processor on an Intel iPSC/860), the
41 * parameters we use are TOTMEM=6200000 (leaving 1.8 MB for OS,
42 * code, BLACS buffer, etc). However, for PVM, we usually set
43 * TOTMEM = 2000000. Some experimenting with the maximum value
44 * of TOTMEM may be required.
45 *
46 * DBLESZ INTEGER, default = 8 bytes.
47 * DBLESZ indicate the length in bytes on the given platform
48 * for a double precision real.
49 * MEM DOUBLE PRECISION array, dimension ( TOTMEM / DBLESZ )
50 *
51 * All arrays used by SCALAPACK routines are allocated from
52 * this array and referenced by pointers. The integer IPA,
53 * for example, is a pointer to the starting element of MEM for
54 * the matrix A.
55 *
56 * =====================================================================
57 *
58 * .. Parameters ..
59  INTEGER block_cyclic_2d, csrc_, ctxt_, dlen_, dt_,
60  $ lld_, mb_, m_, nb_, n_, rsrc_
61  parameter( block_cyclic_2d = 1, dlen_ = 9, dt_ = 1,
62  $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
63  $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
64  INTEGER dblesz, totmem, memsiz, ntests
65  DOUBLE PRECISION padval, zero, one
66  parameter( dblesz = 8, totmem = 2000000,
67  $ memsiz = totmem / dblesz, ntests = 20,
68  $ padval = -9923.0d+0, zero = 0.0d+0,
69  $ one = 1.0d+0 )
70 * ..
71 * .. Local Scalars ..
72  LOGICAL check
73  CHARACTER*6 passed
74  CHARACTER*80 outfile
75  INTEGER i, iam, iaseed, ictxt, iii, imidpad, info, ipa,
76  $ ipostpad, iprepad, ipw, ipwi, ipwr, ipz, j, k,
77  $ kfail, kpass, kskip, ktests, lda, ldz, lwork,
78  $ mycol, myrow, n, nb, ngrids, nmat, nnb, nout,
79  $ np, npcol, nprocs, nprow, nq, worksiz
80  REAL thresh
81  DOUBLE PRECISION anorm, fresid, nops, qresid, tmflops, znorm
82 * ..
83 * .. Local Arrays ..
84  INTEGER desca( dlen_ ), descz( dlen_ ), ierr( 2 ),
85  $ idum( 1 ), nbval( ntests ), nval( ntests ),
86  $ pval( ntests ), qval( ntests )
87  DOUBLE PRECISION ctime( 1 ), mem( memsiz ), wtime( 1 )
88 * ..
89 * .. External Subroutines ..
90  EXTERNAL blacs_barrier, blacs_exit, blacs_get,
91  $ blacs_gridexit, blacs_gridinfo, blacs_gridinit,
92  $ blacs_pinfo, descinit, igsum2d, pdchekpad,
93  $ pdfillpad, pdgemm, pdlahqr, pdlaset, pdmatgen,
95  $ sltimer
96 * ..
97 * .. External Functions ..
98  INTEGER ilcm, numroc
99  DOUBLE PRECISION pdlamch, pdlange, pdlanhs
100  EXTERNAL ilcm, numroc, pdlamch, pdlange, pdlanhs
101 * ..
102 * .. Intrinsic Functions ..
103  INTRINSIC dble, max, min
104 * ..
105 * .. Data statements ..
106  DATA kfail, kpass, kskip, ktests / 4*0 /
107 * ..
108 * .. Executable Statements ..
109 *
110 * Get starting information
111 *
112  CALL blacs_pinfo( iam, nprocs )
113  iaseed = 100
114  CALL pdnepinfo( outfile, nout, nmat, nval, ntests, nnb, nbval,
115  $ ntests, ngrids, pval, ntests, qval, ntests,
116  $ thresh, mem, iam, nprocs )
117  check = ( thresh.GE.0.0e+0 )
118 *
119 * Print headings
120 *
121  IF( iam.EQ.0 ) THEN
122  WRITE( nout, fmt = * )
123  WRITE( nout, fmt = 9995 )
124  WRITE( nout, fmt = 9994 )
125  WRITE( nout, fmt = * )
126  END IF
127 *
128 * Loop over different process grids
129 *
130  DO 30 i = 1, ngrids
131 *
132  nprow = pval( i )
133  npcol = qval( i )
134 *
135 * Make sure grid information is correct
136 *
137  ierr( 1 ) = 0
138  IF( nprow.LT.1 ) THEN
139  IF( iam.EQ.0 )
140  $ WRITE( nout, fmt = 9999 )'GRID', 'nprow', nprow
141  ierr( 1 ) = 1
142  ELSE IF( npcol.LT.1 ) THEN
143  IF( iam.EQ.0 )
144  $ WRITE( nout, fmt = 9999 )'GRID', 'npcol', npcol
145  ierr( 1 ) = 1
146  ELSE IF( nprow*npcol.GT.nprocs ) THEN
147  IF( iam.EQ.0 )
148  $ WRITE( nout, fmt = 9998 )nprow*npcol, nprocs
149  ierr( 1 ) = 1
150  END IF
151 *
152  IF( ierr( 1 ).GT.0 ) THEN
153  IF( iam.EQ.0 )
154  $ WRITE( nout, fmt = 9997 )'grid'
155  kskip = kskip + 1
156  GO TO 30
157  END IF
158 *
159 * Define process grid
160 *
161  CALL blacs_get( -1, 0, ictxt )
162  CALL blacs_gridinit( ictxt, 'Row-major', nprow, npcol )
163  CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
164 *
165 * Go to bottom of process grid loop if this case doesn't use my
166 * process
167 *
168  IF( myrow.GE.nprow .OR. mycol.GE.npcol )
169  $ GO TO 30
170 *
171  DO 20 j = 1, nmat
172 *
173  n = nval( j )
174 *
175 * Make sure matrix information is correct
176 *
177  ierr( 1 ) = 0
178  IF( n.LT.1 ) THEN
179  IF( iam.EQ.0 )
180  $ WRITE( nout, fmt = 9999 )'MATRIX', 'N', n
181  ierr( 1 ) = 1
182  END IF
183 *
184 * Check all processes for an error
185 *
186  CALL igsum2d( ictxt, 'All', ' ', 1, 1, ierr, 1, -1, 0 )
187 *
188  IF( ierr( 1 ).GT.0 ) THEN
189  IF( iam.EQ.0 )
190  $ WRITE( nout, fmt = 9997 )'matrix'
191  kskip = kskip + 1
192  GO TO 20
193  END IF
194 *
195  DO 10 k = 1, nnb
196 *
197  nb = nbval( k )
198 *
199 * Make sure nb is legal
200 *
201  ierr( 1 ) = 0
202  IF( nb.LT.6 ) THEN
203  ierr( 1 ) = 1
204  IF( iam.EQ.0 )
205  $ WRITE( nout, fmt = 9999 )'NB', 'NB', nb
206  END IF
207 *
208 * Check all processes for an error
209 *
210  CALL igsum2d( ictxt, 'All', ' ', 1, 1, ierr, 1, -1, 0 )
211 *
212  IF( ierr( 1 ).GT.0 ) THEN
213  IF( iam.EQ.0 )
214  $ WRITE( nout, fmt = 9997 )'NB'
215  kskip = kskip + 1
216  GO TO 10
217  END IF
218 *
219 * Padding constants
220 *
221  np = numroc( n, nb, myrow, 0, nprow )
222  nq = numroc( n, nb, mycol, 0, npcol )
223  IF( check ) THEN
224  iprepad = max( nb, np )
225  imidpad = nb
226  ipostpad = max( nb, nq )
227  iprepad = iprepad + 1000
228  imidpad = imidpad + 1000
229  ipostpad = ipostpad + 1000
230  ELSE
231  iprepad = 0
232  imidpad = 0
233  ipostpad = 0
234  END IF
235 *
236 * Initialize the array descriptor for the matrix A
237 *
238  CALL descinit( desca, n, n, nb, nb, 0, 0, ictxt,
239  $ max( 1, np )+imidpad, ierr( 1 ) )
240 *
241 * Initialize the array descriptor for the matrix Z
242 *
243  CALL descinit( descz, n, n, nb, nb, 0, 0, ictxt,
244  $ max( 1, np )+imidpad, ierr( 2 ) )
245 *
246  lda = desca( lld_ )
247  ldz = descz( lld_ )
248 *
249 * Check all processes for an error
250 *
251  CALL igsum2d( ictxt, 'All', ' ', 2, 1, ierr, 2, -1, 0 )
252 *
253  IF( ierr( 1 ).LT.0 .OR. ierr( 2 ).LT.0 ) THEN
254  IF( iam.EQ.0 )
255  $ WRITE( nout, fmt = 9997 )'descriptor'
256  kskip = kskip + 1
257  GO TO 10
258  END IF
259 *
260 * Assign pointers into MEM for SCALAPACK arrays, A is
261 * allocated starting at position MEM( IPREPAD+1 )
262 *
263  ipa = iprepad + 1
264  ipz = ipa + desca( lld_ )*nq + ipostpad + iprepad
265  ipwr = ipz + descz( lld_ )*nq + ipostpad + iprepad
266  ipwi = ipwr + n + ipostpad + iprepad
267  ipw = ipwi + n + ipostpad + iprepad
268  iii = n / nb
269  IF( iii*nb.LT.n )
270  $ iii = iii + 1
271  iii = 7*iii / ilcm( nprow, npcol )
272 *
273 *
274  lwork = 3*n + max( 2*max( lda, ldz )+2*nq, iii )
275  lwork = lwork + max(2*n, (8*ilcm(nprow,npcol)+2)**2 )
276 *
277  IF( check ) THEN
278 *
279 * Figure the amount of workspace required by the
280 * checking routines PDNEPFCHK and PDLANHS
281 *
282  worksiz = lwork + max( np*desca( nb_ ),
283  $ desca( mb_ )*nq ) + ipostpad
284 *
285  ELSE
286 *
287  worksiz = lwork + ipostpad
288 *
289  END IF
290 *
291 * Check for adequate memory for problem size
292 *
293  ierr( 1 ) = 0
294  IF( ipw+worksiz.GT.memsiz ) THEN
295  IF( iam.EQ.0 )
296  $ WRITE( nout, fmt = 9996 )'Schur reduction',
297  $ ( ipw+worksiz )*dblesz
298  ierr( 1 ) = 1
299  END IF
300 *
301 * Check all processes for an error
302 *
303  CALL igsum2d( ictxt, 'All', ' ', 1, 1, ierr, 1, -1, 0 )
304 *
305  IF( ierr( 1 ).GT.0 ) THEN
306  IF( iam.EQ.0 )
307  $ WRITE( nout, fmt = 9997 )'MEMORY'
308  kskip = kskip + 1
309  GO TO 10
310  END IF
311 *
312 * Generate matrix Z = In
313 *
314  CALL pdlaset( 'All', n, n, zero, one, mem( ipz ), 1, 1,
315  $ descz )
316 *
317 * Generate matrix A upper Hessenberg
318 *
319  CALL pdmatgen( ictxt, 'No transpose', 'No transpose',
320  $ desca( m_ ), desca( n_ ), desca( mb_ ),
321  $ desca( nb_ ), mem( ipa ), desca( lld_ ),
322  $ desca( rsrc_ ), desca( csrc_ ), iaseed, 0,
323  $ np, 0, nq, myrow, mycol, nprow, npcol )
324  CALL pdlaset( 'Lower', max( 0, n-2 ), max( 0, n-2 ),
325  $ zero, zero, mem( ipa ), min( n, 3 ), 1,
326  $ desca )
327 *
328 * Calculate inf-norm of A for residual error-checking
329 *
330  IF( check ) THEN
331  CALL pdfillpad( ictxt, np, nq, mem( ipa-iprepad ),
332  $ desca( lld_ ), iprepad, ipostpad,
333  $ padval )
334  CALL pdfillpad( ictxt, np, nq, mem( ipz-iprepad ),
335  $ descz( lld_ ), iprepad, ipostpad,
336  $ padval )
337  CALL pdfillpad( ictxt, worksiz-ipostpad, 1,
338  $ mem( ipw-iprepad ), worksiz-ipostpad,
339  $ iprepad, ipostpad, padval )
340  anorm = pdlanhs( 'I', n, mem( ipa ), 1, 1, desca,
341  $ mem( ipw ) )
342  CALL pdchekpad( ictxt, 'PDLANHS', np, nq,
343  $ mem( ipa-iprepad ), desca( lld_ ),
344  $ iprepad, ipostpad, padval )
345  CALL pdchekpad( ictxt, 'PDLANHS', worksiz-ipostpad, 1,
346  $ mem( ipw-iprepad ), worksiz-ipostpad,
347  $ iprepad, ipostpad, padval )
348 *
349  CALL pdfillpad( ictxt, n, 1, mem( ipwr-iprepad ), n,
350  $ iprepad, ipostpad, padval )
351  CALL pdfillpad( ictxt, n, 1, mem( ipwi-iprepad ), n,
352  $ iprepad, ipostpad, padval )
353  CALL pdfillpad( ictxt, lwork, 1, mem( ipw-iprepad ),
354  $ lwork, iprepad, ipostpad, padval )
355 *
356  END IF
357 *
358  CALL slboot( )
359  CALL blacs_barrier( ictxt, 'All' )
360  CALL sltimer( 1 )
361 *
362 * Perform NEP factorization
363 *
364  CALL pdlahqr( .true., .true., n, 1, n, mem( ipa ), desca,
365  $ mem( ipwr ), mem( ipwi ), 1, n, mem( ipz ),
366  $ descz, mem( ipw ), lwork, idum, 0, info )
367 *
368  CALL sltimer( 1 )
369 *
370  IF( info.NE.0 ) THEN
371  IF( iam.EQ.0 )
372  $ WRITE( nout, fmt = * )'PDLAHQR INFO=', info
373  kfail = kfail + 1
374  GO TO 10
375  END IF
376 *
377  IF( check ) THEN
378 *
379 * Check for memory overwrite in NEP factorization
380 *
381  CALL pdchekpad( ictxt, 'PDLAHQR (A)', np, nq,
382  $ mem( ipa-iprepad ), desca( lld_ ),
383  $ iprepad, ipostpad, padval )
384  CALL pdchekpad( ictxt, 'PDLAHQR (Z)', np, nq,
385  $ mem( ipz-iprepad ), descz( lld_ ),
386  $ iprepad, ipostpad, padval )
387  CALL pdchekpad( ictxt, 'PDLAHQR (WR)', n, 1,
388  $ mem( ipwr-iprepad ), n, iprepad,
389  $ ipostpad, padval )
390  CALL pdchekpad( ictxt, 'PDLAHQR (WI)', n, 1,
391  $ mem( ipwi-iprepad ), n, iprepad,
392  $ ipostpad, padval )
393  CALL pdchekpad( ictxt, 'PDLAHQR (WORK)', lwork, 1,
394  $ mem( ipw-iprepad ), lwork, iprepad,
395  $ ipostpad, padval )
396 *
397  CALL pdfillpad( ictxt, worksiz-ipostpad, 1,
398  $ mem( ipw-iprepad ), worksiz-ipostpad,
399  $ iprepad, ipostpad, padval )
400 *
401 * Compute || Z * H * Z**T - H0 || / ( N*|| H0 ||*EPS )
402 *
403  CALL pdnepfchk( n, mem( ipa ), 1, 1, desca, iaseed,
404  $ mem( ipz ), 1, 1, descz, anorm,
405  $ fresid, mem( ipw ) )
406 *
407  CALL pdchekpad( ictxt, 'PDNEPFCHK (A)', np, nq,
408  $ mem( ipa-iprepad ), desca( lld_ ),
409  $ iprepad, ipostpad, padval )
410  CALL pdchekpad( ictxt, 'PDNEPFCHK (Z)', np, nq,
411  $ mem( ipz-iprepad ), descz( lld_ ),
412  $ iprepad, ipostpad, padval )
413  CALL pdchekpad( ictxt, 'PDNEPFCHK (WORK)',
414  $ worksiz-ipostpad, 1,
415  $ mem( ipw-iprepad ), worksiz-ipostpad,
416  $ iprepad, ipostpad, padval )
417 *
418 * Compute || (Z**T)*Z - In ||_1
419 *
420  CALL pdlaset( 'All', n, n, zero, one, mem( ipa ), 1,
421  $ 1, desca )
422  CALL pdgemm( 'Transpose', 'No transpose', n, n, n,
423  $ -one, mem( ipz ), 1, 1, descz,
424  $ mem( ipz ), 1, 1, descz, one, mem( ipa ),
425  $ 1, 1, desca )
426  znorm = pdlange( '1', n, n, mem( ipa ), 1, 1, desca,
427  $ mem( ipw ) )
428  qresid = znorm / ( dble( n )*pdlamch( ictxt, 'P' ) )
429 *
430 * Test residual and detect NaN result
431 *
432  IF( ( fresid.LE.thresh ) .AND.
433  $ ( ( fresid-fresid ).EQ.0.0d+0 ) .AND.
434  $ ( qresid.LE.thresh ) .AND.
435  $ ( ( qresid-qresid ).EQ.0.0d+0 ) ) THEN
436  kpass = kpass + 1
437  passed = 'PASSED'
438  ELSE
439  kfail = kfail + 1
440  passed = 'FAILED'
441  IF( iam.EQ.0 ) THEN
442  WRITE( nout, fmt = 9986 )fresid
443  WRITE( nout, fmt = 9985 )qresid
444  END IF
445  END IF
446 *
447  ELSE
448 *
449 * Don't perform the checking, only timing
450 *
451  kpass = kpass + 1
452  fresid = fresid - fresid
453  qresid = qresid - qresid
454  passed = 'BYPASS'
455 *
456  END IF
457 *
458 * Gather maximum of all CPU and WALL clock timings
459 *
460  CALL slcombine( ictxt, 'All', '>', 'W', 1, 1, wtime )
461  CALL slcombine( ictxt, 'All', '>', 'C', 1, 1, ctime )
462 *
463 * Print results
464 *
465  IF( myrow.EQ.0 .AND. mycol.EQ.0 ) THEN
466 *
467 * 18 N^3 flops for PxLAHQR
468 *
469  nops = 18.0d+0*dble( n )**3
470 *
471 * Calculate total megaflops -- factorization only,
472 * -- for WALL and CPU time, and print output
473 *
474 * Print WALL time if machine supports it
475 *
476  IF( wtime( 1 ).GT.0.0d+0 ) THEN
477  tmflops = nops / ( wtime( 1 )*1.0d+6 )
478  ELSE
479  tmflops = 0.0d+0
480  END IF
481  IF( wtime( 1 ).GE.0.0d+0 )
482  $ WRITE( nout, fmt = 9993 )'WALL', n, nb, nprow,
483  $ npcol, wtime( 1 ), tmflops, passed
484 *
485 * Print CPU time if machine supports it
486 *
487  IF( ctime( 1 ).GT.0.0d+0 ) THEN
488  tmflops = nops / ( ctime( 1 )*1.0d+6 )
489  ELSE
490  tmflops = 0.0d+0
491  END IF
492 *
493  IF( ctime( 1 ).GE.0.0d+0 )
494  $ WRITE( nout, fmt = 9993 )'CPU ', n, nb, nprow,
495  $ npcol, ctime( 1 ), tmflops, passed
496  END IF
497 *
498  10 CONTINUE
499 *
500  20 CONTINUE
501 *
502  CALL blacs_gridexit( ictxt )
503 *
504  30 CONTINUE
505 *
506 * Print ending messages and close output file
507 *
508  IF( iam.EQ.0 ) THEN
509  ktests = kpass + kfail + kskip
510  WRITE( nout, fmt = * )
511  WRITE( nout, fmt = 9992 )ktests
512  IF( check ) THEN
513  WRITE( nout, fmt = 9991 )kpass
514  WRITE( nout, fmt = 9989 )kfail
515  ELSE
516  WRITE( nout, fmt = 9990 )kpass
517  END IF
518  WRITE( nout, fmt = 9988 )kskip
519  WRITE( nout, fmt = * )
520  WRITE( nout, fmt = * )
521  WRITE( nout, fmt = 9987 )
522  IF( nout.NE.6 .AND. nout.NE.0 )
523  $ CLOSE ( nout )
524  END IF
525 *
526  CALL blacs_exit( 0 )
527 *
528  9999 FORMAT( 'ILLEGAL ', a6, ': ', a5, ' = ', i3,
529  $ '; It should be at least 1' )
530  9998 FORMAT( 'ILLEGAL GRID: nprow*npcol = ', i4, '. It can be at most',
531  $ i4 )
532  9997 FORMAT( 'Bad ', a6, ' parameters: going on to next test case.' )
533  9996 FORMAT( 'Unable to perform ', a, ': need TOTMEM of at least',
534  $ i11 )
535  9995 FORMAT( 'TIME N NB P Q NEP Time MFLOPS CHECK' )
536  9994 FORMAT( '---- ----- --- ---- ---- -------- -------- ------' )
537  9993 FORMAT( a4, 1x, i5, 1x, i3, 1x, i4, 1x, i4, 1x, f8.2, 1x, f8.2,
538  $ 1x, a6 )
539  9992 FORMAT( 'Finished ', i6, ' tests, with the following results:' )
540  9991 FORMAT( i5, ' tests completed and passed residual checks.' )
541  9990 FORMAT( i5, ' tests completed without checking.' )
542  9989 FORMAT( i5, ' tests completed and failed residual checks.' )
543  9988 FORMAT( i5, ' tests skipped because of illegal input values.' )
544  9987 FORMAT( 'END OF TESTS.' )
545  9986 FORMAT( '||H - Q*S*Q^T|| / (||H|| * N * eps) = ', g25.7 )
546  9985 FORMAT( '||Q^T*Q - I|| / ( N * eps ) ', g25.7 )
547 *
548  stop
549 *
550 * End of PDNEPDRIVER
551 *
552  END
max
#define max(A, B)
Definition: pcgemr.c:180
pdlanhs
double precision function pdlanhs(NORM, N, A, IA, JA, DESCA, WORK)
Definition: pdlanhs.f:3
pdnepinfo
subroutine pdnepinfo(SUMMRY, NOUT, NMAT, NVAL, LDNVAL, NNB, NBVAL, LDNBVAL, NGRIDS, PVAL, LDPVAL, QVAL, LDQVAL, THRESH, WORK, IAM, NPROCS)
Definition: pdnepinfo.f:4
pdchekpad
subroutine pdchekpad(ICTXT, MESS, M, N, A, LDA, IPRE, IPOST, CHKVAL)
Definition: pdchekpad.f:3
ilcm
integer function ilcm(M, N)
Definition: ilcm.f:2
sltimer
subroutine sltimer(I)
Definition: sltimer.f:47
pdnepfchk
subroutine pdnepfchk(N, A, IA, JA, DESCA, IASEED, Z, IZ, JZ, DESCZ, ANORM, FRESID, WORK)
Definition: pdnepfchk.f:3
descinit
subroutine descinit(DESC, M, N, MB, NB, IRSRC, ICSRC, ICTXT, LLD, INFO)
Definition: descinit.f:3
slboot
subroutine slboot()
Definition: sltimer.f:2
pdmatgen
subroutine pdmatgen(ICTXT, AFORM, DIAG, M, N, MB, NB, A, LDA, IAROW, IACOL, ISEED, IROFF, IRNUM, ICOFF, ICNUM, MYROW, MYCOL, NPROW, NPCOL)
Definition: pdmatgen.f:4
pdlange
double precision function pdlange(NORM, M, N, A, IA, JA, DESCA, WORK)
Definition: pdlange.f:3
pdlaset
subroutine pdlaset(UPLO, M, N, ALPHA, BETA, A, IA, JA, DESCA)
Definition: pdblastst.f:6862
numroc
integer function numroc(N, NB, IPROC, ISRCPROC, NPROCS)
Definition: numroc.f:2
pdfillpad
subroutine pdfillpad(ICTXT, M, N, A, LDA, IPRE, IPOST, CHKVAL)
Definition: pdfillpad.f:2
pdnepdriver
program pdnepdriver
Definition: pdnepdriver.f:1
pdlamch
double precision function pdlamch(ICTXT, CMACH)
Definition: pdblastst.f:6769
slcombine
subroutine slcombine(ICTXT, SCOPE, OP, TIMETYPE, N, IBEG, TIMES)
Definition: sltimer.f:267
pdlahqr
subroutine pdlahqr(WANTT, WANTZ, N, ILO, IHI, A, DESCA, WR, WI, ILOZ, IHIZ, Z, DESCZ, WORK, LWORK, IWORK, ILWORK, INFO)
Definition: pdlahqr.f:4
min
#define min(A, B)
Definition: pcgemr.c:181