Logo AND Algorithmique Numérique Distribuée

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