1 /* Copyright (c) 2009-2010, 2012-2015. The SimGrid Team.
2 * All rights reserved. */
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. */
9 #include "simgrid/instr.h"
11 #define DATATOSENT 100000
13 int main(int argc, char *argv[])
18 MPI_Init(&argc, &argv);
19 MPI_Comm_size(MPI_COMM_WORLD, &numprocs);
20 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
23 * A: 0(isend -> wait) with 1(recv)
24 * B: 0(send) with 1(irecv -> wait)
25 * C: 0(N * isend -> N * wait) with 1(N * recv)
26 * D: 0(N * isend -> N * waitany) with 1(N * recv)
27 * E: 0(N*send) with 1(N*irecv, N*wait)
28 * F: 0(N*send) with 1(N*irecv, N*waitany)
29 * G: 0(N* isend -> waitall) with 1(N*recv)
30 * H: 0(N*send) with 1(N*irecv, waitall)
31 * I: 0(2*N*send, 2*N*Irecv, Waitall) with
32 * 1(N*irecv, waitall, N*isend, N*waitany) with
33 * 2(N*irecv, N*waitany, N*isend, waitall)
34 * J: 0(N*isend, N*test, N*wait) with (N*irecv, N*test, N*wait)
39 /////////////////////////////////////////
40 ////////////////// RANK 0
41 ///////////////////////////////////
45 MPI_Request req[2 * N];
46 MPI_Status sta[2 * N];
47 int *r = (int *) malloc(sizeof(int) * DATATOSENT);
50 TRACE_smpi_set_category("A");
51 MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &request);
52 MPI_Wait(&request, &status);
53 MPI_Barrier(MPI_COMM_WORLD);
56 TRACE_smpi_set_category("B");
57 MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
58 MPI_Barrier(MPI_COMM_WORLD);
61 TRACE_smpi_set_category("C");
62 for (int i = 0; i < N; i++) {
63 MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
65 for (int i = 0; i < N; i++) {
66 MPI_Wait(&req[i], &sta[i]);
68 MPI_Barrier(MPI_COMM_WORLD);
70 TRACE_smpi_set_category("D");
71 for (int i = 0; i < N; i++) {
72 MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
74 for (int i = 0; i < N; i++) {
76 MPI_Waitany(N, req, &completed, sta);
78 MPI_Barrier(MPI_COMM_WORLD);
80 TRACE_smpi_set_category("E");
81 for (int i = 0; i < N; i++) {
82 MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
84 MPI_Barrier(MPI_COMM_WORLD);
86 TRACE_smpi_set_category("F");
87 for (int i = 0; i < N; i++) {
88 MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
90 MPI_Barrier(MPI_COMM_WORLD);
92 TRACE_smpi_set_category("G");
93 for (int i = 0; i < N; i++) {
94 MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
96 MPI_Waitall(N, req, sta);
97 MPI_Barrier(MPI_COMM_WORLD);
99 TRACE_smpi_set_category("H");
100 for (int i = 0; i < N; i++) {
101 MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
103 MPI_Barrier(MPI_COMM_WORLD);
105 TRACE_smpi_set_category("I");
106 for (int i = 0; i < 2 * N; i++) {
108 MPI_Send(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD);
110 MPI_Send(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD);
113 MPI_Barrier(MPI_COMM_WORLD);
114 for (int i = 0; i < 2 * N; i++) {
116 MPI_Irecv(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
118 MPI_Irecv(r, DATATOSENT, MPI_INT, 2, tag, MPI_COMM_WORLD, &req[i]);
121 MPI_Waitall(2 * N, req, sta);
122 MPI_Barrier(MPI_COMM_WORLD);
124 TRACE_smpi_set_category("J");
125 for (int i = 0; i < N; i++) {
126 MPI_Isend(r, DATATOSENT, MPI_INT, 1, tag, MPI_COMM_WORLD, &req[i]);
128 for (int i = 0; i < N; i++) {
130 MPI_Test(&req[i], &flag, &sta[i]);
132 for (int i = 0; i < N; i++) {
133 MPI_Wait(&req[i], &sta[i]);
136 /////////////////////////////////////////
137 ////////////////// RANK 1
138 ///////////////////////////////////
139 } else if (rank == 1) {
144 int *r = (int *) malloc(sizeof(int) * DATATOSENT);
146 TRACE_smpi_set_category("A");
147 MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
148 MPI_Barrier(MPI_COMM_WORLD);
150 TRACE_smpi_set_category("B");
151 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &request);
152 MPI_Wait(&request, &status);
153 MPI_Barrier(MPI_COMM_WORLD);
155 TRACE_smpi_set_category("C");
156 for (int i = 0; i < N; i++) {
157 MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
159 MPI_Barrier(MPI_COMM_WORLD);
161 TRACE_smpi_set_category("D");
162 for (int i = 0; i < N; i++) {
163 MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
165 MPI_Barrier(MPI_COMM_WORLD);
167 TRACE_smpi_set_category("E");
168 for (int i = 0; i < N; i++) {
169 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
171 for (int i = 0; i < N; i++) {
172 MPI_Wait(&req[i], &sta[i]);
174 MPI_Barrier(MPI_COMM_WORLD);
176 TRACE_smpi_set_category("F");
177 for (int i = 0; i < N; i++) {
178 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
180 for (int i = 0; i < N; i++) {
182 MPI_Waitany(N, req, &completed, sta);
184 MPI_Barrier(MPI_COMM_WORLD);
186 TRACE_smpi_set_category("G");
187 for (int i = 0; i < N; i++) {
188 MPI_Recv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &sta[i]);
190 MPI_Barrier(MPI_COMM_WORLD);
192 TRACE_smpi_set_category("H");
193 for (int i = 0; i < N; i++) {
194 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
196 MPI_Waitall(N, req, sta);
197 MPI_Barrier(MPI_COMM_WORLD);
199 TRACE_smpi_set_category("I");
200 for (int i = 0; i < N; i++) {
201 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
203 MPI_Waitall(N, req, sta);
205 MPI_Barrier(MPI_COMM_WORLD);
206 for (int i = 0; i < N; i++) {
207 MPI_Isend(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
209 MPI_Waitall(N, req, sta);
210 MPI_Barrier(MPI_COMM_WORLD);
212 TRACE_smpi_set_category("J");
213 for (int i = 0; i < N; i++) {
214 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
216 for (int i = 0; i < N; i++) {
218 MPI_Test(&req[i], &flag, &sta[i]);
220 for (int i = 0; i < N; i++) {
221 MPI_Wait(&req[i], &sta[i]);
224 /////////////////////////////////////////
225 ////////////////// RANK 2
226 ///////////////////////////////////
227 } else if (rank == 2) {
230 int *r = (int *) malloc(sizeof(int) * DATATOSENT);
232 MPI_Barrier(MPI_COMM_WORLD);
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 TRACE_smpi_set_category("I");
241 for (int i = 0; i < N; i++) {
242 MPI_Irecv(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD, &req[i]);
244 for (int i = 0; i < N; i++) {
246 MPI_Waitany(N, req, &completed, sta);
248 MPI_Barrier(MPI_COMM_WORLD);
250 for (int i = 0; i < N; i++) {
251 MPI_Send(r, DATATOSENT, MPI_INT, 0, tag, MPI_COMM_WORLD);
253 MPI_Barrier(MPI_COMM_WORLD);
256 MPI_Barrier(MPI_COMM_WORLD);