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
|
/*
* Copyright (C) by Argonne National Laboratory
* See COPYRIGHT in top-level directory
*/
#include "mpi.h"
#include <stdio.h>
#include <string.h>
#include "mpitest.h"
/*
static char MTest_descrip[] = "Receive partial datatypes and check that\
MPI_Getelements gives the correct version";
*/
int main(int argc, char *argv[])
{
int errs = 0;
MPI_Datatype outtype, oldtypes[2];
MPI_Aint offsets[2];
int blklens[2];
MPI_Comm comm;
int size, rank, src, dest, tag;
MTest_Init(&argc, &argv);
comm = MPI_COMM_WORLD;
MPI_Comm_rank(comm, &rank);
MPI_Comm_size(comm, &size);
if (size < 2) {
errs++;
printf("This test requires at least 2 processes\n");
MPI_Abort(MPI_COMM_WORLD, 1);
}
src = 0;
dest = 1;
if (rank == src) {
int buf[128], position, cnt;
MTEST_VG_MEM_INIT(buf, 128 * sizeof(buf[0]));
/* sender */
/* Create a datatype and send it (multiple of sizeof(int)) */
/* Create a send struct type */
oldtypes[0] = MPI_INT;
oldtypes[1] = MPI_CHAR;
blklens[0] = 1;
blklens[1] = 4 * sizeof(int);
offsets[0] = 0;
offsets[1] = sizeof(int);
MPI_Type_create_struct(2, blklens, offsets, oldtypes, &outtype);
MPI_Type_commit(&outtype);
buf[0] = 4 * sizeof(int);
/* printf("About to send to %d\n", dest); */
MPI_Send(buf, 1, outtype, dest, 0, comm);
MPI_Type_free(&outtype);
/* Create a datatype and send it (not a multiple of sizeof(int)) */
/* Create a send struct type */
oldtypes[0] = MPI_INT;
oldtypes[1] = MPI_CHAR;
blklens[0] = 1;
blklens[1] = 4 * sizeof(int) + 1;
offsets[0] = 0;
offsets[1] = sizeof(int);
MPI_Type_create_struct(2, blklens, offsets, oldtypes, &outtype);
MPI_Type_commit(&outtype);
buf[0] = 4 * sizeof(int) + 1;
MPI_Send(buf, 1, outtype, dest, 1, comm);
MPI_Type_free(&outtype);
/* Pack data and send as packed */
position = 0;
cnt = 7;
MPI_Pack(&cnt, 1, MPI_INT, buf, 128 * sizeof(int), &position, comm);
MPI_Pack((void *) "message", 7, MPI_CHAR, buf, 128 * sizeof(int), &position, comm);
MPI_Send(buf, position, MPI_PACKED, dest, 2, comm);
} else if (rank == dest) {
MPI_Status status;
int buf[128], i, elms, count;
/* Receiver */
/* Create a receive struct type */
oldtypes[0] = MPI_INT;
oldtypes[1] = MPI_CHAR;
blklens[0] = 1;
blklens[1] = 256;
offsets[0] = 0;
offsets[1] = sizeof(int);
MPI_Type_create_struct(2, blklens, offsets, oldtypes, &outtype);
MPI_Type_commit(&outtype);
for (i = 0; i < 3; i++) {
tag = i;
/* printf("about to receive tag %d from %d\n", i, src); */
MPI_Recv(buf, 1, outtype, src, tag, comm, &status);
MPI_Get_elements(&status, outtype, &elms);
if (elms != buf[0] + 1) {
errs++;
printf("For test %d, Get elements gave %d but should be %d\n", i, elms, buf[0] + 1);
}
MPI_Get_count(&status, outtype, &count);
if (count != MPI_UNDEFINED) {
errs++;
printf("For partial send, Get_count did not return MPI_UNDEFINED\n");
}
}
MPI_Type_free(&outtype);
}
MTest_Finalize(errs);
return MTestReturnValue(errs);
}
|