Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
using categories to trace platform utilization on smpi tracing example
[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 100000000
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, &req[i]);
139         }else{
140           MPI_Irecv(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD, &req[i]);
141         }
142       }
143       MPI_Waitall (2*N, req, sta);
144     }
145     MPI_Barrier (MPI_COMM_WORLD);
146
147     if (J){
148       TRACE_smpi_set_category ("J");
149       for (i = 0; i < N; i++){
150         MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
151       }
152       for (i = 0; i < N; i++){
153         int flag;
154         MPI_Test (&req[i], &flag, &sta[i]);
155       }
156       for (i = 0; i < N; i++){
157         MPI_Wait (&req[i], &sta[i]);
158       }
159     }
160
161 /////////////////////////////////////////
162 ////////////////// RANK 1
163 ///////////////////////////////////
164   }else if (rank == 1){
165     MPI_Request request;
166     MPI_Status status;
167     MPI_Request req[N];
168     MPI_Status sta[N];
169     int *r = (int*)malloc(sizeof(int)*DATATOSENT);
170
171     if (A){
172       TRACE_smpi_set_category ("A");
173       MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
174     }
175     MPI_Barrier (MPI_COMM_WORLD);
176
177     if (B){
178       TRACE_smpi_set_category ("B");
179       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &request);
180       MPI_Wait (&request, &status);
181     }
182     MPI_Barrier (MPI_COMM_WORLD);
183
184     if (C){
185       TRACE_smpi_set_category ("C");
186       for (i = 0; i < N; i++){
187         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
188       }
189     }
190     MPI_Barrier (MPI_COMM_WORLD);
191
192     if (D){
193       TRACE_smpi_set_category ("D");
194       for (i = 0; i < N; i++){
195         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
196       }
197     }
198     MPI_Barrier (MPI_COMM_WORLD);
199
200     if (E){
201       TRACE_smpi_set_category ("E");
202       for (i = 0; i < N; i++){
203         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
204       }
205       for (i = 0; i < N; i++){
206         MPI_Wait (&req[i], &sta[i]);
207       }
208     }
209     MPI_Barrier (MPI_COMM_WORLD);
210
211     if (F){
212       TRACE_smpi_set_category ("F");
213       for (i = 0; i < N; i++){
214         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
215       }
216       for (i = 0; i < N; i++){
217         int completed;
218         MPI_Waitany (N, req, &completed, sta);
219       }
220     }
221     MPI_Barrier (MPI_COMM_WORLD);
222
223     if (G){
224       TRACE_smpi_set_category ("G");
225       for (i = 0; i < N; i++){
226         MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
227       }
228     }
229     MPI_Barrier (MPI_COMM_WORLD);
230
231     if (H){
232       TRACE_smpi_set_category ("H");
233       for (i = 0; i < N; i++){
234         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
235       }
236       MPI_Waitall (N, req, sta);
237     }
238     MPI_Barrier (MPI_COMM_WORLD);
239
240     if (I){
241       TRACE_smpi_set_category ("I");
242       for (i = 0; i < N; i++){
243         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
244       }
245       MPI_Waitall (N, req, sta);
246
247       MPI_Barrier (MPI_COMM_WORLD);
248       for (i = 0; i < N; i++){
249         MPI_Isend(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
250       }
251       MPI_Waitall (N, req, sta);
252 //      for (i = 0; i < N; i++){
253 //        MPI_Wait (&req[i], &sta[i]);
254 //      }
255     }
256     MPI_Barrier (MPI_COMM_WORLD);
257
258     if (J){
259       TRACE_smpi_set_category ("J");
260       for (i = 0; i < N; i++){
261         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
262       }
263       for (i = 0; i < N; i++){
264         int flag;
265         MPI_Test (&req[i], &flag, &sta[i]);
266       }
267       for (i = 0; i < N; i++){
268         MPI_Wait (&req[i], &sta[i]);
269       }
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     MPI_Barrier (MPI_COMM_WORLD);
284     if (B) {} 
285     MPI_Barrier (MPI_COMM_WORLD);
286     if (C) {} 
287     MPI_Barrier (MPI_COMM_WORLD);
288     if (D) {} 
289     MPI_Barrier (MPI_COMM_WORLD);
290     if (E) {} 
291     MPI_Barrier (MPI_COMM_WORLD);
292     if (F) {} 
293     MPI_Barrier (MPI_COMM_WORLD);
294     if (G) {} 
295     MPI_Barrier (MPI_COMM_WORLD);
296     if (H) {} 
297     MPI_Barrier (MPI_COMM_WORLD);
298     if (I){
299       TRACE_smpi_set_category ("I");
300       for (i = 0; i < N; i++){
301         MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
302       }
303       for (i = 0; i < N; i++){
304         int completed;
305         MPI_Waitany (N, req, &completed, sta);
306       }
307       MPI_Barrier (MPI_COMM_WORLD);
308
309       for (i = 0; i < N; i++){
310         MPI_Send(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD);
311       }
312     }
313     MPI_Barrier (MPI_COMM_WORLD);
314     if (J){}
315   }
316   MPI_Barrier (MPI_COMM_WORLD);
317   MPI_Finalize();
318   return 0;
319 }