Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Basic implementation of MPI_Cancel
[simgrid.git] / teshsuite / smpi / mpich3-test / pt2pt / scancel.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2003 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6 #include "mpi.h"
7 #include <stdio.h>
8 #include <stdlib.h>
9 #include "mpitest.h"
10
11 /*
12 static char MTEST_Descrip[] = "Test of various send cancel calls";
13 */
14
15 int main(int argc, char *argv[])
16 {
17     int errs = 0;
18     int rank, size, source, dest;
19     MPI_Comm comm;
20     MPI_Status status;
21     MPI_Request req;
22     static int bufsizes[4] = { 1, 100, 10000, 1000000 };
23     char *buf;
24 #ifdef TEST_IRSEND
25     int veryPicky = 0;          /* Set to 1 to test "quality of implementation" in
26                                  * a tricky part of cancel */
27 #endif
28     int cs, flag, n;
29
30     MTest_Init(&argc, &argv);
31
32     comm = MPI_COMM_WORLD;
33     MPI_Comm_rank(comm, &rank);
34     MPI_Comm_size(comm, &size);
35
36     source = 0;
37     dest = size - 1;
38
39     MTestPrintfMsg(1, "Starting scancel test\n");
40     for (cs = 0; cs < 4; cs++) {
41         if (rank == 0) {
42             n = bufsizes[cs];
43             buf = (char *) malloc(n);
44             if (!buf) {
45                 fprintf(stderr, "Unable to allocate %d bytes\n", n);
46                 MPI_Abort(MPI_COMM_WORLD, 1);
47             }
48             MTestPrintfMsg(1, "(%d) About to create isend and cancel\n", cs);
49             MPI_Isend(buf, n, MPI_CHAR, dest, cs + n + 1, comm, &req);
50             MPI_Cancel(&req);
51             MPI_Wait(&req, &status);
52             MTestPrintfMsg(1, "Completed wait on isend\n");
53             MPI_Test_cancelled(&status, &flag);
54             if (!flag) {
55                 errs++;
56                 printf("Failed to cancel an Isend request\n");
57                 fflush(stdout);
58             }
59             else {
60                 n = 0;
61             }
62             /* Send the size, zero for successfully cancelled */
63             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
64             /* Send the tag so the message can be received */
65             n = cs + n + 1;
66             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
67             free(buf);
68         }
69         else if (rank == dest) {
70             int nn, tag;
71             char *btemp;
72             MPI_Recv(&nn, 1, MPI_INT, 0, 123, comm, &status);
73             MPI_Recv(&tag, 1, MPI_INT, 0, 123, comm, &status);
74             if (nn > 0) {
75                 /* If the message was not cancelled, receive it here */
76                 btemp = (char *) malloc(nn);
77                 if (!btemp) {
78                     fprintf(stderr, "Unable to allocate %d bytes\n", nn);
79                     MPI_Abort(MPI_COMM_WORLD, 1);
80                 }
81                 MPI_Recv(btemp, nn, MPI_CHAR, 0, tag, comm, &status);
82                 free(btemp);
83             }
84         }
85         MPI_Barrier(comm);
86
87         if (rank == 0) {
88             char *bsendbuf;
89             int bsendbufsize;
90             int bf, bs;
91             n = bufsizes[cs];
92             buf = (char *) malloc(n);
93             if (!buf) {
94                 fprintf(stderr, "Unable to allocate %d bytes\n", n);
95                 MPI_Abort(MPI_COMM_WORLD, 1);
96             }
97             bsendbufsize = n + MPI_BSEND_OVERHEAD;
98             bsendbuf = (char *) malloc(bsendbufsize);
99             if (!bsendbuf) {
100                 fprintf(stderr, "Unable to allocate %d bytes for bsend\n", n);
101                 MPI_Abort(MPI_COMM_WORLD, 1);
102             }
103             MPI_Buffer_attach(bsendbuf, bsendbufsize);
104             MTestPrintfMsg(1, "About to create and cancel ibsend\n");
105             MPI_Ibsend(buf, n, MPI_CHAR, dest, cs + n + 2, comm, &req);
106             MPI_Cancel(&req);
107             MPI_Wait(&req, &status);
108             MPI_Test_cancelled(&status, &flag);
109             if (!flag) {
110                 errs++;
111                 printf("Failed to cancel an Ibsend request\n");
112                 fflush(stdout);
113             }
114             else {
115                 n = 0;
116             }
117             /* Send the size, zero for successfully cancelled */
118             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
119             /* Send the tag so the message can be received */
120             n = cs + n + 2;
121             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
122             free(buf);
123             MPI_Buffer_detach(&bf, &bs);
124             free(bsendbuf);
125         }
126         else if (rank == dest) {
127             int nn, tag;
128             char *btemp;
129             MPI_Recv(&nn, 1, MPI_INT, 0, 123, comm, &status);
130             MPI_Recv(&tag, 1, MPI_INT, 0, 123, comm, &status);
131             if (nn > 0) {
132                 /* If the message was not cancelled, receive it here */
133                 btemp = (char *) malloc(nn);
134                 if (!btemp) {
135                     fprintf(stderr, "Unable to allocate %d bytes\n", nn);
136                     MPI_Abort(MPI_COMM_WORLD, 1);
137                 }
138                 MPI_Recv(btemp, nn, MPI_CHAR, 0, tag, comm, &status);
139                 free(btemp);
140             }
141         }
142         MPI_Barrier(comm);
143
144         /* Because this test is erroneous, we do not perform it unless
145          * TEST_IRSEND is defined.  */
146 #ifdef TEST_IRSEND
147         /* We avoid ready send to self because an implementation
148          * is free to detect the error in delivering a message to
149          * itself without a pending receive; we could also check
150          * for an error return from the MPI_Irsend */
151         if (rank == 0 && dest != rank) {
152             n = bufsizes[cs];
153             buf = (char *) malloc(n);
154             if (!buf) {
155                 fprintf(stderr, "Unable to allocate %d bytes\n", n);
156                 MPI_Abort(MPI_COMM_WORLD, 1);
157             }
158             MTestPrintfMsg(1, "About to create and cancel irsend\n");
159             MPI_Irsend(buf, n, MPI_CHAR, dest, cs + n + 3, comm, &req);
160             MPI_Cancel(&req);
161             MPI_Wait(&req, &status);
162             MPI_Test_cancelled(&status, &flag);
163             /* This can be pretty ugly.  The standard is clear (Section 3.8)
164              * that either a sent message is received or the
165              * sent message is successfully cancelled.  Since this message
166              * can never be received, the cancel must complete
167              * successfully.
168              *
169              * However, since there is no matching receive, this
170              * program is erroneous.  In this case, we can't really
171              * flag this as an error */
172             if (!flag && veryPicky) {
173                 errs++;
174                 printf("Failed to cancel an Irsend request\n");
175                 fflush(stdout);
176             }
177             if (flag) {
178                 n = 0;
179             }
180             /* Send the size, zero for successfully cancelled */
181             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
182             /* Send the tag so the message can be received */
183             n = cs + n + 3;
184             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
185             free(buf);
186         }
187         else if (rank == dest) {
188             int n, tag;
189             char *btemp;
190             MPI_Recv(&n, 1, MPI_INT, 0, 123, comm, &status);
191             MPI_Recv(&tag, 1, MPI_INT, 0, 123, comm, &status);
192             if (n > 0) {
193                 /* If the message was not cancelled, receive it here */
194                 btemp = (char *) malloc(n);
195                 if (!btemp) {
196                     fprintf(stderr, "Unable to allocate %d bytes\n", n);
197                     MPI_Abort(MPI_COMM_WORLD, 1);
198                 }
199                 MPI_Recv(btemp, n, MPI_CHAR, 0, tag, comm, &status);
200                 free(btemp);
201             }
202         }
203         MPI_Barrier(comm);
204 #endif
205
206         if (rank == 0) {
207             n = bufsizes[cs];
208             buf = (char *) malloc(n);
209             if (!buf) {
210                 fprintf(stderr, "Unable to allocate %d bytes\n", n);
211                 MPI_Abort(MPI_COMM_WORLD, 1);
212             }
213             MTestPrintfMsg(1, "About to create and cancel issend\n");
214             MPI_Issend(buf, n, MPI_CHAR, dest, cs + n + 4, comm, &req);
215             MPI_Cancel(&req);
216             MPI_Wait(&req, &status);
217             MPI_Test_cancelled(&status, &flag);
218             if (!flag) {
219                 errs++;
220                 printf("Failed to cancel an Issend request\n");
221                 fflush(stdout);
222             }
223             else {
224                 n = 0;
225             }
226             /* Send the size, zero for successfully cancelled */
227             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
228             /* Send the tag so the message can be received */
229             n = cs + n + 4;
230             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
231             free(buf);
232         }
233         else if (rank == dest) {
234             int nn, tag;
235             char *btemp;
236             MPI_Recv(&nn, 1, MPI_INT, 0, 123, comm, &status);
237             MPI_Recv(&tag, 1, MPI_INT, 0, 123, comm, &status);
238             if (nn > 0) {
239                 /* If the message was not cancelled, receive it here */
240                 btemp = (char *) malloc(nn);
241                 if (!btemp) {
242                     fprintf(stderr, "Unable to allocate %d bytes\n", nn);
243                     MPI_Abort(MPI_COMM_WORLD, 1);
244                 }
245                 MPI_Recv(btemp, nn, MPI_CHAR, 0, tag, comm, &status);
246                 free(btemp);
247             }
248         }
249         MPI_Barrier(comm);
250     }
251
252     MTest_Finalize(errs);
253     MPI_Finalize();
254     return 0;
255 }