Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Reduce the size of partial shared malloc tests.
[simgrid.git] / teshsuite / smpi / mpich3-test / pt2pt / pscancel.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     int cs, flag, n;
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
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     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             MPI_Send_init(buf, n, MPI_CHAR, dest, cs + n + 1, comm, &req);
48             MPI_Start(&req);
49             MPI_Cancel(&req);
50             MPI_Wait(&req, &status);
51             MPI_Test_cancelled(&status, &flag);
52             if (!flag) {
53                 errs++;
54                 printf("Failed to cancel a persistent send request\n");
55                 fflush(stdout);
56             }
57             else {
58                 n = 0;
59             }
60             MPI_Request_free(&req);
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
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             MPI_Bsend_init(buf, n, MPI_CHAR, dest, cs + n + 2, comm, &req);
104             MPI_Start(&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 a persistent bsend request\n");
111                 fflush(stdout);
112             }
113             else {
114                 n = 0;
115             }
116             MPI_Request_free(&req);
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             MPI_Rsend_init(buf, n, MPI_CHAR, dest, cs + n + 3, comm, &req);
159             MPI_Start(&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 a persistent rsend request\n");
175                 fflush(stdout);
176             }
177             if (flag) {
178                 n = 0;
179             }
180             MPI_Request_free(&req);
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             MPI_Ssend_init(buf, n, MPI_CHAR, dest, cs + n + 4, comm, &req);
215             MPI_Start(&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 a persistent ssend request\n");
222                 fflush(stdout);
223             }
224             else {
225                 n = 0;
226             }
227             MPI_Request_free(&req);
228             /* Send the size, zero for successfully cancelled */
229             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
230             /* Send the tag so the message can be received */
231             n = cs + n + 4;
232             MPI_Send(&n, 1, MPI_INT, dest, 123, comm);
233             free(buf);
234         }
235         else if (rank == dest) {
236             int nn, tag;
237             char *btemp;
238             MPI_Recv(&nn, 1, MPI_INT, 0, 123, comm, &status);
239             MPI_Recv(&tag, 1, MPI_INT, 0, 123, comm, &status);
240             if (nn > 0) {
241                 /* If the message was not cancelled, receive it here */
242                 btemp = (char *) malloc(nn);
243                 if (!btemp) {
244                     fprintf(stderr, "Unable to allocate %d bytes\n", nn);
245                     MPI_Abort(MPI_COMM_WORLD, 1);
246                 }
247                 MPI_Recv(btemp, nn, MPI_CHAR, 0, tag, comm, &status);
248                 free(btemp);
249             }
250         }
251         MPI_Barrier(comm);
252     }
253
254     MTest_Finalize(errs);
255     MPI_Finalize();
256     return 0;
257 }