Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright notices
[simgrid.git] / examples / smpi / tracing / smpi_traced.c
1 /* Copyright (c) 2009-2010, 2012-2015. 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 <stdio.h>
8 #include "mpi.h"
9 #include "simgrid/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     free(r);
163 /////////////////////////////////////////
164 ////////////////// RANK 1
165 ///////////////////////////////////
166   } else if (rank == 1) {
167     MPI_Request request;
168     MPI_Status status;
169     MPI_Request req[N];
170     MPI_Status sta[N];
171     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
172
173     if (A) {
174       TRACE_smpi_set_category("A");
175       MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
176     }
177     MPI_Barrier(MPI_COMM_WORLD);
178
179     if (B) {
180       TRACE_smpi_set_category("B");
181       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &request);
182       MPI_Wait(&request, &status);
183     }
184     MPI_Barrier(MPI_COMM_WORLD);
185
186     if (C) {
187       TRACE_smpi_set_category("C");
188       for (i = 0; i < N; i++) {
189         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
190       }
191     }
192     MPI_Barrier(MPI_COMM_WORLD);
193
194     if (D) {
195       TRACE_smpi_set_category("D");
196       for (i = 0; i < N; i++) {
197         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
198       }
199     }
200     MPI_Barrier(MPI_COMM_WORLD);
201
202     if (E) {
203       TRACE_smpi_set_category("E");
204       for (i = 0; i < N; i++) {
205         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
206       }
207       for (i = 0; i < N; i++) {
208         MPI_Wait(&req[i], &sta[i]);
209       }
210     }
211     MPI_Barrier(MPI_COMM_WORLD);
212
213     if (F) {
214       TRACE_smpi_set_category("F");
215       for (i = 0; i < N; i++) {
216         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
217       }
218       for (i = 0; i < N; i++) {
219         int completed;
220         MPI_Waitany(N, req, &completed, sta);
221       }
222     }
223     MPI_Barrier(MPI_COMM_WORLD);
224
225     if (G) {
226       TRACE_smpi_set_category("G");
227       for (i = 0; i < N; i++) {
228         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
229       }
230     }
231     MPI_Barrier(MPI_COMM_WORLD);
232
233     if (H) {
234       TRACE_smpi_set_category("H");
235       for (i = 0; i < N; i++) {
236         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
237       }
238       MPI_Waitall(N, req, sta);
239     }
240     MPI_Barrier(MPI_COMM_WORLD);
241
242     if (I) {
243       TRACE_smpi_set_category("I");
244       for (i = 0; i < N; i++) {
245         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
246       }
247       MPI_Waitall(N, req, sta);
248
249       MPI_Barrier(MPI_COMM_WORLD);
250       for (i = 0; i < N; i++) {
251         MPI_Isend(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
252       }
253       MPI_Waitall(N, req, sta);
254 //      for (i = 0; i < N; i++){
255 //        MPI_Wait (&req[i], &sta[i]);
256 //      }
257     }
258     MPI_Barrier(MPI_COMM_WORLD);
259
260     if (J) {
261       TRACE_smpi_set_category("J");
262       for (i = 0; i < N; i++) {
263         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
264       }
265       for (i = 0; i < N; i++) {
266         int flag;
267         MPI_Test(&req[i], &flag, &sta[i]);
268       }
269       for (i = 0; i < N; i++) {
270         MPI_Wait(&req[i], &sta[i]);
271       }
272     }
273     free(r);
274 /////////////////////////////////////////
275 ////////////////// RANK 2
276 ///////////////////////////////////
277   } else if (rank == 2) {
278 //    MPI_Request request;
279 //    MPI_Status status;
280     MPI_Request req[N];
281     MPI_Status sta[N];
282     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
283
284     if (A) {
285     }
286     MPI_Barrier(MPI_COMM_WORLD);
287     if (B) {
288     }
289     MPI_Barrier(MPI_COMM_WORLD);
290     if (C) {
291     }
292     MPI_Barrier(MPI_COMM_WORLD);
293     if (D) {
294     }
295     MPI_Barrier(MPI_COMM_WORLD);
296     if (E) {
297     }
298     MPI_Barrier(MPI_COMM_WORLD);
299     if (F) {
300     }
301     MPI_Barrier(MPI_COMM_WORLD);
302     if (G) {
303     }
304     MPI_Barrier(MPI_COMM_WORLD);
305     if (H) {
306     }
307     MPI_Barrier(MPI_COMM_WORLD);
308     if (I) {
309       TRACE_smpi_set_category("I");
310       for (i = 0; i < N; i++) {
311         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
312       }
313       for (i = 0; i < N; i++) {
314         int completed;
315         MPI_Waitany(N, req, &completed, sta);
316       }
317       MPI_Barrier(MPI_COMM_WORLD);
318
319       for (i = 0; i < N; i++) {
320         MPI_Send(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD);
321       }
322     }
323     MPI_Barrier(MPI_COMM_WORLD);
324     if (J) {
325     }
326     free(r);
327   }
328   MPI_Barrier(MPI_COMM_WORLD);
329   MPI_Finalize();
330   return 0;
331 }