Logo AND Algorithmique Numérique Distribuée

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