Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Fix warnings about unused variables in mpich3-test.
[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             {
59                 n = 0;
60             }
61             MPI_Request_free( &req );
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         {
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             {
77                 /* If the message was not cancelled, receive it here */
78                 btemp = (char*)malloc( nn );
79                 if (!btemp)
80                 {
81                     fprintf( stderr, "Unable to allocate %d bytes\n", nn);
82                     MPI_Abort( MPI_COMM_WORLD, 1 );
83                 }
84                 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
85                 free(btemp);
86             }
87         }
88         MPI_Barrier( comm );
89
90         if (rank == 0) {
91             char *bsendbuf;
92             int bsendbufsize;
93             int bf, bs;
94             n = bufsizes[cs];
95             buf = (char *)malloc( n );
96             if (!buf) {
97                 fprintf( stderr, "Unable to allocate %d bytes\n", n );
98                 MPI_Abort( MPI_COMM_WORLD, 1 );
99             }
100             bsendbufsize = n + MPI_BSEND_OVERHEAD;
101             bsendbuf = (char *)malloc( bsendbufsize );
102             if (!bsendbuf) {
103                 fprintf( stderr, "Unable to allocate %d bytes for bsend\n", n );
104                 MPI_Abort( MPI_COMM_WORLD, 1 );
105             }
106             MPI_Buffer_attach( bsendbuf, bsendbufsize );
107             MPI_Bsend_init( buf, n, MPI_CHAR, dest, cs+n+2, comm, &req );
108             MPI_Start( &req );
109             MPI_Cancel( &req );
110             MPI_Wait( &req, &status );
111             MPI_Test_cancelled( &status, &flag );
112             if (!flag) {
113                 errs ++;
114                 printf( "Failed to cancel a persistent bsend request\n" );
115                 fflush(stdout);
116             }
117             else
118             {
119                 n = 0;
120             }
121             MPI_Request_free( &req );
122             /* Send the size, zero for successfully cancelled */
123             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
124             /* Send the tag so the message can be received */
125             n = cs+n+2;
126             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
127             free( buf );
128             MPI_Buffer_detach( &bf, &bs );
129             free( bsendbuf );
130         }
131         else if (rank == dest)
132         {
133             int nn, tag;
134             char *btemp;
135             MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
136             MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
137             if (nn > 0)
138             {
139                 /* If the message was not cancelled, receive it here */
140                 btemp = (char*)malloc( nn );
141                 if (!btemp)
142                 {
143                     fprintf( stderr, "Unable to allocate %d bytes\n", nn);
144                     MPI_Abort( MPI_COMM_WORLD, 1 );
145                 }
146                 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
147                 free(btemp);
148             }
149         }
150         MPI_Barrier( comm );
151
152         /* Because this test is erroneous, we do not perform it unless
153            TEST_IRSEND is defined.  */
154 #ifdef TEST_IRSEND
155         /* We avoid ready send to self because an implementation
156            is free to detect the error in delivering a message to
157            itself without a pending receive; we could also check
158            for an error return from the MPI_Irsend */
159         if (rank == 0 && dest != rank) {
160             n = bufsizes[cs];
161             buf = (char *)malloc( n );
162             if (!buf) {
163                 fprintf( stderr, "Unable to allocate %d bytes\n", n );
164                 MPI_Abort( MPI_COMM_WORLD, 1 );
165             }
166             MPI_Rsend_init( buf, n, MPI_CHAR, dest, cs+n+3, comm, &req );
167             MPI_Start( &req );
168             MPI_Cancel( &req );
169             MPI_Wait( &req, &status );
170             MPI_Test_cancelled( &status, &flag );
171             /* This can be pretty ugly.  The standard is clear (Section 3.8)
172                that either a sent message is received or the 
173                sent message is successfully cancelled.  Since this message
174                can never be received, the cancel must complete
175                successfully.  
176
177                However, since there is no matching receive, this
178                program is erroneous.  In this case, we can't really
179                flag this as an error */
180             if (!flag && veryPicky) {
181                 errs ++;
182                 printf( "Failed to cancel a persistent rsend request\n" );
183                 fflush(stdout);
184             }
185             if (flag)
186             {
187                 n = 0;
188             }
189             MPI_Request_free( &req );
190             /* Send the size, zero for successfully cancelled */
191             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
192             /* Send the tag so the message can be received */
193             n = cs+n+3;
194             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
195             free( buf );
196         }
197         else if (rank == dest)
198         {
199             int n, tag;
200             char *btemp;
201             MPI_Recv( &n, 1, MPI_INT, 0, 123, comm, &status );
202             MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
203             if (n > 0)
204             {
205                 /* If the message was not cancelled, receive it here */
206                 btemp = (char*)malloc( n );
207                 if (!btemp)
208                 {
209                     fprintf( stderr, "Unable to allocate %d bytes\n", n);
210                     MPI_Abort( MPI_COMM_WORLD, 1 );
211                 }
212                 MPI_Recv( btemp, n, MPI_CHAR, 0, tag, comm, &status );
213                 free(btemp);
214             }
215         }
216         MPI_Barrier( comm );
217 #endif
218
219         if (rank == 0) {
220             n = bufsizes[cs];
221             buf = (char *)malloc( n );
222             if (!buf) {
223                 fprintf( stderr, "Unable to allocate %d bytes\n", n );
224                 MPI_Abort( MPI_COMM_WORLD, 1 );
225             }
226             MPI_Ssend_init( buf, n, MPI_CHAR, dest, cs+n+4, comm, &req );
227             MPI_Start( &req );
228             MPI_Cancel( &req );
229             MPI_Wait( &req, &status );
230             MPI_Test_cancelled( &status, &flag );
231             if (!flag) {
232                 errs ++;
233                 printf( "Failed to cancel a persistent ssend request\n" );
234                 fflush(stdout);
235             }
236             else
237             {
238                 n = 0;
239             }
240             MPI_Request_free( &req );
241             /* Send the size, zero for successfully cancelled */
242             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
243             /* Send the tag so the message can be received */
244             n = cs+n+4;
245             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
246             free( buf );
247         }
248         else if (rank == dest)
249         {
250             int nn, tag;
251             char *btemp;
252             MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
253             MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
254             if (nn > 0)
255             {
256                 /* If the message was not cancelled, receive it here */
257                 btemp = (char*)malloc( nn );
258                 if (!btemp)
259                 {
260                     fprintf( stderr, "Unable to allocate %d bytes\n", nn);
261                     MPI_Abort( MPI_COMM_WORLD, 1 );
262                 }
263                 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
264                 free(btemp);
265             }
266         }
267         MPI_Barrier( comm );
268     }
269
270     MTest_Finalize( errs );
271     MPI_Finalize();
272     return 0;
273 }