Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge tag 'v3_9_90' into hypervisor
[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             {
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         {
72             int nn, tag;
73             char *btemp;
74             MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
75             MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
76             if (nn > 0)
77             {
78                 /* If the message was not cancelled, receive it here */
79                 btemp = (char*)malloc( nn );
80                 if (!btemp)
81                 {
82                     fprintf( stderr, "Unable to allocate %d bytes\n", nn );
83                     MPI_Abort( MPI_COMM_WORLD, 1 );
84                 }
85                 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
86                 free(btemp);
87             }
88         }
89         MPI_Barrier( comm );
90
91         if (rank == 0) {
92             char *bsendbuf;
93             int bsendbufsize;
94             int bf, bs;
95             n = bufsizes[cs];
96             buf = (char *)malloc( n );
97             if (!buf) {
98                 fprintf( stderr, "Unable to allocate %d bytes\n", n );
99                 MPI_Abort( MPI_COMM_WORLD, 1 );
100             }
101             bsendbufsize = n + MPI_BSEND_OVERHEAD;
102             bsendbuf = (char *)malloc( bsendbufsize );
103             if (!bsendbuf) {
104                 fprintf( stderr, "Unable to allocate %d bytes for bsend\n", n );
105                 MPI_Abort( MPI_COMM_WORLD, 1 );
106             }
107             MPI_Buffer_attach( bsendbuf, bsendbufsize );
108             MTestPrintfMsg( 1, "About to create and cancel ibsend\n" );
109             MPI_Ibsend( buf, n, MPI_CHAR, dest, cs+n+2, comm, &req );
110             MPI_Cancel( &req );
111             MPI_Wait( &req, &status );
112             MPI_Test_cancelled( &status, &flag );
113             if (!flag) {
114                 errs ++;
115                 printf( "Failed to cancel an Ibsend request\n" );
116                 fflush(stdout);
117             }
118             else
119             {
120                 n = 0;
121             }
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             MTestPrintfMsg( 1, "About to create and cancel irsend\n" );
167             MPI_Irsend( buf, n, MPI_CHAR, dest, cs+n+3, comm, &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 an Irsend request\n" );
183                 fflush(stdout);
184             }
185             if (flag)
186             {
187                 n = 0;
188             }
189             /* Send the size, zero for successfully cancelled */
190             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
191             /* Send the tag so the message can be received */
192             n = cs+n+3;
193             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
194             free( buf );
195         }
196         else if (rank == dest)
197         {
198             int n, tag;
199             char *btemp;
200             MPI_Recv( &n, 1, MPI_INT, 0, 123, comm, &status );
201             MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
202             if (n > 0)
203             {
204                 /* If the message was not cancelled, receive it here */
205                 btemp = (char*)malloc( n );
206                 if (!btemp)
207                 {
208                     fprintf( stderr, "Unable to allocate %d bytes\n", n);
209                     MPI_Abort( MPI_COMM_WORLD, 1 );
210                 }
211                 MPI_Recv( btemp, n, MPI_CHAR, 0, tag, comm, &status );
212                 free(btemp);
213             }
214         }
215         MPI_Barrier( comm );
216 #endif
217
218         if (rank == 0) {
219             n = bufsizes[cs];
220             buf = (char *)malloc( n );
221             if (!buf) {
222                 fprintf( stderr, "Unable to allocate %d bytes\n", n );
223                 MPI_Abort( MPI_COMM_WORLD, 1 );
224             }
225             MTestPrintfMsg( 1, "About to create and cancel issend\n" );
226             MPI_Issend( buf, n, MPI_CHAR, dest, cs+n+4, comm, &req );
227             MPI_Cancel( &req );
228             MPI_Wait( &req, &status );
229             MPI_Test_cancelled( &status, &flag );
230             if (!flag) {
231                 errs ++;
232                 printf( "Failed to cancel an Issend request\n" );
233                 fflush(stdout);
234             }
235             else
236             {
237                 n = 0;
238             }
239             /* Send the size, zero for successfully cancelled */
240             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
241             /* Send the tag so the message can be received */
242             n = cs+n+4;
243             MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
244             free( buf );
245         }
246         else if (rank == dest)
247         {
248             int nn, tag;
249             char *btemp;
250             MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
251             MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
252             if (nn > 0)
253             {
254                 /* If the message was not cancelled, receive it here */
255                 btemp = (char*)malloc( nn );
256                 if (!btemp)
257                 {
258                     fprintf( stderr, "Unable to allocate %d bytes\n", nn);
259                     MPI_Abort( MPI_COMM_WORLD, 1 );
260                 }
261                 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
262                 free(btemp);
263             }
264         }
265         MPI_Barrier( comm );
266     }
267
268     MTest_Finalize( errs );
269     MPI_Finalize();
270     return 0;
271 }