SCALAPACK 2.2.2
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ pshrdinfo()

subroutine pshrdinfo ( character*( * )  summry,
integer  nout,
integer  nmat,
integer, dimension( ldnval )  nval,
integer, dimension( ldnval )  nvlo,
integer, dimension( ldnval )  nvhi,
integer  ldnval,
integer  nnb,
integer, dimension( ldnbval )  nbval,
integer  ldnbval,
integer  ngrids,
integer, dimension( ldpval )  pval,
integer  ldpval,
integer, dimension( ldqval )  qval,
integer  ldqval,
real  thresh,
integer, dimension( * )  work,
integer  iam,
integer  nprocs 
)

Definition at line 1 of file pshrdinfo.f.

5*
6* -- ScaLAPACK routine (version 1.7) --
7* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
8* and University of California, Berkeley.
9* May 1, 1997
10*
11* .. Scalar Arguments ..
12 INTEGER IAM, LDNBVAL, LDNVAL, LDPVAL, LDQVAL,
13 $ NGRIDS, NMAT, NNB, NOUT, NPROCS
14 REAL THRESH
15* ..
16* .. Array Arguments ..
17 CHARACTER*( * ) SUMMRY
18 INTEGER NBVAL( LDNBVAL ), NVAL( LDNVAL ),
19 $ NVHI( LDNVAL ), NVLO( LDNVAL ),
20 $ PVAL( LDPVAL ), QVAL( LDQVAL ), WORK( * )
21* ..
22*
23* Purpose
24* =======
25*
26* PSHRDINFO get the needed startup information for the Hessenberg
27* reduction tests and transmits it to all processes.
28*
29* Arguments
30* =========
31*
32* SUMMRY (global output) CHARACTER*(*)
33* Name of output (summary) file (if any). Only defined for
34* process 0.
35*
36* NOUT (global output) INTEGER
37* The unit number for output file. NOUT = 6, output to screen,
38* NOUT = 0, output to stderr. Only defined for process 0.
39*
40* NMAT (global output) INTEGER
41* The number of different values that can be used for
42* N, IHI & ILO.
43*
44* NVAL (global output) INTEGER array, dimension (LDNVAL)
45* The values of N (number of rows & columns in matrix).
46*
47* NVLO (global output) INTEGER array, dimension (LDNVAL)
48* The values of ILO.
49*
50* NVHI (global output) INTEGER array, dimension (LDNVAL)
51* The values of IHI.
52*
53* LDNVAL (global input) INTEGER
54* The maximum number of different values that can be used for
55* N, ILO and IHI. LDNVAL >= NMAT.
56*
57* NNB (global output) INTEGER
58* The number of different values that can be used for NB.
59*
60* NBVAL (global output) INTEGER array, dimension (LDNBVAL)
61* The values of NB (blocksize) to run the code with.
62*
63* LDNBVAL (global input) INTEGER
64* The maximum number of different values that can be used for
65* NB, LDNBVAL >= NNB.
66*
67* NGRIDS (global output) INTEGER
68* The number of different values that can be used for P & Q.
69*
70* PVAL (global output) INTEGER array, dimension (LDPVAL)
71* The values of P (number of process rows) to run the code
72* with.
73*
74* LDPVAL (global input) INTEGER
75* The maximum number of different values that can be used for
76* P, LDPVAL >= NGRIDS.
77*
78* QVAL (global output) INTEGER array, dimension (LDQVAL)
79* The values of Q (number of process columns) to run the code
80* with.
81*
82* LDQVAL (global input) INTEGER
83* The maximum number of different values that can be used for
84* Q, LDQVAL >= NGRIDS.
85*
86* THRESH (global output) REAL
87* Indicates what error checks shall be run and printed out:
88* = 0 : Perform no error checking
89* > 0 : report all residuals greater than THRESH.
90*
91* WORK (local workspace) INTEGER array, dimension >=
92* 3*LDNVAL+LDNBVAL+2*LDPVAL. Used to pack all input arrays
93* in order to send info in one message.
94*
95* IAM (local input) INTEGER
96* My process number.
97*
98* NPROCS (global input) INTEGER
99* The total number of processes.
100*
101* Note
102* ====
103*
104* For packing the information we assumed that the length in bytes of an
105* integer is equal to the length in bytes of a real single precision.
106*
107* =====================================================================
108*
109* .. Parameters ..
110 INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
111 $ LLD_, MB_, M_, NB_, N_, RSRC_
112 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
113 $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
114 $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
115 INTEGER NIN
116 parameter( nin = 11 )
117* ..
118* .. Local Scalars ..
119 CHARACTER*79 USRINFO
120 INTEGER I, ICTXT
121 REAL EPS
122* ..
123* .. External Subroutines ..
124 EXTERNAL blacs_abort, blacs_get, blacs_gridexit,
125 $ blacs_gridinit, blacs_setup, icopy, igebr2d,
126 $ igebs2d, sgebr2d, sgebs2d
127* ..
128* .. External Functions ..
129 REAL PSLAMCH
130 EXTERNAL pslamch
131* ..
132* .. Intrinsic Functions ..
133 INTRINSIC max, min
134* ..
135* .. Executable Statements ..
136*
137* Process 0 reads the input data, broadcasts to other processes and
138* writes needed information to NOUT
139*
140 IF( iam.EQ.0 ) THEN
141*
142* Open file and skip data file header
143*
144 OPEN( unit = nin, file = 'HRD.dat', status = 'OLD' )
145 READ( nin, fmt = * )summry
146 summry = ' '
147*
148* Read in user-supplied info about machine type, compiler, etc.
149*
150 READ( nin, fmt = * ) usrinfo
151*
152* Read name and unit number for summary output file
153*
154 READ( nin, fmt = * ) summry
155 READ( nin, fmt = * ) nout
156 IF( nout.NE.0 .AND. nout.NE.6 )
157 $ OPEN( unit = nout, file = summry, status = 'UNKNOWN' )
158*
159* Read and check the parameter values for the tests.
160*
161* Get number of matrices
162*
163 READ( nin, fmt = * ) nmat
164 IF( nmat.LT.1. .OR. nmat.GT.ldnval ) THEN
165 WRITE( nout, fmt = 9997 ) 'N', ldnval
166 GO TO 20
167 END IF
168*
169* Get values of N, ILO, IHI
170*
171 READ( nin, fmt = * ) ( nval( i ), i = 1, nmat )
172 READ( nin, fmt = * ) ( nvlo( i ), i = 1, nmat )
173 READ( nin, fmt = * ) ( nvhi( i ), i = 1, nmat )
174*
175* Get values of NB
176*
177 READ( nin, fmt = * ) nnb
178 IF( nnb.LT.1 .OR. nnb.GT.ldnbval ) THEN
179 WRITE( nout, fmt = 9997 ) 'NB', ldnbval
180 GO TO 20
181 END IF
182 READ( nin, fmt = * ) ( nbval( i ), i = 1, nnb )
183*
184* Get number of grids
185*
186 READ( nin, fmt = * ) ngrids
187 IF( ngrids.LT.1 .OR. ngrids.GT.ldpval ) THEN
188 WRITE( nout, fmt = 9997 ) 'Grids', ldpval
189 GO TO 20
190 ELSE IF( ngrids.GT.ldqval ) THEN
191 WRITE( nout, fmt = 9997 ) 'Grids', ldqval
192 GO TO 20
193 END IF
194*
195* Get values of P and Q
196*
197 READ( nin, fmt = * ) ( pval( i ), i = 1, ngrids )
198 READ( nin, fmt = * ) ( qval( i ), i = 1, ngrids )
199*
200* Get level of checking
201*
202 READ( nin, fmt = * ) thresh
203*
204* Close input file
205*
206 CLOSE( nin )
207*
208* For pvm only: if virtual machine not set up, allocate it and
209* spawn the correct number of processes.
210*
211 IF( nprocs.LT.1 ) THEN
212 nprocs = 0
213 DO 10 i = 1, ngrids
214 nprocs = max( nprocs, pval( i )*qval( i ) )
215 10 CONTINUE
216 CALL blacs_setup( iam, nprocs )
217 END IF
218*
219* Temporarily define blacs grid to include all processes so
220* information can be broadcast to all processes
221*
222 CALL blacs_get( -1, 0, ictxt )
223 CALL blacs_gridinit( ictxt, 'Row-major', 1, nprocs )
224*
225* Compute machine epsilon
226*
227 eps = pslamch( ictxt, 'eps' )
228*
229* Pack information arrays and broadcast
230*
231 CALL sgebs2d( ictxt, 'All', ' ', 1, 1, thresh, 1 )
232*
233 work( 1 ) = nmat
234 work( 2 ) = nnb
235 work( 3 ) = ngrids
236 CALL igebs2d( ictxt, 'All', ' ', 1, 3, work, 1 )
237*
238 i = 1
239 CALL icopy( nmat, nval, 1, work( i ), 1 )
240 i = i + nmat
241 CALL icopy( nmat, nvlo, 1, work( i ), 1 )
242 i = i + nmat
243 CALL icopy( nmat, nvhi, 1, work( i ), 1 )
244 i = i + nmat
245 CALL icopy( nnb, nbval, 1, work( i ), 1 )
246 i = i + nnb
247 CALL icopy( ngrids, pval, 1, work( i ), 1 )
248 i = i + ngrids
249 CALL icopy( ngrids, qval, 1, work( i ), 1 )
250 i = i + ngrids -1
251 CALL igebs2d( ictxt, 'All', ' ', 1, i, work, 1 )
252*
253* regurgitate input
254*
255 WRITE( nout, fmt = 9999 )
256 $ 'ScaLAPACK Reduction routine to Hessenberg form.'
257 WRITE( nout, fmt = 9999 ) usrinfo
258 WRITE( nout, fmt = * )
259 WRITE( nout, fmt = 9999 )
260 $ 'Tests of the parallel '//
261 $ 'real single precision Hessenberg '
262 WRITE( nout, fmt = 9999 ) 'reduction routines.'
263 WRITE( nout, fmt = 9999 )
264 $ 'The following scaled residual '//
265 $ 'checks will be computed:'
266 WRITE( nout, fmt = 9999 )
267 $ ' ||A - Q H Q''|| / (||A|| * eps * N)'
268 WRITE( nout, fmt = 9999 )
269 $ 'The matrix A is randomly '//
270 $ 'generated for each test.'
271 WRITE( nout, fmt = * )
272 WRITE( nout, fmt = 9999 )
273 $ 'An explanation of the input/output '//
274 $ 'parameters follows:'
275 WRITE( nout, fmt = 9999 )
276 $ 'TIME : Indicates whether WALL or '//
277 $ 'CPU time was used.'
278 WRITE( nout, fmt = 9999 )
279 $ 'N : The number of rows and columns '//
280 $ 'of the matrix A.'
281 WRITE( nout, fmt = 9999 )
282 $ 'NB : The size of the square blocks'//
283 $ ' the matrix A is split into.'
284 WRITE( nout, fmt = 9999 )
285 $ ' on to the next column of processes.'
286 WRITE( nout, fmt = 9999 )
287 $ 'P : The number of process rows.'
288 WRITE( nout, fmt = 9999 )
289 $ 'Q : The number of process columns.'
290 WRITE( nout, fmt = 9999 )
291 $ 'HRD time : Time in seconds to compute HRD '
292 WRITE( nout, fmt = 9999 )
293 $ 'MFLOPS : Rate of execution for HRD ' //
294 $ 'reduction.'
295 WRITE( nout, fmt = * )
296 WRITE( nout, fmt = 9999 )
297 $ 'The following parameter values will be used:'
298 WRITE( nout, fmt = 9995 )
299 $ 'N ', ( nval( i ), i = 1, min( nmat, 10 ) )
300 IF( nmat.GT.10 )
301 $ WRITE( nout, fmt = 9994 ) ( nval( i ), i = 11, nmat )
302 WRITE( nout, fmt = 9995 )
303 $ 'ILO ', ( nvlo( i ), i = 1, min( nmat, 10 ) )
304 IF( nmat.GT.10 )
305 $ WRITE( nout, fmt = 9994 ) ( nvlo( i ), i = 11, nmat )
306 WRITE( nout, fmt = 9995 )
307 $ 'IHI ', ( nvhi( i ), i = 1, min( nmat, 10 ) )
308 IF( nmat.GT.10 )
309 $ WRITE( nout, fmt = 9994 ) ( nvhi( i ), i = 11, nmat )
310 WRITE( nout, fmt = 9995 )
311 $ 'NB ', ( nbval( i ), i = 1, min( nnb, 10 ) )
312 IF( nnb.GT.10 )
313 $ WRITE( nout, fmt = 9994 ) ( nbval( i ), i = 11, nnb )
314 WRITE( nout, fmt = 9995 )
315 $ 'P ', ( pval( i ), i = 1, min( ngrids, 10 ) )
316 IF( ngrids.GT.10 )
317 $ WRITE( nout, fmt = 9994 ) ( pval( i ), i = 11, ngrids )
318 WRITE( nout, fmt = 9995 )
319 $ 'Q ', ( qval( i ), i = 1, min( ngrids, 10 ) )
320 IF( ngrids.GT.10 )
321 $ WRITE( nout, fmt = 9994 ) ( qval( i ), i = 11, ngrids )
322 WRITE( nout, fmt = * )
323 WRITE( nout, fmt = 9996 ) eps
324 WRITE( nout, fmt = 9993 ) thresh
325*
326 ELSE
327*
328* If in pvm, must participate setting up virtual machine
329*
330 IF( nprocs.LT.1 )
331 $ CALL blacs_setup( iam, nprocs )
332*
333* Temporarily define blacs grid to include all processes so
334* all processes have needed startup information
335*
336 CALL blacs_get( -1, 0, ictxt )
337 CALL blacs_gridinit( ictxt, 'Row-major', 1, nprocs )
338*
339* Compute machine epsilon
340*
341 eps = pslamch( ictxt, 'eps' )
342*
343 CALL sgebr2d( ictxt, 'All', ' ', 1, 1, thresh, 1, 0, 0 )
344 CALL igebr2d( ictxt, 'All', ' ', 1, 3, work, 1, 0, 0 )
345 nmat = work( 1 )
346 nnb = work( 2 )
347 ngrids = work( 3 )
348*
349 i = 3*nmat + nnb + 2*ngrids
350 CALL igebr2d( ictxt, 'All', ' ', 1, i, work, 1, 0, 0 )
351*
352 i = 1
353 CALL icopy( nmat, work( i ), 1, nval, 1 )
354 i = i + nmat
355 CALL icopy( nmat, work( i ), 1, nvlo, 1 )
356 i = i + nmat
357 CALL icopy( nmat, work( i ), 1, nvhi, 1 )
358 i = i + nmat
359 CALL icopy( nnb, work( i ), 1, nbval, 1 )
360 i = i + nnb
361 CALL icopy( ngrids, work( i ), 1, pval, 1 )
362 i = i + ngrids
363 CALL icopy( ngrids, work( i ), 1, qval, 1 )
364*
365 END IF
366*
367 CALL blacs_gridexit( ictxt )
368*
369 RETURN
370*
371 20 CONTINUE
372 WRITE( nout, fmt = 9998 )
373 CLOSE( nin )
374 IF( nout.NE.6 .AND. nout.NE.0 )
375 $ CLOSE( nout )
376 CALL blacs_abort( ictxt, 1 )
377*
378 stop
379*
380 9999 FORMAT( a )
381 9998 FORMAT( ' ILLEGAL INPUT IN FILE ', 40a, '. ABORTING RUN.' )
382 9997 FORMAT( ' NUMBER OF VALUES OF ', 5a,
383 $ ' IS LESS THAN 1 OR GREATER ', 'THAN ', i2 )
384 9996 FORMAT( 'Relative machine precision (eps) is taken to be ',
385 $ e18.6 )
386 9995 FORMAT( 2x, a5, ': ', 10i6 )
387 9994 FORMAT( ' ', 10i6 )
388 9993 FORMAT( 'Routines pass computational tests if scaled residual is',
389 $ ' less than ', g14.7 )
390*
391* End of PSHRDINFO
392*
subroutine icopy(n, sx, incx, sy, incy)
Definition pblastst.f:1525
real function pslamch(ictxt, cmach)
Definition pcblastst.f:7455
#define max(A, B)
Definition pcgemr.c:180
#define min(A, B)
Definition pcgemr.c:181
Here is the call graph for this function:
Here is the caller graph for this function: