Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Prepare the out of source build.
[simgrid.git] / examples / smpi / smpi_traced.c
1 /* Copyright (c) 2009, 2010. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "mpi.h"
8 #include <stdio.h>
9 #include "instr/instr.h"
10
11 #define DATATOSENT 100000
12
13 int main(int argc, char *argv[])
14 {
15   int rank, numprocs, i;
16
17   MPI_Init(&argc, &argv);
18   MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
19   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
20
21   //Tests: 
22   //A: 0(isend -> wait) with 1(recv)
23   int A = 1;
24   //B: 0(send) with 1(irecv -> wait)
25   int B = 1;
26   //C: 0(N * isend -> N * wait) with 1(N * recv)
27   int C = 1;
28   //D: 0(N * isend -> N * waitany) with 1(N * recv)
29   int D = 1;
30   //E: 0(N*send) with 1(N*irecv, N*wait)
31   int E = 1;
32   //F: 0(N*send) with 1(N*irecv, N*waitany)
33   int F = 1;
34   //G: 0(N* isend -> waitall) with 1(N*recv)
35   int G = 1;
36   //H: 0(N*send) with 1(N*irecv, waitall)
37   int H = 1;
38   //I: 0(2*N*send, 2*N*Irecv, Waitall) with
39   //   1(N*irecv, waitall, N*isend, N*waitany) with 
40   //   2(N*irecv, N*waitany, N*isend, waitall)
41   int I = 1;
42   //J: 0(N*isend, N*test, N*wait) with (N*irecv, N*test, N*wait)
43   int J = 1;
44
45   int N = 13;
46
47   int tag = 12345;
48 /////////////////////////////////////////
49 ////////////////// RANK 0
50 ///////////////////////////////////
51   if (rank == 0) {
52     MPI_Request request;
53     MPI_Status status;
54     MPI_Request req[2 * N];
55     MPI_Status sta[2 * N];
56     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
57     if (A) {
58       TRACE_smpi_set_category("A");
59       MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &request);
60       MPI_Wait(&request, &status);
61     }
62     MPI_Barrier(MPI_COMM_WORLD);
63
64     if (B) {
65       TRACE_smpi_set_category("B");
66       MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
67     }
68     MPI_Barrier(MPI_COMM_WORLD);
69
70     if (C) {
71       TRACE_smpi_set_category("C");
72       for (i = 0; i < N; i++) {
73         MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
74       }
75       for (i = 0; i < N; i++) {
76         MPI_Wait(&req[i], &sta[i]);
77       }
78     }
79     MPI_Barrier(MPI_COMM_WORLD);
80
81     if (D) {
82       TRACE_smpi_set_category("D");
83       for (i = 0; i < N; i++) {
84         MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
85       }
86       for (i = 0; i < N; i++) {
87         int completed;
88         MPI_Waitany(N, req, &completed, sta);
89       }
90     }
91     MPI_Barrier(MPI_COMM_WORLD);
92
93     if (E) {
94       TRACE_smpi_set_category("E");
95       for (i = 0; i < N; i++) {
96         MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
97       }
98     }
99     MPI_Barrier(MPI_COMM_WORLD);
100
101     if (F) {
102       TRACE_smpi_set_category("F");
103       for (i = 0; i < N; i++) {
104         MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
105       }
106     }
107     MPI_Barrier(MPI_COMM_WORLD);
108
109     if (G) {
110       TRACE_smpi_set_category("G");
111       for (i = 0; i < N; i++) {
112         MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
113       }
114       MPI_Waitall(N, req, sta);
115     }
116     MPI_Barrier(MPI_COMM_WORLD);
117
118     if (H) {
119       TRACE_smpi_set_category("H");
120       for (i = 0; i < N; i++) {
121         MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
122       }
123     }
124     MPI_Barrier(MPI_COMM_WORLD);
125
126     if (I) {
127       TRACE_smpi_set_category("I");
128       for (i = 0; i < 2 * N; i++) {
129         if (i < N) {
130           MPI_Send(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD);
131         } else {
132           MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
133         }
134       }
135       MPI_Barrier(MPI_COMM_WORLD);
136       for (i = 0; i < 2 * N; i++) {
137         if (i < N) {
138           MPI_Irecv(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD,
139                     &req[i]);
140         } else {
141           MPI_Irecv(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD,
142                     &req[i]);
143         }
144       }
145       MPI_Waitall(2 * N, req, sta);
146     }
147     MPI_Barrier(MPI_COMM_WORLD);
148
149     if (J) {
150       TRACE_smpi_set_category("J");
151       for (i = 0; i < N; i++) {
152         MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
153       }
154       for (i = 0; i < N; i++) {
155         int flag;
156         MPI_Test(&req[i], &flag, &sta[i]);
157       }
158       for (i = 0; i < N; i++) {
159         MPI_Wait(&req[i], &sta[i]);
160       }
161     }
162 /////////////////////////////////////////
163 ////////////////// RANK 1
164 ///////////////////////////////////
165   } else if (rank == 1) {
166     MPI_Request request;
167     MPI_Status status;
168     MPI_Request req[N];
169     MPI_Status sta[N];
170     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
171
172     if (A) {
173       TRACE_smpi_set_category("A");
174       MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
175     }
176     MPI_Barrier(MPI_COMM_WORLD);
177
178     if (B) {
179       TRACE_smpi_set_category("B");
180       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &request);
181       MPI_Wait(&request, &status);
182     }
183     MPI_Barrier(MPI_COMM_WORLD);
184
185     if (C) {
186       TRACE_smpi_set_category("C");
187       for (i = 0; i < N; i++) {
188         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
189       }
190     }
191     MPI_Barrier(MPI_COMM_WORLD);
192
193     if (D) {
194       TRACE_smpi_set_category("D");
195       for (i = 0; i < N; i++) {
196         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
197       }
198     }
199     MPI_Barrier(MPI_COMM_WORLD);
200
201     if (E) {
202       TRACE_smpi_set_category("E");
203       for (i = 0; i < N; i++) {
204         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
205       }
206       for (i = 0; i < N; i++) {
207         MPI_Wait(&req[i], &sta[i]);
208       }
209     }
210     MPI_Barrier(MPI_COMM_WORLD);
211
212     if (F) {
213       TRACE_smpi_set_category("F");
214       for (i = 0; i < N; i++) {
215         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
216       }
217       for (i = 0; i < N; i++) {
218         int completed;
219         MPI_Waitany(N, req, &completed, sta);
220       }
221     }
222     MPI_Barrier(MPI_COMM_WORLD);
223
224     if (G) {
225       TRACE_smpi_set_category("G");
226       for (i = 0; i < N; i++) {
227         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
228       }
229     }
230     MPI_Barrier(MPI_COMM_WORLD);
231
232     if (H) {
233       TRACE_smpi_set_category("H");
234       for (i = 0; i < N; i++) {
235         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
236       }
237       MPI_Waitall(N, req, sta);
238     }
239     MPI_Barrier(MPI_COMM_WORLD);
240
241     if (I) {
242       TRACE_smpi_set_category("I");
243       for (i = 0; i < N; i++) {
244         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
245       }
246       MPI_Waitall(N, req, sta);
247
248       MPI_Barrier(MPI_COMM_WORLD);
249       for (i = 0; i < N; i++) {
250         MPI_Isend(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
251       }
252       MPI_Waitall(N, req, sta);
253 //      for (i = 0; i < N; i++){
254 //        MPI_Wait (&req[i], &sta[i]);
255 //      }
256     }
257     MPI_Barrier(MPI_COMM_WORLD);
258
259     if (J) {
260       TRACE_smpi_set_category("J");
261       for (i = 0; i < N; i++) {
262         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
263       }
264       for (i = 0; i < N; i++) {
265         int flag;
266         MPI_Test(&req[i], &flag, &sta[i]);
267       }
268       for (i = 0; i < N; i++) {
269         MPI_Wait(&req[i], &sta[i]);
270       }
271     }
272 /////////////////////////////////////////
273 ////////////////// RANK 2
274 ///////////////////////////////////
275   } else if (rank == 2) {
276 //    MPI_Request request;
277 //    MPI_Status status;
278     MPI_Request req[N];
279     MPI_Status sta[N];
280     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
281
282     if (A) {
283     }
284     MPI_Barrier(MPI_COMM_WORLD);
285     if (B) {
286     }
287     MPI_Barrier(MPI_COMM_WORLD);
288     if (C) {
289     }
290     MPI_Barrier(MPI_COMM_WORLD);
291     if (D) {
292     }
293     MPI_Barrier(MPI_COMM_WORLD);
294     if (E) {
295     }
296     MPI_Barrier(MPI_COMM_WORLD);
297     if (F) {
298     }
299     MPI_Barrier(MPI_COMM_WORLD);
300     if (G) {
301     }
302     MPI_Barrier(MPI_COMM_WORLD);
303     if (H) {
304     }
305     MPI_Barrier(MPI_COMM_WORLD);
306     if (I) {
307       TRACE_smpi_set_category("I");
308       for (i = 0; i < N; i++) {
309         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
310       }
311       for (i = 0; i < N; i++) {
312         int completed;
313         MPI_Waitany(N, req, &completed, sta);
314       }
315       MPI_Barrier(MPI_COMM_WORLD);
316
317       for (i = 0; i < N; i++) {
318         MPI_Send(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD);
319       }
320     }
321     MPI_Barrier(MPI_COMM_WORLD);
322     if (J) {
323     }
324   }
325   MPI_Barrier(MPI_COMM_WORLD);
326   MPI_Finalize();
327   return 0;
328 }