Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'master' of https://framagit.org/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 50
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)
67         retval = -retval;
68
69     free(margin);
70     return retval;
71 }
72
73 /* Here are the tests */
74
75 /* Test packing a vector of individual doubles */
76 /* We don't use restrict in the function args because assignments between
77    restrict pointers is not valid in C and some compilers, such as the
78    IBM xlc compilers, flag that use as an error.*/
79 int TestVecPackDouble(int n, int stride,
80                       double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src);
81 int TestVecPackDouble(int n, int stride,
82                       double *avgTimeUser, double *avgTimeMPI, 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)
93         printf("TestVecPackDouble (USER): ");
94     for (j = 0; j < NTRIALS; j++) {
95         t1 = MPI_Wtime();
96         for (rep = 0; rep < N_REPS; rep++) {
97             i = n;
98             d_dest = dest;
99             d_src = src;
100             while (i--) {
101                 *d_dest++ = *d_src;
102                 d_src += stride;
103             }
104         }
105         t2 = MPI_Wtime() - t1;
106         t[j] = t2;
107         if (verbose)
108             printf("%.3f ", t[j]);
109     }
110     if (verbose)
111         printf("[%.3f]\n", noise(t, NTRIALS));
112     /* If there is too much noise, discard the test */
113     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
114         *avgTimeUser = 0;
115         *avgTimeMPI = 0;
116         if (verbose)
117             printf("Too much noise; discarding measurement\n");
118         return 0;
119     }
120     *avgTimeUser = mean(t, NTRIALS) / N_REPS;
121
122     /* MPI Vector code */
123     MPI_Type_vector(n, 1, stride, MPI_DOUBLE, &vectype);
124     MPI_Type_commit(&vectype);
125
126     if (verbose)
127         printf("TestVecPackDouble (MPI): ");
128     for (j = 0; j < NTRIALS; j++) {
129         t1 = MPI_Wtime();
130         for (rep = 0; rep < N_REPS; rep++) {
131             position = 0;
132             MPI_Pack((void *) src, 1, vectype, dest, n * sizeof(double), &position, MPI_COMM_SELF);
133         }
134         t2 = MPI_Wtime() - t1;
135         t[j] = t2;
136         if (verbose)
137             printf("%.3f ", t[j]);
138     }
139     if (verbose)
140         printf("[%.3f]\n", noise(t, NTRIALS));
141     /* If there is too much noise, discard the test */
142     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
143         *avgTimeUser = 0;
144         *avgTimeMPI = 0;
145         if (verbose)
146             printf("Too much noise; discarding measurement\n");
147     }
148     else {
149         *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
150     }
151
152     MPI_Type_free(&vectype);
153
154     return 0;
155 }
156
157 /* Test unpacking a vector of individual doubles */
158 /* See above for why restrict is not used in the function args */
159 int TestVecUnPackDouble(int n, int stride,
160                         double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src);
161 int TestVecUnPackDouble(int n, int stride,
162                         double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src)
163 {
164     double *restrict d_dest;
165     const double *restrict d_src;
166     register int i, j;
167     int rep, position;
168     double t1, t2, t[NTRIALS];
169     MPI_Datatype vectype;
170
171     /* User code */
172     if (verbose)
173         printf("TestVecUnPackDouble (USER): ");
174     for (j = 0; j < NTRIALS; j++) {
175         t1 = MPI_Wtime();
176         for (rep = 0; rep < N_REPS; rep++) {
177             i = n;
178             d_dest = dest;
179             d_src = src;
180             while (i--) {
181                 *d_dest = *d_src++;
182                 d_dest += stride;
183             }
184         }
185         t2 = MPI_Wtime() - t1;
186         t[j] = t2;
187         if (verbose)
188             printf("%.3f ", t[j]);
189     }
190     if (verbose)
191         printf("[%.3f]\n", noise(t, NTRIALS));
192     /* If there is too much noise, discard the test */
193     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
194         *avgTimeUser = 0;
195         *avgTimeMPI = 0;
196         if (verbose)
197             printf("Too much noise; discarding measurement\n");
198         return 0;
199     }
200     *avgTimeUser = mean(t, NTRIALS) / N_REPS;
201
202     /* MPI Vector code */
203     MPI_Type_vector(n, 1, stride, MPI_DOUBLE, &vectype);
204     MPI_Type_commit(&vectype);
205
206     if (verbose)
207         printf("TestVecUnPackDouble (MPI): ");
208     for (j = 0; j < NTRIALS; j++) {
209         t1 = MPI_Wtime();
210         for (rep = 0; rep < N_REPS; rep++) {
211             position = 0;
212             MPI_Unpack((void *) src, n * sizeof(double),
213                        &position, dest, 1, vectype, MPI_COMM_SELF);
214         }
215         t2 = MPI_Wtime() - t1;
216         t[j] = t2;
217         if (verbose)
218             printf("%.3f ", t[j]);
219     }
220     if (verbose)
221         printf("[%.3f]\n", noise(t, NTRIALS));
222     /* If there is too much noise, discard the test */
223     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
224         *avgTimeUser = 0;
225         *avgTimeMPI = 0;
226         if (verbose)
227             printf("Too much noise; discarding measurement\n");
228     }
229     else {
230         *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
231     }
232
233     MPI_Type_free(&vectype);
234
235     return 0;
236 }
237
238 /* Test packing a vector of 2-individual doubles */
239 /* See above for why restrict is not used in the function args */
240 int TestVecPack2Double(int n, int stride,
241                        double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src);
242 int TestVecPack2Double(int n, int stride,
243                        double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src)
244 {
245     double *restrict d_dest;
246     const double *restrict d_src;
247     register int i, j;
248     int rep, position;
249     double t1, t2, t[NTRIALS];
250     MPI_Datatype vectype;
251
252     /* User code */
253     if (verbose)
254         printf("TestVecPack2Double (USER): ");
255     for (j = 0; j < NTRIALS; j++) {
256         t1 = MPI_Wtime();
257         for (rep = 0; rep < N_REPS; rep++) {
258             i = n;
259             d_dest = dest;
260             d_src = src;
261             while (i--) {
262                 *d_dest++ = d_src[0];
263                 *d_dest++ = d_src[1];
264                 d_src += stride;
265             }
266         }
267         t2 = MPI_Wtime() - t1;
268         t[j] = t2;
269         if (verbose)
270             printf("%.3f ", t[j]);
271     }
272     if (verbose)
273         printf("[%.3f]\n", noise(t, NTRIALS));
274     /* If there is too much noise, discard the test */
275     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
276         *avgTimeUser = 0;
277         *avgTimeMPI = 0;
278         if (verbose)
279             printf("Too much noise; discarding measurement\n");
280         return 0;
281     }
282     *avgTimeUser = mean(t, NTRIALS) / N_REPS;
283
284     /* MPI Vector code */
285     MPI_Type_vector(n, 2, stride, MPI_DOUBLE, &vectype);
286     MPI_Type_commit(&vectype);
287
288     if (verbose)
289         printf("TestVecPack2Double (MPI): ");
290     for (j = 0; j < NTRIALS; j++) {
291         t1 = MPI_Wtime();
292         for (rep = 0; rep < N_REPS; rep++) {
293             position = 0;
294             MPI_Pack((void *) src, 1, vectype, dest, 2 * n * sizeof(double),
295                      &position, MPI_COMM_SELF);
296         }
297         t2 = MPI_Wtime() - t1;
298         t[j] = t2;
299         if (verbose)
300             printf("%.3f ", t[j]);
301     }
302     if (verbose)
303         printf("[%.3f]\n", noise(t, NTRIALS));
304     /* If there is too much noise, discard the test */
305     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
306         *avgTimeUser = 0;
307         *avgTimeMPI = 0;
308         if (verbose)
309             printf("Too much noise; discarding measurement\n");
310     }
311     else {
312         *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
313     }
314     MPI_Type_free(&vectype);
315
316     return 0;
317 }
318
319 /* This creates an indexed type that is like a vector (for simplicity
320    of construction).  There is a possibility that the MPI implementation
321    will recognize and simplify this (e.g., in MPI_Type_commit); if so,
322    let us know and we'll add a version that is not as regular
323 */
324 /* See above for why restrict is not used in the function args */
325 int TestIndexPackDouble(int n, int stride,
326                         double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src);
327 int TestIndexPackDouble(int n, int stride,
328                         double *avgTimeUser, double *avgTimeMPI, double *dest, const double *src)
329 {
330     double *restrict d_dest;
331     const double *restrict d_src;
332     register int i, j;
333     int rep, position;
334     int *restrict displs = 0;
335     double t1, t2, t[NTRIALS];
336     MPI_Datatype indextype;
337
338     displs = (int *) malloc(n * sizeof(int));
339     for (i = 0; i < n; i++)
340         displs[i] = i * stride;
341
342     /* User code */
343     if (verbose)
344         printf("TestIndexPackDouble (USER): ");
345     for (j = 0; j < NTRIALS; j++) {
346         t1 = MPI_Wtime();
347         for (rep = 0; rep < N_REPS; rep++) {
348             i = n;
349             d_dest = dest;
350             d_src = src;
351             for (i = 0; i < n; i++) {
352                 *d_dest++ = d_src[displs[i]];
353             }
354         }
355         t2 = MPI_Wtime() - t1;
356         t[j] = t2;
357         if (verbose)
358             printf("%.3f ", t[j]);
359     }
360     if (verbose)
361         printf("[%.3f]\n", noise(t, NTRIALS));
362     /* If there is too much noise, discard the test */
363     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
364         *avgTimeUser = 0;
365         *avgTimeMPI = 0;
366         if (verbose)
367             printf("Too much noise; discarding measurement\n");
368         return 0;
369     }
370     *avgTimeUser = mean(t, NTRIALS) / N_REPS;
371
372     /* MPI Index code */
373     MPI_Type_create_indexed_block(n, 1, displs, MPI_DOUBLE, &indextype);
374     MPI_Type_commit(&indextype);
375
376     free(displs);
377
378     if (verbose)
379         printf("TestIndexPackDouble (MPI): ");
380     for (j = 0; j < NTRIALS; j++) {
381         t1 = MPI_Wtime();
382         for (rep = 0; rep < N_REPS; rep++) {
383             position = 0;
384             MPI_Pack((void *) src, 1, indextype, dest, n * sizeof(double),
385                      &position, MPI_COMM_SELF);
386         }
387         t2 = MPI_Wtime() - t1;
388         t[j] = t2;
389         if (verbose)
390             printf("%.3f ", t[j]);
391     }
392     if (verbose)
393         printf("[%.3f]\n", noise(t, NTRIALS));
394     /* If there is too much noise, discard the test */
395     if (noise(t, NTRIALS) > VARIANCE_THRESHOLD) {
396         *avgTimeUser = 0;
397         *avgTimeMPI = 0;
398         if (verbose)
399             printf("Too much noise; discarding measurement\n");
400     }
401     else {
402         *avgTimeMPI = mean(t, NTRIALS) / N_REPS;
403     }
404     MPI_Type_free(&indextype);
405
406     return 0;
407 }
408
409 int Report(const char *name, const char *packname, double avgTimeMPI, double avgTimeUser);
410 int Report(const char *name, const char *packname, double avgTimeMPI, double avgTimeUser)
411 {
412     double diffTime, maxTime;
413     int errs = 0;
414
415     /* Move this into a common routine */
416     diffTime = avgTimeMPI - avgTimeUser;
417     if (diffTime < 0)
418         diffTime = -diffTime;
419     if (avgTimeMPI > avgTimeUser)
420         maxTime = avgTimeMPI;
421     else
422         maxTime = avgTimeUser;
423
424     if (verbose) {
425         printf("%-30s:\t%g\t%g\t(%g%%)\n", name,
426                avgTimeMPI, avgTimeUser, 100 * (diffTime / maxTime));
427         fflush(stdout);
428     }
429     if (avgTimeMPI > avgTimeUser && (diffTime > THRESHOLD * maxTime)) {
430         errs++;
431         printf("%s:\tMPI %s code is too slow: MPI %g\t User %g\n",
432                name, packname, avgTimeMPI, avgTimeUser);
433     }
434
435     return errs;
436 }
437
438 /* Finally, here's the main program */
439 int main(int argc, char *argv[])
440 {
441     int n, stride, errs = 0;
442     void *dest, *src;
443     double avgTimeUser, avgTimeMPI;
444
445     MPI_Init(&argc, &argv);
446     if (getenv("MPITEST_VERBOSE"))
447         verbose = 1;
448
449     n = 30000;
450     stride = 4;
451     dest = (void *) malloc(n * sizeof(double));
452     src = (void *) malloc(n * ((1 + stride) * sizeof(double)));
453     /* Touch the source and destination arrays */
454     memset(src, 0, n * (1 + stride) * sizeof(double));
455     memset(dest, 0, n * sizeof(double));
456
457     TestVecPackDouble(n, stride, &avgTimeUser, &avgTimeMPI, dest, src);
458     errs += Report("VecPackDouble", "Pack", avgTimeMPI, avgTimeUser);
459
460     TestVecUnPackDouble(n, stride, &avgTimeUser, &avgTimeMPI, src, dest);
461     errs += Report("VecUnPackDouble", "Unpack", avgTimeMPI, avgTimeUser);
462
463     TestIndexPackDouble(n, stride, &avgTimeUser, &avgTimeMPI, dest, src);
464     errs += Report("VecIndexDouble", "Pack", avgTimeMPI, avgTimeUser);
465
466     free(dest);
467     free(src);
468
469     dest = (void *) malloc(2 * n * sizeof(double));
470     src = (void *) malloc((1 + n) * ((1 + stride) * sizeof(double)));
471     memset(dest, 0, 2 * n * sizeof(double));
472     memset(src, 0, (1 + n) * (1 + stride) * sizeof(double));
473     TestVecPack2Double(n, stride, &avgTimeUser, &avgTimeMPI, dest, src);
474     errs += Report("VecPack2Double", "Pack", avgTimeMPI, avgTimeUser);
475
476     free(dest);
477     free(src);
478
479
480
481     if (errs == 0) {
482         printf(" No Errors\n");
483     }
484     else {
485         printf(" Found %d performance problems\n", errs);
486     }
487
488     fflush(stdout);
489     MPI_Finalize();
490
491     return 0;
492 }