Logo AND Algorithmique Numérique Distribuée

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