OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
pzbrddriver.f
Go to the documentation of this file.
1 PROGRAM pzbrddriver
2*
3* -- ScaLAPACK testing driver (version 1.7) --
4* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
5* and University of California, Berkeley.
6* March 13, 2000
7*
8* Purpose
9* =======
10*
11* PZBRDDRIVER is the main test program for the COMPLEX*16
12* ScaLAPACK BRD (bidiagonal reduction) routines.
13*
14* The program must be driven by a short data file. An annotated
15* example of a data file can be obtained by deleting the first 3
16* characters from the following 13 lines:
17* 'ScaLAPACK BRD computation input file'
18* 'PVM machine'
19* 'BRD.out' output file name
20* 6 device out
21* 3 number of problems sizes
22* 16 20 18 values of M
23* 16 18 20 values of N
24* 3 number of NB's
25* 2 3 5 values of NB
26* 7 number of process grids (ordered pairs of P & Q)
27* 1 2 1 4 2 3 8 values of P
28* 1 2 4 1 3 2 1 values of Q
29* 1.0 threshold
30*
31*
32* Internal Parameters
33* ===================
34*
35* TOTMEM INTEGER, default = 2000000
36* TOTMEM is a machine-specific parameter indicating the
37* maximum amount of available memory in bytes.
38* The user should customize TOTMEM to his platform. Remember
39* to leave room in memory for the operating system, the BLACS
40* buffer, etc. For example, on a system with 8 MB of memory
41* per process (e.g., one processor on an Intel iPSC/860), the
42* parameters we use are TOTMEM=6200000 (leaving 1.8 MB for OS,
43* code, BLACS buffer, etc). However, for PVM, we usually set
44* TOTMEM = 2000000. Some experimenting with the maximum value
45* of TOTMEM may be required.
46*
47* INTGSZ INTEGER, default = 4 bytes.
48* ZPLXSZ INTEGER, default = 16 bytes.
49* INTGSZ and ZPLXSZ indicate the length in bytes on the
50* given platform for an integer and a double precision
51* complex.
52* MEM COMPLEX*16 array, dimension ( TOTMEM / ZPLXSZ )
53*
54* All arrays used by SCALAPACK routines are allocated from
55* this array and referenced by pointers. The integer IPA,
56* for example, is a pointer to the starting element of MEM for
57* the matrix A.
58*
59* =====================================================================
60*
61* .. Parameters ..
62 INTEGER block_cyclic_2d, csrc_, ctxt_, dlen_, dtype_,
63 $ lld_, mb_, m_, nb_, n_, rsrc_
64 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
65 $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
66 $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
67 INTEGER memsiz, ntests, totmem, zplxsz, dblesz
68 COMPLEX*16 padval
69 parameter( totmem = 2000000, zplxsz = 16, dblesz = 8,
70 $ memsiz = totmem / zplxsz, ntests = 20,
71 $ padval = ( -9923.0d+0, -9923.0d+0 ) )
72* ..
73* .. Local Scalars ..
74 LOGICAL check
75 CHARACTER*6 passed
76 CHARACTER*80 OUTFILE
77 INTEGER i, iam, iaseed, ictxt, imidpad, info, ipa, ipd,
78 $ ipe, ipostpad, iprepad, iptp, iptq, ipw, j, k,
79 $ kfail, kpass, kskip, ktests, lwork, m, MAXMN,
80 $ minmn, mnp, mnq, mp, mycol, myrow, n, nb,
81 $ ndiag, ngrids, NMAT, nnb, noffd, nout, npcol,
82 $ nprocs, nprow, nq, workbrd, worksiz
83 REAL thresh
84 DOUBLE PRECISION anorm, FRESID, nops, tmflops
85* ..
86* .. Local Arrays ..
87 INTEGER desca( dlen_ ), ierr( 1 ), nbval( ntests ),
88 $ mval( ntests ), nval( ntests ),
89 $ pval( ntests ), qval( ntests )
90 DOUBLE PRECISION ctime( 1 ), WTIME( 1 )
91 COMPLEX*16 mem( memsiz )
92* ..
93* .. External Subroutines ..
94 EXTERNAL blacs_barrier, blacs_exit, blacs_get,
96 $ blacs_pinfo, descinit, igsum2d, pzchekpad,
100* ..
101* .. External Functions ..
102 INTEGER iceil, numroc
103 DOUBLE PRECISION pzlange
104 EXTERNAL iceil, numroc, pzlange
105* ..
106* .. Intrinsic Functions ..
107 INTRINSIC dble, max, min
108* ..
109* .. Data statements ..
110 DATA ktests, kpass, kfail, kskip / 4*0 /
111* ..
112* .. Executable Statements ..
113*
114* Get starting information
115*
116 CALL blacs_pinfo( iam, nprocs )
117 iaseed = 100
118 CALL pzbrdinfo( outfile, nout, nmat, mval, ntests, nval, ntests,
119 $ nnb, nbval, ntests, ngrids, pval, ntests, qval,
120 $ ntests, thresh, mem, iam, nprocs )
121 check = ( thresh.GE.0.0e+0 )
122*
123* Print headings
124*
125 IF( iam.EQ.0 ) THEN
126 WRITE( nout, fmt = * )
127 WRITE( nout, fmt = 9995 )
128 WRITE( nout, fmt = 9994 )
129 WRITE( nout, fmt = * )
130 END IF
131*
132* Loop over different process grids
133*
134 DO 30 i = 1, ngrids
135*
136 nprow = pval( i )
137 npcol = qval( i )
138*
139* Make sure grid information is correct
140*
141 ierr( 1 ) = 0
142 IF( nprow.LT.1 ) THEN
143 IF( iam.EQ.0 )
144 $ WRITE( nout, fmt = 9999 ) 'GRID', 'nprow', nprow
145 ierr( 1 ) = 1
146 ELSE IF( npcol.LT.1 ) THEN
147 IF( iam.EQ.0 )
148 $ WRITE( nout, fmt = 9999 ) 'GRID', 'npcol', npcol
149 ierr( 1 ) = 1
150 ELSE IF( nprow*npcol.GT.nprocs ) THEN
151 IF( iam.EQ.0 )
152 $ WRITE( nout, fmt = 9998 ) nprow*npcol, nprocs
153 ierr( 1 ) = 1
154 END IF
155*
156 IF( ierr( 1 ).GT.0 ) THEN
157 IF( iam.EQ.0 )
158 $ WRITE( nout, fmt = 9997 ) 'grid'
159 kskip = kskip + 1
160 GO TO 30
161 END IF
162*
163* Define process grid
164*
165 CALL blacs_get( -1, 0, ictxt )
166 CALL blacs_gridinit( ictxt, 'Row-major', nprow, npcol )
167 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
168*
169 IF( myrow.GE.nprow .OR. mycol.GE.npcol )
170 $ GO TO 30
171*
172* Go to bottom of loop if this case doesn't use my process
173*
174 DO 20 j = 1, nmat
175*
176 m = mval( j )
177 n = nval( j )
178*
179* Make sure matrix information is correct
180*
181 ierr( 1 ) = 0
182 IF( m.LT.1 ) THEN
183 IF( iam.EQ.0 )
184 $ WRITE( nout, fmt = 9999 ) 'MATRIX', 'M', m
185 ierr( 1 ) = 1
186 ELSE IF( n.LT.1 ) THEN
187 IF( iam.EQ.0 )
188 $ WRITE( nout, fmt = 9999 ) 'MATRIX', 'N', n
189 ierr( 1 ) = 1
190 END IF
191*
192* Make sure no one had error
193*
194 CALL igsum2d( ictxt, 'All', ' ', 1, 1, ierr, 1, -1, 0 )
195*
196 IF( ierr( 1 ).GT.0 ) THEN
197 IF( iam.EQ.0 )
198 $ WRITE( nout, fmt = 9997 ) 'matrix'
199 kskip = kskip + 1
200 GO TO 20
201 END IF
202*
203* Loop over different blocking sizes
204*
205 DO 10 k = 1, nnb
206*
207 nb = nbval( k )
208*
209* Make sure nb is legal
210*
211 ierr( 1 ) = 0
212 IF( nb.LT.1 ) THEN
213 ierr( 1 ) = 1
214 IF( iam.EQ.0 )
215 $ WRITE( nout, fmt = 9999 ) 'nb', 'nb', NB
216 END IF
217*
218* Check all processes for an error
219*
220 CALL IGSUM2D( ICTXT, 'all', ' ', 1, 1, IERR, 1, -1, 0 )
221*
222.GT. IF( IERR( 1 )0 ) THEN
223.EQ. IF( IAM0 )
224 $ WRITE( NOUT, FMT = 9997 ) 'nb'
225 KSKIP = KSKIP + 1
226 GO TO 10
227 END IF
228*
229* Padding constants
230*
231 MP = NUMROC( M, NB, MYROW, 0, NPROW )
232 NQ = NUMROC( N, NB, MYCOL, 0, NPCOL )
233 MNP = NUMROC( MIN( M, N ), NB, MYROW, 0, NPROW )
234 MNQ = NUMROC( MIN( M, N ), NB, MYCOL, 0, NPCOL )
235 IF( CHECK ) THEN
236 IPREPAD = MAX( NB, MP )
237 IMIDPAD = NB
238 IPOSTPAD = MAX( NB, NQ )
239 ELSE
240 IPREPAD = 0
241 IMIDPAD = 0
242 IPOSTPAD = 0
243 END IF
244*
245* Initialize the array descriptor for the matrix A
246*
247 CALL DESCINIT( DESCA, M, N, NB, NB, 0, 0, ICTXT,
248 $ MAX( 1, MP )+IMIDPAD, IERR( 1 ) )
249*
250 CALL IGSUM2D( ICTXT, 'all', ' ', 1, 1, IERR, 1, -1, 0 )
251*
252.LT. IF( IERR( 1 )0 ) THEN
253.EQ. IF( IAM0 )
254 $ WRITE( NOUT, FMT = 9997 ) 'descriptor'
255 KSKIP = KSKIP + 1
256 GO TO 10
257 END IF
258*
259* Assign pointers into MEM for SCALAPACK arrays, A is
260* allocated starting at position MEM( IPREPAD+1 )
261*
262.GE. IF( MN ) THEN
263 NDIAG = MNQ
264 NOFFD = MNP
265 NDIAG = ICEIL( DBLESZ*NDIAG, ZPLXSZ )
266 NOFFD = ICEIL( DBLESZ*NOFFD, ZPLXSZ )
267 ELSE
268 NDIAG = MNP
269 NOFFD = NUMROC( MIN( M, N )-1, NB, MYCOL, 0, NPCOL )
270 NDIAG = ICEIL( DBLESZ*NDIAG, ZPLXSZ )
271 NOFFD = ICEIL( DBLESZ*NOFFD, ZPLXSZ )
272 END IF
273*
274 IPA = IPREPAD + 1
275 IPD = IPA + DESCA( LLD_ )*NQ + IPOSTPAD + IPREPAD
276 IPE = IPD + NDIAG + IPOSTPAD + IPREPAD
277 IPTQ = IPE + NOFFD + IPOSTPAD + IPREPAD
278 IPTP = IPTQ + MNQ + IPOSTPAD + IPREPAD
279 IPW = IPTP + MNP + IPOSTPAD + IPREPAD
280*
281* Calculate the amount of workspace required for the
282* reduction
283*
284 LWORK = NB*( MP+NQ+1 ) + NQ
285 WORKBRD = LWORK + IPOSTPAD
286 WORKSIZ = WORKBRD
287*
288* Figure the amount of workspace required by the check
289*
290 IF( CHECK ) THEN
291 WORKSIZ = MAX( LWORK, 2*NB*( MP+NQ+NB ) ) + IPOSTPAD
292 END IF
293*
294* Check for adequate memory for problem size
295*
296 IERR( 1 ) = 0
297.GT. IF( IPW+WORKSIZMEMSIZ ) THEN
298.EQ. IF( IAM0 )
299 $ WRITE( NOUT, FMT = 9996 ) 'bidiagonal reduction',
300 $ ( IPW+WORKSIZ )*ZPLXSZ
301 IERR( 1 ) = 1
302 END IF
303*
304* Check all processes for an error
305*
306 CALL IGSUM2D( ICTXT, 'all', ' ', 1, 1, IERR, 1, -1, 0 )
307*
308.GT. IF( IERR( 1 )0 ) THEN
309.EQ. IF( IAM0 )
310 $ WRITE( NOUT, FMT = 9997 ) 'memory'
311 KSKIP = KSKIP + 1
312 GO TO 10
313 END IF
314*
315* Generate the matrix A
316*
317 CALL PZMATGEN( ICTXT, 'no', 'no', DESCA( M_ ),
318 $ DESCA( N_ ), DESCA( MB_ ), DESCA( NB_ ),
319 $ MEM( IPA ), DESCA( LLD_ ), DESCA( RSRC_ ),
320 $ DESCA( CSRC_ ), IASEED, 0, MP, 0, NQ,
321 $ MYROW, MYCOL, NPROW, NPCOL )
322*
323* Need Infinity-norm of A for checking
324*
325 IF( CHECK ) THEN
326 CALL PZFILLPAD( ICTXT, MP, NQ, MEM( IPA-IPREPAD ),
327 $ DESCA( LLD_ ), IPREPAD, IPOSTPAD,
328 $ PADVAL )
329 CALL PZFILLPAD( ICTXT, NDIAG, 1, MEM( IPD-IPREPAD ),
330 $ NDIAG, IPREPAD, IPOSTPAD, PADVAL )
331 CALL PZFILLPAD( ICTXT, NOFFD, 1, MEM( IPE-IPREPAD ),
332 $ NOFFD, IPREPAD, IPOSTPAD, PADVAL )
333 CALL PZFILLPAD( ICTXT, MNQ, 1, MEM( IPTQ-IPREPAD ),
334 $ MNQ, IPREPAD, IPOSTPAD, PADVAL )
335 CALL PZFILLPAD( ICTXT, MNP, 1, MEM( IPTP-IPREPAD ),
336 $ MNP, IPREPAD, IPOSTPAD, PADVAL )
337 CALL PZFILLPAD( ICTXT, WORKSIZ-IPOSTPAD, 1,
338 $ MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
339 $ IPREPAD, IPOSTPAD, PADVAL )
340 ANORM = PZLANGE( 'i', M, N, MEM( IPA ), 1, 1, DESCA,
341 $ MEM( IPW ) )
342 CALL PZCHEKPAD( ICTXT, 'pzlange', MP, NQ,
343 $ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
344 $ IPREPAD, IPOSTPAD, PADVAL )
345 CALL PZCHEKPAD( ICTXT, 'pzlange', WORKSIZ-IPOSTPAD,
346 $ 1, MEM( IPW-IPREPAD ),
347 $ WORKSIZ-IPOSTPAD, IPREPAD, IPOSTPAD,
348 $ PADVAL )
349 CALL PZFILLPAD( ICTXT, WORKBRD-IPOSTPAD, 1,
350 $ MEM( IPW-IPREPAD ), WORKBRD-IPOSTPAD,
351 $ IPREPAD, IPOSTPAD, PADVAL )
352 END IF
353*
354 CALL SLBOOT()
355 CALL BLACS_BARRIER( ICTXT, 'all' )
356 CALL SLTIMER( 1 )
357*
358* Reduce to bidiagonal form
359*
360 CALL PZGEBRD( M, N, MEM( IPA ), 1, 1, DESCA, MEM( IPD ),
361 $ MEM( IPE ), MEM( IPTQ ), MEM( IPTP ),
362 $ MEM( IPW ), LWORK, INFO )
363*
364 CALL SLTIMER( 1 )
365*
366 IF( CHECK ) THEN
367*
368* Check for memory overwrite
369*
370 CALL PZCHEKPAD( ICTXT, 'pzgebrd', MP, NQ,
371 $ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
372 $ IPREPAD, IPOSTPAD, PADVAL )
373 CALL PZCHEKPAD( ICTXT, 'pzgebrd', NDIAG, 1,
374 $ MEM( IPD-IPREPAD ), NDIAG, IPREPAD,
375 $ IPOSTPAD, PADVAL )
376 CALL PZCHEKPAD( ICTXT, 'pzgebrd', NOFFD, 1,
377 $ MEM( IPE-IPREPAD ), NOFFD, IPREPAD,
378 $ IPOSTPAD, PADVAL )
379 CALL PZCHEKPAD( ICTXT, 'pzgebrd', MNQ, 1,
380 $ MEM( IPTQ-IPREPAD ), MNQ, IPREPAD,
381 $ IPOSTPAD, PADVAL )
382 CALL PZCHEKPAD( ICTXT, 'pzgebrd', MNP, 1,
383 $ MEM( IPTP-IPREPAD ), MNP, IPREPAD,
384 $ IPOSTPAD, PADVAL )
385 CALL PZCHEKPAD( ICTXT, 'pzgebrd', WORKBRD-IPOSTPAD,
386 $ 1, MEM( IPW-IPREPAD ),
387 $ WORKBRD-IPOSTPAD, IPREPAD,
388 $ IPOSTPAD, PADVAL )
389 CALL PZFILLPAD( ICTXT, WORKSIZ-IPOSTPAD, 1,
390 $ MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
391 $ IPREPAD, IPOSTPAD, PADVAL )
392*
393* Compute fctres = ||A-Q*B*P|| / (||A|| * N * eps)
394*
395 CALL PZGEBDRV( M, N, MEM( IPA ), 1, 1, DESCA,
396 $ MEM( IPD ), MEM( IPE ), MEM( IPTQ ),
397 $ MEM( IPTP ), MEM( IPW ), IERR( 1 ) )
398 CALL PZLAFCHK( 'no', 'no', M, N, MEM( IPA ), 1, 1,
399 $ DESCA, IASEED, ANORM, FRESID,
400 $ MEM( IPW ) )
401*
402* Check for memory overwrite
403*
404 CALL PZCHEKPAD( ICTXT, 'pzgebdrv', MP, NQ,
405 $ MEM( IPA-IPREPAD ), DESCA( LLD_ ),
406 $ IPREPAD, IPOSTPAD, PADVAL )
407 CALL PZCHEKPAD( ICTXT, 'pzgebdrv', NDIAG, 1,
408 $ MEM( IPD-IPREPAD ), NDIAG, IPREPAD,
409 $ IPOSTPAD, PADVAL )
410 CALL PZCHEKPAD( ICTXT, 'pzgebdrv', NOFFD, 1,
411 $ MEM( IPE-IPREPAD ), NOFFD, IPREPAD,
412 $ IPOSTPAD, PADVAL )
413 CALL PZCHEKPAD( ICTXT, 'pzgebdrv', WORKSIZ-IPOSTPAD,
414 $ 1, MEM( IPW-IPREPAD ),
415 $ WORKSIZ-IPOSTPAD, IPREPAD,
416 $ IPOSTPAD, PADVAL )
417*
418* Test residual and detect NaN result
419*
420.LE..AND..EQ. IF( FRESIDTHRESH FRESID-FRESID0.0D+0
421.AND..EQ. $ IERR( 1 )0 ) THEN
422 KPASS = KPASS + 1
423 PASSED = 'passed'
424 ELSE
425.EQ..AND..EQ. IF( MYROW0 MYCOL0 )
426 $ WRITE( NOUT, FMT = 9986 ) FRESID
427*
428 KFAIL = KFAIL + 1
429 PASSED = 'failed'
430 END IF
431*
432.EQ..AND..EQ..AND..NE. IF( MYROW0 MYCOL0 IERR( 1 )0 )
433 $ WRITE( NOUT, FMT = * )
434 $ 'd or e copies incorrect ...'
435 ELSE
436*
437* Don't perform the checking, only the timing operation
438*
439 KPASS = KPASS + 1
440 FRESID = FRESID - FRESID
441 PASSED = 'bypass'
442*
443 END IF
444*
445* Gather maximum of all CPU and WALL clock timings
446*
447 CALL SLCOMBINE( ICTXT, 'all', '>', 'w', 1, 1, WTIME )
448 CALL SLCOMBINE( ICTXT, 'all', '>', 'c', 1, 1, CTIME )
449*
450* Print results
451*
452.EQ..AND..EQ. IF( MYROW0 MYCOL0 ) THEN
453*
454* BRD requires 32/3 N^3 floating point operations
455*
456 MAXMN = MAX( M, N )
457 MINMN = MIN( M, N )
458 NOPS = 16.0D+0 * DBLE( MINMN ) * DBLE( MINMN ) *
459 $ ( DBLE( MAXMN ) - DBLE( MINMN ) / 3.0D+0 )
460 NOPS = NOPS / 1.0D+6
461*
462* Print WALL time
463*
464.GT. IF( WTIME( 1 )0.0D+0 ) THEN
465 TMFLOPS = NOPS / WTIME( 1 )
466 ELSE
467 TMFLOPS = 0.0D+0
468 END IF
469.GE. IF( WTIME( 1 )0.0D+0 )
470 $ WRITE( NOUT, FMT = 9993 ) 'wall', M, N, NB, NPROW,
471 $ NPCOL, WTIME( 1 ), TMFLOPS, FRESID, PASSED
472*
473* Print CPU time
474*
475.GT. IF( CTIME( 1 )0.0D+0 ) THEN
476 TMFLOPS = NOPS / CTIME( 1 )
477 ELSE
478 TMFLOPS = 0.0D+0
479 END IF
480.GE. IF( CTIME( 1 )0.0D+0 )
481 $ WRITE( NOUT, FMT = 9993 ) 'cpu ', M, N, NB, NPROW,
482 $ NPCOL, CTIME( 1 ), TMFLOPS, FRESID, PASSED
483 END IF
484 10 CONTINUE
485 20 CONTINUE
486*
487 CALL BLACS_GRIDEXIT( ICTXT )
488 30 CONTINUE
489*
490* Print ending messages and close output file
491*
492.EQ. IF( IAM0 ) THEN
493 KTESTS = KPASS + KFAIL + KSKIP
494 WRITE( NOUT, FMT = * )
495 WRITE( NOUT, FMT = 9992 ) KTESTS
496 IF( CHECK ) THEN
497 WRITE( NOUT, FMT = 9991 ) KPASS
498 WRITE( NOUT, FMT = 9989 ) KFAIL
499 ELSE
500 WRITE( NOUT, FMT = 9990 ) KPASS
501 END IF
502 WRITE( NOUT, FMT = 9988 ) KSKIP
503 WRITE( NOUT, FMT = * )
504 WRITE( NOUT, FMT = * )
505 WRITE( NOUT, FMT = 9987 )
506.NE..AND..NE. IF( NOUT6 NOUT0 ) CLOSE ( NOUT )
507 END IF
508*
509 CALL BLACS_EXIT( 0 )
510*
511 9999 FORMAT( 'illegal ', A6, ': ', A5, ' = ', I3,
512 $ '; it should be at least 1' )
513 9998 FORMAT( 'illegal grid: nprow*npcol = ', i4, '. It can be at most',
514 $ i4 )
515 9997 FORMAT( 'Bad ', a6, ' parameters: going on to next test case.' )
516 9996 FORMAT( 'Unable to perform ', a, ': need TOTMEM of at least',
517 $ i11 )
518 9995 FORMAT( 'TIME M N NB P Q BRD Time ',
519 $ ' MFLOPS Residual CHECK' )
520 9994 FORMAT( '---- ------ ------ --- ----- ----- --------- ',
521 $ '----------- -------- ------' )
522 9993 FORMAT( a4, 1x, i6, 1x, i6, 1x, i3, 1x, i5, 1x, i5, 1x, f9.2, 1x,
523 $ f11.2, 1x, f8.2, 1x, a6 )
524 9992 FORMAT( 'Finished', i4, ' tests, with the following results:' )
525 9991 FORMAT( i5, ' tests completed and passed residual checks.' )
526 9990 FORMAT( i5, ' tests completed without checking.' )
527 9989 FORMAT( i5, ' tests completed and failed residual checks.' )
528 9988 FORMAT( i5, ' tests skipped because of illegal input values.' )
529 9987 FORMAT( 'END OF TESTS.' )
530 9986 FORMAT( '||A - Q*B*P|| / (||A|| * N * eps) = ', g25.7 )
531*
532 stop
533*
534* End of PZBRDDRIVER
535*
536 END
subroutine pzlafchk(aform, diag, m, n, a, ia, ja, desca, iaseed, anorm, fresid, work)
Definition pzlafchk.f:3
subroutine pzmatgen(ictxt, aform, diag, m, n, mb, nb, a, lda, iarow, iacol, iseed, iroff, irnum, icoff, icnum, myrow, mycol, nprow, npcol)
Definition pzmatgen.f:4
integer function iceil(inum, idenom)
Definition iceil.f:2
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine blacs_gridinit(cntxt, c, nprow, npcol)
Definition mpi.f:745
double precision function pzlange(norm, m, n, a, ia, ja, desca, work)
Definition mpi.f:1287
subroutine descinit(desc, m, n, mb, nb, irsrc, icsrc, ictxt, lld, info)
Definition mpi.f:777
subroutine blacs_gridexit(cntxt)
Definition mpi.f:762
subroutine blacs_gridinfo(cntxt, nprow, npcol, myrow, mycol)
Definition mpi.f:754
integer function numroc(n, nb, iproc, isrcproc, nprocs)
Definition mpi.f:786
program pzbrddriver
Definition pzbrddriver.f:1
subroutine pzbrdinfo(summry, nout, nmat, mval, ldmval, nval, ldnval, nnb, nbval, ldnbval, ngrids, pval, ldpval, qval, ldqval, thresh, work, iam, nprocs)
Definition pzbrdinfo.f:5
subroutine pzchekpad(ictxt, mess, m, n, a, lda, ipre, ipost, chkval)
Definition pzchekpad.f:3
subroutine pzfillpad(ictxt, m, n, a, lda, ipre, ipost, chkval)
Definition pzfillpad.f:2
subroutine pzgebdrv(m, n, a, ia, ja, desca, d, e, tauq, taup, work, info)
Definition pzgebdrv.f:3
subroutine pzgebrd(m, n, a, ia, ja, desca, d, e, tauq, taup, work, lwork, info)
Definition pzgebrd.f:3
subroutine slboot()
Definition sltimer.f:2
subroutine sltimer(i)
Definition sltimer.f:47
subroutine slcombine(ictxt, scope, op, timetype, n, ibeg, times)
Definition sltimer.f:267