This source file includes following definitions.
- comp_double
- get_tick
- timing_delay
- op_send
- op_send_pingpong
- op_coll
- op_a2a
- op_put
- op_get
- do_bench
- main
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23 #include <stdlib.h>
24 #include <stdio.h>
25 #include <sys/time.h>
26 #include <time.h>
27 #include <string.h>
28 #include "mpi.h"
29
30 #define NB_ITER 1000
31 #define FULL_NB_ITER (size_world * NB_ITER)
32 #define MAX_SIZE (1024 * 1024 * 1.4)
33 #define NB_OPS 6
34
35 static int rank_world = -1;
36 static int size_world = 0;
37 static int to = -1;
38 static int from = -1;
39 static MPI_Win win = MPI_WIN_NULL;
40
41
42 static int comp_double(const void*_a, const void*_b)
43 {
44 const double*a = _a;
45 const double*b = _b;
46 if(*a < *b)
47 return -1;
48 else if(*a > *b)
49 return 1;
50 else
51 return 0;
52 }
53
54
55 static inline void get_tick(struct timespec*t)
56 {
57 #if defined(__bg__)
58 # define CLOCK_TYPE CLOCK_REALTIME
59 #elif defined(CLOCK_MONOTONIC_RAW)
60 # define CLOCK_TYPE CLOCK_MONOTONIC_RAW
61 #elif defined(CLOCK_MONOTONIC)
62 # define CLOCK_TYPE CLOCK_MONOTONIC
63 #endif
64 #if defined(CLOCK_TYPE)
65 clock_gettime(CLOCK_TYPE, t);
66 #else
67 struct timeval tv;
68 gettimeofday(&tv, NULL);
69 t->tv_sec = tv.tv_sec;
70 t->tv_nsec = tv.tv_usec * 1000;
71 #endif
72 }
73 static inline double timing_delay(const struct timespec*const t1, const struct timespec*const t2)
74 {
75 const double delay = 1000000.0 * (t2->tv_sec - t1->tv_sec) + (t2->tv_nsec - t1->tv_nsec) / 1000.0;
76 return delay;
77 }
78
79
80 static inline void op_send(double*res, void*sbuf, int size, int tagno, void*rbuf) {
81 MPI_Request request;
82 struct timespec start, end;
83
84
85 MPI_Irecv(rbuf, size, MPI_BYTE, from, tagno, MPI_COMM_WORLD, &request);
86
87
88
89
90 if( 0 == rank_world ) {
91 MPI_Send(NULL, 0, MPI_BYTE, from, 100, MPI_COMM_WORLD);
92 MPI_Recv(NULL, 0, MPI_BYTE, to, 100, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
93 } else {
94 MPI_Recv(NULL, 0, MPI_BYTE, to, 100, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
95 MPI_Send(NULL, 0, MPI_BYTE, from, 100, MPI_COMM_WORLD);
96 }
97
98
99 get_tick(&start);
100 MPI_Send(sbuf, size, MPI_BYTE, to, tagno, MPI_COMM_WORLD);
101 get_tick(&end);
102
103 MPI_Wait(&request, MPI_STATUS_IGNORE);
104 *res = timing_delay(&start, &end);
105 }
106
107 static inline void op_send_pingpong(double*res, void*sbuf, int size, int tagno, void*rbuf) {
108 struct timespec start, end;
109
110 MPI_Barrier(MPI_COMM_WORLD);
111
112
113 if(rank_world % 2) {
114 MPI_Recv(rbuf, size, MPI_BYTE, from, tagno, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
115 MPI_Send(sbuf, size, MPI_BYTE, from, tagno, MPI_COMM_WORLD);
116 MPI_Barrier(MPI_COMM_WORLD);
117 get_tick(&start);
118 MPI_Send(sbuf, size, MPI_BYTE, from, tagno, MPI_COMM_WORLD);
119 MPI_Recv(rbuf, size, MPI_BYTE, from, tagno, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
120 get_tick(&end);
121 } else {
122 get_tick(&start);
123 MPI_Send(sbuf, size, MPI_BYTE, to, tagno, MPI_COMM_WORLD);
124 MPI_Recv(rbuf, size, MPI_BYTE, to, tagno, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
125 get_tick(&end);
126 MPI_Barrier(MPI_COMM_WORLD);
127 MPI_Recv(rbuf, size, MPI_BYTE, to, tagno, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
128 MPI_Send(sbuf, size, MPI_BYTE, to, tagno, MPI_COMM_WORLD);
129 }
130
131 *res = timing_delay(&start, &end) / 2;
132 }
133
134 static inline void op_coll(double*res, void*buff, int size, int tagno, void*rbuf) {
135 struct timespec start, end;
136 MPI_Barrier(MPI_COMM_WORLD);
137
138
139 get_tick(&start);
140 MPI_Bcast(buff, size, MPI_BYTE, 0, MPI_COMM_WORLD);
141 get_tick(&end);
142
143 *res = timing_delay(&start, &end);
144 }
145
146 static inline void op_a2a(double*res, void*sbuf, int size, int tagno, void*rbuf) {
147 struct timespec start, end;
148 MPI_Barrier(MPI_COMM_WORLD);
149
150
151 get_tick(&start);
152 MPI_Alltoall(sbuf, size, MPI_BYTE, rbuf, size, MPI_BYTE, MPI_COMM_WORLD);
153 get_tick(&end);
154
155 *res = timing_delay(&start, &end);
156 }
157
158 static inline void op_put(double*res, void*sbuf, int size, int tagno, void*rbuf) {
159 struct timespec start, end;
160
161 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, to, 0, win);
162
163
164 get_tick(&start);
165 MPI_Put(sbuf, size, MPI_BYTE, to, 0, size, MPI_BYTE, win);
166 MPI_Win_unlock(to, win);
167 get_tick(&end);
168
169 *res = timing_delay(&start, &end);
170 }
171
172 static inline void op_get(double*res, void*rbuf, int size, int tagno, void*sbuf) {
173 struct timespec start, end;
174
175 MPI_Win_lock(MPI_LOCK_SHARED, to, 0, win);
176
177
178 get_tick(&start);
179 MPI_Get(rbuf, size, MPI_BYTE, to, 0, size, MPI_BYTE, win);
180 MPI_Win_unlock(to, win);
181 get_tick(&end);
182
183 *res = timing_delay(&start, &end);
184 }
185
186 static inline void do_bench(int size, char*sbuf, double*results,
187 void(*op)(double*, void*, int, int, void*)) {
188 int iter;
189 int tagno = 201;
190 char*rbuf = sbuf ? sbuf + size : NULL;
191
192 if(op == op_put || op == op_get){
193 win = MPI_WIN_NULL;
194 MPI_Win_create(rbuf, size, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &win);
195 }
196
197 for( iter = 0; iter < NB_ITER; ++iter ) {
198 op(&results[iter], sbuf, size, tagno, rbuf);
199 MPI_Barrier(MPI_COMM_WORLD);
200 }
201
202 if(op == op_put || op == op_get){
203 MPI_Win_free(&win);
204 win = MPI_WIN_NULL;
205 }
206 }
207
208 int main(int argc, char* argv[])
209 {
210 int size, iter, nop;
211 char*sbuf = NULL;
212 double results[NB_ITER];
213 void(*op)(double*, void*, int, int, void*);
214 char name[255];
215 MPI_Init(&argc, &argv);
216 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
217 MPI_Comm_size(MPI_COMM_WORLD, &size_world);
218 to = (rank_world + 1) % size_world;
219 from = (rank_world + size_world - 1) % size_world;
220
221 double full_res[FULL_NB_ITER];
222
223 for( nop = 0; nop < NB_OPS; ++nop ) {
224 switch(nop) {
225 case 0:
226 op = op_send;
227 sprintf(name, "MPI_Send");
228 break;
229 case 1:
230 op = op_coll;
231 sprintf(name, "MPI_Bcast");
232 break;
233 case 2:
234 op = op_a2a;
235 sprintf(name, "MPI_Alltoall");
236 break;
237 case 3:
238 op = op_send_pingpong;
239 sprintf(name, "MPI_Send_pp");
240 break;
241 case 4:
242 op = op_put;
243 sprintf(name, "MPI_Put");
244 break;
245 case 5:
246 op = op_get;
247 sprintf(name, "MPI_Get");
248 break;
249 }
250
251 if( 0 == rank_world )
252 printf("# %s%%%d\n# size \t| latency \t| 10^6 B/s \t| MB/s \t| median \t| q1 \t| q3 \t| d1 \t| d9 \t| avg \t| max\n", name, size_world);
253
254 for(size = 0; size < MAX_SIZE; size = ((int)(size * 1.4) > size) ? (size * 1.4) : (size + 1)) {
255
256 if( 0 != size ) {
257 sbuf = (char *)realloc(sbuf, (size_world + 1) * size);
258 }
259
260 do_bench(size, sbuf, results, op);
261
262 MPI_Gather(results, NB_ITER, MPI_DOUBLE, full_res, NB_ITER, MPI_DOUBLE, 0, MPI_COMM_WORLD);
263
264 if( 0 == rank_world ) {
265 qsort(full_res, FULL_NB_ITER, sizeof(double), &comp_double);
266 const double min_lat = full_res[0];
267 const double max_lat = full_res[FULL_NB_ITER - 1];
268 const double med_lat = full_res[(FULL_NB_ITER - 1) / 2];
269 const double q1_lat = full_res[(FULL_NB_ITER - 1) / 4];
270 const double q3_lat = full_res[ 3 * (FULL_NB_ITER - 1) / 4];
271 const double d1_lat = full_res[(FULL_NB_ITER - 1) / 10];
272 const double d9_lat = full_res[ 9 * (FULL_NB_ITER - 1) / 10];
273 double avg_lat = 0.0;
274 for( iter = 0; iter < FULL_NB_ITER; iter++ ){
275 avg_lat += full_res[iter];
276 }
277 avg_lat /= FULL_NB_ITER;
278 const double bw_million_byte = size / min_lat;
279 const double bw_mbyte = bw_million_byte / 1.048576;
280
281 printf("%9lld\t%9.3lf\t%9.3f\t%9.3f\t%9.3lf\t%9.3lf\t%9.3lf\t%9.3lf\t%9.3lf\t%9.3lf\t%9.3lf",
282 (long long)size, min_lat, bw_million_byte, bw_mbyte,
283 med_lat, q1_lat, q3_lat, d1_lat, d9_lat,
284 avg_lat, max_lat);
285 printf("\n");
286 }
287 }
288 free(sbuf);
289 sbuf = NULL;
290 }
291
292 MPI_Finalize();
293 return EXIT_SUCCESS;
294 }