Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
dtd update for new cluster routing options
[simgrid.git] / teshsuite / smpi / mpich3-test / coll / allgatherv4.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *
4  *  (C) 2003 by Argonne National Laboratory.
5  *      See COPYRIGHT in top-level directory.
6  */
7
8 #include "mpi.h"
9 #include "mpitest.h"
10 #include <stdio.h>
11 #include <stdlib.h>
12 #ifdef HAVE_SYS_TIME_H
13 #include <sys/time.h>
14 #endif
15 #include <time.h>
16 #include <math.h>
17 #include <assert.h>
18
19 /* FIXME: What is this test supposed to accomplish? */
20
21 #define START_BUF (1)
22 #define LARGE_BUF (256 * 1024)
23
24 /* FIXME: MAX_BUF is too large */
25 #define MAX_BUF   (32 * 1024 * 1024)
26 #define LOOPS 10
27
28 __thread char * sbuf, * rbuf;
29 __thread int * recvcounts, * displs;
30 int errs = 0;
31
32 /* #define dprintf printf */
33 #define dprintf(...)
34
35 typedef enum {
36     REGULAR,
37     BCAST,
38     SPIKE,
39     HALF_FULL,
40     LINEAR_DECREASE,
41     BELL_CURVE
42 } test_t;
43
44 void comm_tests(MPI_Comm comm);
45 double run_test(long long msg_size, MPI_Comm comm, test_t test_type, double * max_time);
46
47 int main(int argc, char ** argv)
48 {
49     int comm_size, comm_rank;
50     MPI_Comm comm;
51
52     MTest_Init(&argc, &argv);
53     MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
54     MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
55
56     if (LARGE_BUF * comm_size > MAX_BUF)
57         goto fn_exit;
58
59     sbuf = (void *) calloc(MAX_BUF, 1);
60     rbuf = (void *) calloc(MAX_BUF, 1);
61
62     srand(time(NULL));
63
64     recvcounts = (void *) malloc(comm_size * sizeof(int));
65     displs = (void *) malloc(comm_size * sizeof(int));
66     if (!recvcounts || !displs || !sbuf || !rbuf) {
67         fprintf(stderr, "Unable to allocate memory:\n");
68         if (!sbuf) fprintf(stderr,"\tsbuf of %d bytes\n", MAX_BUF );
69         if (!rbuf) fprintf(stderr,"\trbuf of %d bytes\n", MAX_BUF );
70         if (!recvcounts) fprintf(stderr,"\trecvcounts of %zd bytes\n", comm_size * sizeof(int) );
71         if (!displs) fprintf(stderr,"\tdispls of %zd bytes\n", comm_size * sizeof(int) );
72         fflush(stderr);
73         MPI_Abort(MPI_COMM_WORLD, -1);
74         exit(-1);
75     }
76
77     if (!comm_rank) {
78         dprintf("Message Range: (%d, %d); System size: %d\n", START_BUF, LARGE_BUF, comm_size);
79         fflush(stdout);
80     }
81
82
83     /* COMM_WORLD tests */
84     if (!comm_rank) {
85         dprintf("\n\n==========================================================\n");
86         dprintf("                         MPI_COMM_WORLD\n");
87         dprintf("==========================================================\n");
88     }
89     comm_tests(MPI_COMM_WORLD);
90
91     /* non-COMM_WORLD tests */
92     if (!comm_rank) {
93         dprintf("\n\n==========================================================\n");
94         dprintf("                         non-COMM_WORLD\n");
95         dprintf("==========================================================\n");
96     }
97     MPI_Comm_split(MPI_COMM_WORLD, (comm_rank == comm_size - 1) ? 0 : 1, 0, &comm);
98     if (comm_rank < comm_size - 1)
99         comm_tests(comm);
100     MPI_Comm_free(&comm);
101
102     /* Randomized communicator tests */
103     if (!comm_rank) {
104         dprintf("\n\n==========================================================\n");
105         dprintf("                         Randomized Communicator\n");
106         dprintf("==========================================================\n");
107     }
108     MPI_Comm_split(MPI_COMM_WORLD, 0, rand(), &comm);
109     comm_tests(comm);
110     MPI_Comm_free(&comm);
111
112     //free(sbuf);
113     //free(rbuf);
114     free(recvcounts);
115     free(displs);
116
117 fn_exit:
118     MTest_Finalize(errs);
119     MPI_Finalize();
120
121     return 0;
122 }
123
124 void comm_tests(MPI_Comm comm)
125 {
126     int comm_size, comm_rank;
127     double rtime, max_time;
128     long long msg_size;
129
130     MPI_Comm_size(comm, &comm_size);
131     MPI_Comm_rank(comm, &comm_rank);
132
133     for (msg_size = START_BUF; msg_size <= LARGE_BUF; msg_size *= 2) {
134         if (!comm_rank) {
135             dprintf("\n====> MSG_SIZE: %d\n", (int) msg_size);
136             fflush(stdout);
137         }
138
139         rtime = run_test(msg_size, comm, REGULAR, &max_time);
140         if (!comm_rank) {
141             dprintf("REGULAR:\tAVG: %.3f\tMAX: %.3f\n", rtime, max_time);
142             fflush(stdout);
143         }
144
145         rtime = run_test(msg_size, comm, BCAST, &max_time);
146         if (!comm_rank) {
147             dprintf("BCAST:\tAVG: %.3f\tMAX: %.3f\n", rtime, max_time);
148             fflush(stdout);
149         }
150
151         rtime = run_test(msg_size, comm, SPIKE, &max_time);
152         if (!comm_rank) {
153             dprintf("SPIKE:\tAVG: %.3f\tMAX: %.3f\n", rtime, max_time);
154             fflush(stdout);
155         }
156
157         rtime = run_test(msg_size, comm, HALF_FULL, &max_time);
158         if (!comm_rank) {
159             dprintf("HALF_FULL:\tAVG: %.3f\tMAX: %.3f\n", rtime, max_time);
160             fflush(stdout);
161         }
162
163         rtime = run_test(msg_size, comm, LINEAR_DECREASE, &max_time);
164         if (!comm_rank) {
165             dprintf("LINEAR_DECREASE:\tAVG: %.3f\tMAX: %.3f\n", rtime, max_time);
166             fflush(stdout);
167         }
168
169         rtime = run_test(msg_size, comm, BELL_CURVE, &max_time);
170         if (!comm_rank) {
171             dprintf("BELL_CURVE:\tAVG: %.3f\tMAX: %.3f\n", rtime, max_time);
172             fflush(stdout);
173         }
174     }
175 }
176
177 double run_test(long long msg_size, MPI_Comm comm, test_t test_type, 
178                 double * max_time)
179 {
180     int i, j;
181     int comm_size, comm_rank;
182     double start, end;
183     double total_time, avg_time;
184     MPI_Aint tmp;
185
186     MPI_Comm_size(comm, &comm_size);
187     MPI_Comm_rank(comm, &comm_rank);
188
189     displs[0] = 0;
190     for (i = 0; i < comm_size; i++) {
191         if (test_type == REGULAR)
192             recvcounts[i] = msg_size;
193         else if (test_type == BCAST)
194             recvcounts[i] = (!i) ? msg_size : 0;
195         else if (test_type == SPIKE)
196             recvcounts[i] = (!i) ? (msg_size / 2) : (msg_size / (2 * (comm_size - 1)));
197         else if (test_type == HALF_FULL)
198             recvcounts[i] = (i < (comm_size / 2)) ? (2 * msg_size) : 0;
199         else if (test_type == LINEAR_DECREASE) {
200             tmp = 2 * msg_size * (comm_size - 1 - i) / (comm_size - 1);
201             if (tmp != (int)tmp) {
202                 fprintf( stderr, "Integer overflow in variable tmp\n" );
203                 MPI_Abort( MPI_COMM_WORLD, 1 );
204                 exit(1);
205             }
206             recvcounts[i] = (int) tmp;
207
208             /* If the maximum message size is too large, don't run */
209             if (tmp > MAX_BUF) return 0;
210         }
211         else if (test_type == BELL_CURVE) {
212             for (j = 0; j < i; j++) {
213                 if (i - 1 + j >= comm_size) continue;
214                 tmp = msg_size * comm_size / (log(comm_size) * i);
215                 recvcounts[i - 1 + j] = (int) tmp;
216                 displs[i - 1 + j] = 0;
217
218                 /* If the maximum message size is too large, don't run */
219                 if (tmp > MAX_BUF) return 0;
220             }
221         }
222
223         if (i < comm_size - 1)
224             displs[i+1] = displs[i] + recvcounts[i];
225     }
226
227     /* Test that:
228        1: sbuf is large enough
229        2: rbuf is large enough
230        3: There were no failures (e.g., tmp nowhere > rbuf size 
231     */
232     MPI_Barrier(comm);
233     start = MPI_Wtime();
234     for (i = 0; i < LOOPS; i++) {
235         MPI_Allgatherv(sbuf, recvcounts[comm_rank], MPI_CHAR,
236                        rbuf, recvcounts, displs, MPI_CHAR, comm);
237     }
238     end = MPI_Wtime();
239     MPI_Barrier(comm);
240
241     /* Convert to microseconds (why?) */
242     total_time = 1.0e6 * (end - start);
243     MPI_Reduce(&total_time, &avg_time, 1, MPI_DOUBLE, MPI_SUM, 0, comm);
244     MPI_Reduce(&total_time, max_time, 1, MPI_DOUBLE, MPI_MAX, 0, comm);
245
246     return (avg_time / (LOOPS * comm_size));
247 }