1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
3 * (C) 2003 by Argonne National Laboratory.
4 * See COPYRIGHT in top-level directory.
12 static char MTEST_Descrip[] = "Test of various send cancel calls";
15 int main( int argc, char *argv[] )
18 int rank, size, source, dest;
22 static int bufsizes[4] = { 1, 100, 10000, 1000000 };
26 int veryPicky = 0; /* Set to 1 to test "quality of implementation" in
27 a tricky part of cancel */
30 MTest_Init( &argc, &argv );
32 comm = MPI_COMM_WORLD;
33 MPI_Comm_rank( comm, &rank );
34 MPI_Comm_size( comm, &size );
39 for (cs=0; cs<4; cs++) {
42 buf = (char *)malloc( n );
44 fprintf( stderr, "Unable to allocate %d bytes\n", n );
45 MPI_Abort( MPI_COMM_WORLD, 1 );
47 MPI_Send_init( buf, n, MPI_CHAR, dest, cs+n+1, comm, &req );
50 MPI_Wait( &req, &status );
51 MPI_Test_cancelled( &status, &flag );
54 printf( "Failed to cancel a persistent send request\n" );
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 */
66 MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
69 else if (rank == dest)
73 MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
74 MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
77 /* If the message was not cancelled, receive it here */
78 btemp = (char*)malloc( nn );
81 fprintf( stderr, "Unable to allocate %d bytes\n", nn);
82 MPI_Abort( MPI_COMM_WORLD, 1 );
84 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
95 buf = (char *)malloc( n );
97 fprintf( stderr, "Unable to allocate %d bytes\n", n );
98 MPI_Abort( MPI_COMM_WORLD, 1 );
100 bsendbufsize = n + MPI_BSEND_OVERHEAD;
101 bsendbuf = (char *)malloc( bsendbufsize );
103 fprintf( stderr, "Unable to allocate %d bytes for bsend\n", n );
104 MPI_Abort( MPI_COMM_WORLD, 1 );
106 MPI_Buffer_attach( bsendbuf, bsendbufsize );
107 MPI_Bsend_init( buf, n, MPI_CHAR, dest, cs+n+2, comm, &req );
110 MPI_Wait( &req, &status );
111 MPI_Test_cancelled( &status, &flag );
114 printf( "Failed to cancel a persistent bsend request\n" );
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 */
126 MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
128 MPI_Buffer_detach( &bf, &bs );
131 else if (rank == dest)
135 MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
136 MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
139 /* If the message was not cancelled, receive it here */
140 btemp = (char*)malloc( nn );
143 fprintf( stderr, "Unable to allocate %d bytes\n", nn);
144 MPI_Abort( MPI_COMM_WORLD, 1 );
146 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
152 /* Because this test is erroneous, we do not perform it unless
153 TEST_IRSEND is defined. */
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) {
161 buf = (char *)malloc( n );
163 fprintf( stderr, "Unable to allocate %d bytes\n", n );
164 MPI_Abort( MPI_COMM_WORLD, 1 );
166 MPI_Rsend_init( buf, n, MPI_CHAR, dest, cs+n+3, comm, &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
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) {
182 printf( "Failed to cancel a persistent rsend request\n" );
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 */
194 MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
197 else if (rank == dest)
201 MPI_Recv( &n, 1, MPI_INT, 0, 123, comm, &status );
202 MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
205 /* If the message was not cancelled, receive it here */
206 btemp = (char*)malloc( n );
209 fprintf( stderr, "Unable to allocate %d bytes\n", n);
210 MPI_Abort( MPI_COMM_WORLD, 1 );
212 MPI_Recv( btemp, n, MPI_CHAR, 0, tag, comm, &status );
221 buf = (char *)malloc( n );
223 fprintf( stderr, "Unable to allocate %d bytes\n", n );
224 MPI_Abort( MPI_COMM_WORLD, 1 );
226 MPI_Ssend_init( buf, n, MPI_CHAR, dest, cs+n+4, comm, &req );
229 MPI_Wait( &req, &status );
230 MPI_Test_cancelled( &status, &flag );
233 printf( "Failed to cancel a persistent ssend request\n" );
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 */
245 MPI_Send( &n, 1, MPI_INT, dest, 123, comm );
248 else if (rank == dest)
252 MPI_Recv( &nn, 1, MPI_INT, 0, 123, comm, &status );
253 MPI_Recv( &tag, 1, MPI_INT, 0, 123, comm, &status );
256 /* If the message was not cancelled, receive it here */
257 btemp = (char*)malloc( nn );
260 fprintf( stderr, "Unable to allocate %d bytes\n", nn);
261 MPI_Abort( MPI_COMM_WORLD, 1 );
263 MPI_Recv( btemp, nn, MPI_CHAR, 0, tag, comm, &status );
270 MTest_Finalize( errs );