Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines with new year.
[simgrid.git] / examples / smpi / trace / trace.c
1 /* Copyright (c) 2009-2019. 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 #include <unistd.h>
11
12 #define DATATOSENT 100000
13
14 int main(int argc, char *argv[])
15 {
16   int rank;
17   int numprocs;
18
19   MPI_Init(&argc, &argv);
20   MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
21   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
22   sleep(rank);
23   /** Tests:
24   * A: 0(isend -> wait) with 1(recv)
25   * B: 0(send) with 1(irecv -> wait)
26   * C: 0(N * isend -> N * wait) with 1(N * recv)
27   * D: 0(N * isend -> N * waitany) with 1(N * recv)
28   * E: 0(N*send) with 1(N*irecv, N*wait)
29   * F: 0(N*send) with 1(N*irecv, N*waitany)
30   * G: 0(N* isend -> waitall) with 1(N*recv)
31   * H: 0(N*send) with 1(N*irecv, waitall)
32   * I: 0(2*N*send, 2*N*Irecv, Waitall) with
33   *    1(N*irecv, waitall, N*isend, N*waitany) with
34   *    2(N*irecv, N*waitany, N*isend, waitall)
35   * J: 0(N*isend, N*test, N*wait) with (N*irecv, N*test, N*wait)
36   s*/
37   int N = 13;
38
39   int tag = 12345;
40 /////////////////////////////////////////
41 ////////////////// RANK 0
42 ///////////////////////////////////
43   if (rank == 0) {
44     MPI_Request request;
45     MPI_Status status;
46     MPI_Request req[2 * N];
47     MPI_Status sta[2 * N];
48     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
49
50     /** Test A */
51     TRACE_smpi_set_category("A");
52     MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &request);
53     MPI_Wait(&request, &status);
54     MPI_Barrier(MPI_COMM_WORLD);
55
56     /** Test B */
57     TRACE_smpi_set_category("B");
58     MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
59     MPI_Barrier(MPI_COMM_WORLD);
60
61     /** Test C */
62     TRACE_smpi_set_category("C");
63     for (int i = 0; i < N; i++) {
64       MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
65     }
66     for (int i = 0; i < N; i++) {
67       MPI_Wait(&req[i], &sta[i]);
68     }
69     MPI_Barrier(MPI_COMM_WORLD);
70
71     TRACE_smpi_set_category("D");
72     for (int i = 0; i < N; i++) {
73       MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
74     }
75     for (int i = 0; i < N; i++) {
76       int completed;
77       MPI_Waitany(N, req, &completed, sta);
78     }
79     MPI_Barrier(MPI_COMM_WORLD);
80
81     TRACE_smpi_set_category("E");
82     for (int i = 0; i < N; i++) {
83       MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
84     }
85     MPI_Barrier(MPI_COMM_WORLD);
86
87     TRACE_smpi_set_category("F");
88     for (int i = 0; i < N; i++) {
89       MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
90     }
91     MPI_Barrier(MPI_COMM_WORLD);
92
93     TRACE_smpi_set_category("G");
94     for (int i = 0; i < N; i++) {
95       MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
96     }
97     MPI_Waitall(N, req, sta);
98     MPI_Barrier(MPI_COMM_WORLD);
99
100     TRACE_smpi_set_category("H");
101     for (int i = 0; i < N; i++) {
102       MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
103     }
104     MPI_Barrier(MPI_COMM_WORLD);
105
106     TRACE_smpi_set_category("I");
107     for (int i = 0; i < 2 * N; i++) {
108       if (i < N) {
109         MPI_Send(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD);
110       } else {
111         MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
112       }
113     }
114     MPI_Barrier(MPI_COMM_WORLD);
115     for (int i = 0; i < 2 * N; i++) {
116       if (i < N) {
117         MPI_Irecv(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
118       } else {
119         MPI_Irecv(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD, &req[i]);
120       }
121     }
122     MPI_Waitall(2 * N, req, sta);
123     MPI_Barrier(MPI_COMM_WORLD);
124
125     TRACE_smpi_set_category("J");
126     for (int i = 0; i < N; i++) {
127       MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
128     }
129     for (int i = 0; i < N; i++) {
130       int flag;
131       MPI_Test(&req[i], &flag, &sta[i]);
132     }
133     for (int i = 0; i < N; i++) {
134       MPI_Wait(&req[i], &sta[i]);
135     }
136     free(r);
137 /////////////////////////////////////////
138 ////////////////// RANK 1
139 ///////////////////////////////////
140   } else if (rank == 1) {
141     MPI_Request request;
142     MPI_Status status;
143     MPI_Request req[N];
144     MPI_Status sta[N];
145     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
146
147     TRACE_smpi_set_category("A");
148     MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
149     MPI_Barrier(MPI_COMM_WORLD);
150
151     TRACE_smpi_set_category("B");
152     MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &request);
153     MPI_Wait(&request, &status);
154     MPI_Barrier(MPI_COMM_WORLD);
155
156     TRACE_smpi_set_category("C");
157     for (int i = 0; i < N; i++) {
158       MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
159     }
160     MPI_Barrier(MPI_COMM_WORLD);
161
162     TRACE_smpi_set_category("D");
163     for (int i = 0; i < N; i++) {
164       MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
165     }
166     MPI_Barrier(MPI_COMM_WORLD);
167
168     TRACE_smpi_set_category("E");
169     for (int i = 0; i < N; i++) {
170       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
171     }
172     for (int i = 0; i < N; i++) {
173       MPI_Wait(&req[i], &sta[i]);
174     }
175     MPI_Barrier(MPI_COMM_WORLD);
176
177     TRACE_smpi_set_category("F");
178     for (int i = 0; i < N; i++) {
179       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
180     }
181     for (int i = 0; i < N; i++) {
182       int completed;
183       MPI_Waitany(N, req, &completed, sta);
184     }
185     MPI_Barrier(MPI_COMM_WORLD);
186
187     TRACE_smpi_set_category("G");
188     for (int i = 0; i < N; i++) {
189       MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
190     }
191     MPI_Barrier(MPI_COMM_WORLD);
192
193     TRACE_smpi_set_category("H");
194     for (int i = 0; i < N; i++) {
195       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
196     }
197     MPI_Waitall(N, req, sta);
198     MPI_Barrier(MPI_COMM_WORLD);
199
200     TRACE_smpi_set_category("I");
201     for (int i = 0; i < N; i++) {
202       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
203     }
204     MPI_Waitall(N, req, sta);
205
206     MPI_Barrier(MPI_COMM_WORLD);
207     for (int i = 0; i < N; i++) {
208       MPI_Isend(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
209     }
210     MPI_Waitall(N, req, sta);
211     MPI_Barrier(MPI_COMM_WORLD);
212
213     TRACE_smpi_set_category("J");
214     for (int i = 0; i < N; i++) {
215       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
216     }
217     for (int i = 0; i < N; i++) {
218       int flag;
219       MPI_Test(&req[i], &flag, &sta[i]);
220     }
221     for (int i = 0; i < N; i++) {
222       MPI_Wait(&req[i], &sta[i]);
223     }
224     free(r);
225 /////////////////////////////////////////
226 ////////////////// RANK 2
227 ///////////////////////////////////
228   } else if (rank == 2) {
229     MPI_Request req[N];
230     MPI_Status sta[N];
231     int *r = (int *) malloc(sizeof(int) * DATATOSENT);
232
233     MPI_Barrier(MPI_COMM_WORLD);
234     MPI_Barrier(MPI_COMM_WORLD);
235     MPI_Barrier(MPI_COMM_WORLD);
236     MPI_Barrier(MPI_COMM_WORLD);
237     MPI_Barrier(MPI_COMM_WORLD);
238     MPI_Barrier(MPI_COMM_WORLD);
239     MPI_Barrier(MPI_COMM_WORLD);
240     MPI_Barrier(MPI_COMM_WORLD);
241     TRACE_smpi_set_category("I");
242     for (int i = 0; i < N; i++) {
243       MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
244     }
245     for (int i = 0; i < N; i++) {
246       int completed;
247       MPI_Waitany(N, req, &completed, sta);
248     }
249     MPI_Barrier(MPI_COMM_WORLD);
250
251     for (int i = 0; i < N; i++) {
252       MPI_Send(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD);
253     }
254     MPI_Barrier(MPI_COMM_WORLD);
255     free(r);
256   }
257   MPI_Barrier(MPI_COMM_WORLD);
258   MPI_Finalize();
259   return 0;
260 }