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