Logo AND Algorithmique Numérique Distribuée

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