Logo AND Algorithmique Numérique Distribuée

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