Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' of scm.gforge.inria.fr:/gitroot/simgrid/simgrid
[simgrid.git] / teshsuite / smpi / mpich3-test / perf / dtpack.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2008 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 /*
7  * This code may be used to test the performance of some of the 
8  * noncontiguous datatype operations, including vector and indexed
9  * pack and unpack operations.  To simplify the use of this code for 
10  * tuning an MPI implementation, it uses no communication, just the
11  * MPI_Pack and MPI_Unpack routines.  In addition, the individual tests are
12  * in separate routines, making it easier to compare the compiler-generated
13  * code for the user (manual) pack/unpack with the code used by 
14  * the MPI implementation.  Further, to be fair to the MPI implementation,
15  * the routines are passed the source and destination buffers; this ensures
16  * that the compiler can't optimize for statically allocated buffers.
17  */
18 #include <stdio.h>
19 #include <stdlib.h>
20 #include <string.h>
21 #include "mpi.h"
22
23 /* Needed for restrict and const definitions */
24 #include "mpitestconf.h"
25
26 static int verbose = 0;
27
28 #define N_REPS 1000
29 #define THRESHOLD 0.10
30 #define VARIANCE_THRESHOLD ((THRESHOLD * THRESHOLD) / 2)
31 #define NTRIALS 10
32
33 double mean(double *list, int count);
34 double mean(double *list, int count)
35 {
36         double retval;
37         int i;
38
39         retval = 0;
40         for (i = 0; i < count; i++)
41                 retval += list[i];
42         retval /= count;
43
44         return retval;
45 }
46
47 double noise(double *list, int count);
48 double noise(double *list, int count)
49 {
50         double *margin, retval;
51         int i;
52
53         if (!(margin = malloc(count * sizeof(double)))) {
54                 printf("Unable to allocate memory\n");
55                 return -1;
56         }
57
58         for (i = 0; i < count; i++)
59                 margin[i] = list[i] / mean(list, count);
60
61         retval = 0;
62         for (i = 0; i < count; i++) {
63                 retval += ((margin[i] - 1) * (margin[i] - 1));
64         }
65         retval /= count;
66         if (retval < 0) retval = -retval;
67
68         return retval;
69 }
70
71 /* Here are the tests */
72
73 /* Test packing a vector of individual doubles */
74 /* We don't use restrict in the function args because assignments between
75    restrict pointers is not valid in C and some compilers, such as the 
76    IBM xlc compilers, flag that use as an error.*/
77 int TestVecPackDouble( int n, int stride, 
78                        double *avgTimeUser, double *avgTimeMPI,
79                        double *dest, const double *src );
80 int TestVecPackDouble( int n, int stride, 
81                        double *avgTimeUser, double *avgTimeMPI,
82                        double *dest, const double *src )
83 {
84         double *restrict d_dest;
85         const double *restrict d_src;
86         register int i, j;
87         int          rep, position;
88         double       t1, t2, t[NTRIALS];
89         MPI_Datatype vectype;
90
91         /* User code */
92         if (verbose) printf("TestVecPackDouble (USER): ");
93         for (j = 0; j < NTRIALS; j++) {
94                 t1 = MPI_Wtime();
95                 for (rep=0; rep<N_REPS; rep++) {
96                         i = n;
97                         d_dest = dest;
98                         d_src  = src;
99                         while (i--) {
100                                 *d_dest++ = *d_src;
101                                 d_src += stride;
102                         }
103                 }
104                 t2 = MPI_Wtime() - t1;
105                 t[j] = t2;
106                 if (verbose) printf("%.3f ", t[j]);
107         }
108         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
109         /* If there is too much noise, discard the test */
110         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
111                 *avgTimeUser = 0;
112                 *avgTimeMPI = 0;
113                 if (verbose)
114                         printf("Too much noise; discarding measurement\n");
115                 return 0;
116         }
117         *avgTimeUser = mean(t, NTRIALS) / N_REPS;
118
119         /* MPI Vector code */
120         MPI_Type_vector( n, 1, stride, MPI_DOUBLE, &vectype );
121         MPI_Type_commit( &vectype );
122
123         if (verbose) printf("TestVecPackDouble (MPI): ");
124         for (j = 0; j < NTRIALS; j++) {
125                 t1 = MPI_Wtime();
126                 for (rep=0; rep<N_REPS; rep++) {
127                         position = 0;
128                         MPI_Pack( (void *)src, 1, vectype, dest, n*sizeof(double),
129                                   &position, MPI_COMM_SELF );
130                 }
131                 t2 = MPI_Wtime() - t1;
132                 t[j] = t2;
133                 if (verbose) printf("%.3f ", t[j]);
134         }
135         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
136         /* If there is too much noise, discard the test */
137         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
138                 *avgTimeUser = 0;
139                 *avgTimeMPI = 0;
140                 if (verbose)
141                         printf("Too much noise; discarding measurement\n");
142         }
143         else {
144             *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
145         }
146
147         MPI_Type_free( &vectype );
148
149         return 0;
150 }
151
152 /* Test unpacking a vector of individual doubles */
153 /* See above for why restrict is not used in the function args */
154 int TestVecUnPackDouble( int n, int stride, 
155                        double *avgTimeUser, double *avgTimeMPI,
156                        double *dest, const double *src );
157 int TestVecUnPackDouble( int n, int stride, 
158                        double *avgTimeUser, double *avgTimeMPI,
159                        double *dest, const double *src )
160 {
161         double *restrict d_dest;
162         const double *restrict d_src;
163         register int i, j;
164         int          rep, position;
165         double       t1, t2, t[NTRIALS];
166         MPI_Datatype vectype;
167
168         /* User code */
169         if (verbose) printf("TestVecUnPackDouble (USER): ");
170         for (j = 0; j < NTRIALS; j++) {
171                 t1 = MPI_Wtime();
172                 for (rep=0; rep<N_REPS; rep++) {
173                         i = n;
174                         d_dest = dest;
175                         d_src  = src;
176                         while (i--) {
177                                 *d_dest = *d_src++;
178                                 d_dest += stride;
179                         }
180                 }
181                 t2 = MPI_Wtime() - t1;
182                 t[j] = t2;
183                 if (verbose) printf("%.3f ", t[j]);
184         }
185         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
186         /* If there is too much noise, discard the test */
187         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
188                 *avgTimeUser = 0;
189                 *avgTimeMPI = 0;
190                 if (verbose)
191                         printf("Too much noise; discarding measurement\n");
192                 return 0;
193         }
194         *avgTimeUser = mean(t, NTRIALS) / N_REPS;
195     
196         /* MPI Vector code */
197         MPI_Type_vector( n, 1, stride, MPI_DOUBLE, &vectype );
198         MPI_Type_commit( &vectype );
199
200         if (verbose) printf("TestVecUnPackDouble (MPI): ");
201         for (j = 0; j < NTRIALS; j++) {
202                 t1 = MPI_Wtime();
203                 for (rep=0; rep<N_REPS; rep++) {
204                         position = 0;
205                         MPI_Unpack( (void *)src, n*sizeof(double), 
206                                     &position, dest, 1, vectype, MPI_COMM_SELF );
207                 }
208                 t2 = MPI_Wtime() - t1;
209                 t[j] = t2;
210                 if (verbose) printf("%.3f ", t[j]);
211         }
212         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
213         /* If there is too much noise, discard the test */
214         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
215                 *avgTimeUser = 0;
216                 *avgTimeMPI = 0;
217                 if (verbose)
218                         printf("Too much noise; discarding measurement\n");
219         }
220         else {
221             *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
222         }
223
224         MPI_Type_free( &vectype );
225
226         return 0;
227 }
228
229 /* Test packing a vector of 2-individual doubles */
230 /* See above for why restrict is not used in the function args */
231 int TestVecPack2Double( int n, int stride, 
232                         double *avgTimeUser, double *avgTimeMPI,
233                         double *dest, const double *src );
234 int TestVecPack2Double( int n, int stride, 
235                         double *avgTimeUser, double *avgTimeMPI,
236                         double *dest, const double *src )
237 {
238         double *restrict d_dest;
239         const double *restrict d_src;
240         register int i, j;
241         int          rep, position;
242         double       t1, t2, t[NTRIALS];
243         MPI_Datatype vectype;
244
245         /* User code */
246         if (verbose) printf("TestVecPack2Double (USER): ");
247         for (j = 0; j < NTRIALS; j++) {
248                 t1 = MPI_Wtime();
249                 for (rep=0; rep<N_REPS; rep++) {
250                         i = n;
251                         d_dest = dest;
252                         d_src  = src;
253                         while (i--) {
254                                 *d_dest++ = d_src[0];
255                                 *d_dest++ = d_src[1];
256                                 d_src += stride;
257                         }
258                 }
259                 t2 = MPI_Wtime() - t1;
260                 t[j] = t2;
261                 if (verbose) printf("%.3f ", t[j]);
262         }
263         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
264         /* If there is too much noise, discard the test */
265         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
266                 *avgTimeUser = 0;
267                 *avgTimeMPI = 0;
268                 if (verbose)
269                         printf("Too much noise; discarding measurement\n");
270                 return 0;
271         }
272         *avgTimeUser = mean(t, NTRIALS) / N_REPS;
273     
274         /* MPI Vector code */
275         MPI_Type_vector( n, 2, stride, MPI_DOUBLE, &vectype );
276         MPI_Type_commit( &vectype );
277     
278         if (verbose) printf("TestVecPack2Double (MPI): ");
279         for (j = 0; j < NTRIALS; j++) {
280                 t1 = MPI_Wtime();
281                 for (rep=0; rep<N_REPS; rep++) {
282                         position = 0;
283                         MPI_Pack( (void *)src, 1, vectype, dest, 2*n*sizeof(double),
284                                   &position, MPI_COMM_SELF );
285                 }
286                 t2 = MPI_Wtime() - t1;
287                 t[j] = t2;
288                 if (verbose) printf("%.3f ", t[j]);
289         }
290         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
291         /* If there is too much noise, discard the test */
292         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
293                 *avgTimeUser = 0;
294                 *avgTimeMPI = 0;
295                 if (verbose)
296                         printf("Too much noise; discarding measurement\n");
297         }
298         else {
299             *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
300         }
301         MPI_Type_free( &vectype );
302
303         return 0;
304 }
305
306 /* This creates an indexed type that is like a vector (for simplicity
307    of construction).  There is a possibility that the MPI implementation 
308    will recognize and simplify this (e.g., in MPI_Type_commit); if so,
309    let us know and we'll add a version that is not as regular 
310 */
311 /* See above for why restrict is not used in the function args */
312 int TestIndexPackDouble( int n, int stride, 
313                          double *avgTimeUser, double *avgTimeMPI,
314                          double *dest, const double *src );
315 int TestIndexPackDouble( int n, int stride, 
316                          double *avgTimeUser, double *avgTimeMPI,
317                          double *dest, const double *src )
318 {
319         double *restrict d_dest;
320         const double *restrict d_src;
321         register int i, j;
322         int          rep, position;
323         int          *restrict displs = 0;
324         double       t1, t2, t[NTRIALS];
325         MPI_Datatype indextype;
326
327         displs = (int *)malloc( n * sizeof(int) );
328         for (i=0; i<n; i++) displs[i] = i * stride;
329
330         /* User code */
331         if (verbose) printf("TestIndexPackDouble (USER): ");
332         for (j = 0; j < NTRIALS; j++) {
333                 t1 = MPI_Wtime();
334                 for (rep=0; rep<N_REPS; rep++) {
335                         i = n;
336                         d_dest = dest;
337                         d_src  = src;
338                         for (i=0; i<n; i++) {
339                                 *d_dest++ = d_src[displs[i]];
340                         }
341                 }
342                 t2 = MPI_Wtime() - t1;
343                 t[j] = t2;
344                 if (verbose) printf("%.3f ", t[j]);
345         }
346         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
347         /* If there is too much noise, discard the test */
348         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
349                 *avgTimeUser = 0;
350                 *avgTimeMPI = 0;
351                 if (verbose)
352                         printf("Too much noise; discarding measurement\n");
353                 return 0;
354         }
355         *avgTimeUser = mean(t, NTRIALS) / N_REPS;
356     
357         /* MPI Index code */
358         MPI_Type_create_indexed_block( n, 1, displs, MPI_DOUBLE, &indextype );
359         MPI_Type_commit( &indextype );
360
361         free( displs );
362     
363         if (verbose) printf("TestIndexPackDouble (MPI): ");
364         for (j = 0; j < NTRIALS; j++) {
365                 t1 = MPI_Wtime();
366                 for (rep=0; rep<N_REPS; rep++) {
367                         position = 0;
368                         MPI_Pack( (void *)src, 1, indextype, dest, n*sizeof(double),
369                                   &position, MPI_COMM_SELF );
370                 }
371                 t2 = MPI_Wtime() - t1;
372                 t[j] = t2;
373                 if (verbose) printf("%.3f ", t[j]);
374         }
375         if (verbose) printf("[%.3f]\n", noise(t, NTRIALS));
376         /* If there is too much noise, discard the test */
377         if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
378                 *avgTimeUser = 0;
379                 *avgTimeMPI = 0;
380                 if (verbose)
381                         printf("Too much noise; discarding measurement\n");
382         }
383         else {
384             *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
385         }
386         MPI_Type_free( &indextype );
387
388         return 0;
389 }
390
391 int Report( const char *name, const char *packname, 
392             double avgTimeMPI, double avgTimeUser );
393 int Report( const char *name, const char *packname, 
394             double avgTimeMPI, double avgTimeUser )
395 {
396         double diffTime, maxTime;
397         int errs=0;
398
399         /* Move this into a common routine */
400         diffTime = avgTimeMPI - avgTimeUser;
401         if (diffTime < 0) diffTime = - diffTime;
402         if (avgTimeMPI > avgTimeUser) maxTime = avgTimeMPI;
403         else                          maxTime = avgTimeUser;
404
405         if (verbose) {
406                 printf( "%-30s:\t%g\t%g\t(%g%%)\n", name, 
407                         avgTimeMPI, avgTimeUser,
408                         100 * (diffTime / maxTime) );
409                 fflush(stdout);
410         }
411         if (avgTimeMPI > avgTimeUser && (diffTime > THRESHOLD * maxTime)) {
412                 errs++;
413                 printf( "%s:\tMPI %s code is too slow: MPI %g\t User %g\n",
414                         name, packname, avgTimeMPI, avgTimeUser );
415         }
416
417         return errs;
418 }
419
420 /* Finally, here's the main program */
421 int main( int argc, char *argv[] )
422 {
423     int n, stride, err, errs = 0;
424     void *dest, *src;
425     double avgTimeUser, avgTimeMPI;
426
427     MPI_Init( &argc, &argv );
428     if (getenv("MPITEST_VERBOSE")) verbose = 1;
429
430     n      = 30000;
431     stride = 4;
432     dest = (void *)malloc( n * sizeof(double) );
433     src  = (void *)malloc( n * ((1+stride)*sizeof(double)) );
434     /* Touch the source and destination arrays */
435     memset( src, 0, n * (1+stride)*sizeof(double) );
436     memset( dest, 0, n * sizeof(double) );
437
438     err = TestVecPackDouble( n, stride, &avgTimeUser, &avgTimeMPI,
439                              dest, src );
440     errs += Report( "VecPackDouble", "Pack", avgTimeMPI, avgTimeUser );
441
442     err = TestVecUnPackDouble( n, stride, &avgTimeUser, &avgTimeMPI,
443                                src, dest );
444     errs += Report( "VecUnPackDouble", "Unpack", avgTimeMPI, avgTimeUser );
445
446     err = TestIndexPackDouble( n, stride, &avgTimeUser, &avgTimeMPI,
447                              dest, src );
448     errs += Report( "VecIndexDouble", "Pack", avgTimeMPI, avgTimeUser );
449
450     free(dest);
451     free(src);
452     
453     dest = (void *)malloc( 2*n * sizeof(double) );
454     src  = (void *)malloc( (1 + n) * ((1+stride)*sizeof(double)) );
455     memset( dest, 0, 2*n * sizeof(double) );
456     memset( src, 0, (1+n) * (1+stride)*sizeof(double) );
457     err = TestVecPack2Double( n, stride, &avgTimeUser, &avgTimeMPI,
458                               dest, src );
459     errs += Report( "VecPack2Double", "Pack", avgTimeMPI, avgTimeUser );
460
461     free(dest);
462     free(src);
463     
464
465
466     if (errs == 0) {
467         printf( " No Errors\n" );
468     }
469     else {
470         printf( " Found %d performance problems\n", errs );
471     }
472
473     fflush(stdout);
474     MPI_Finalize();
475
476     return 0;
477 }