Logo AND Algorithmique Numérique Distribuée

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